Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecInelasticModel Class Reference

#include <G4MicroElecInelasticModel.hh>

+ Inheritance diagram for G4MicroElecInelasticModel:

Public Member Functions

 G4MicroElecInelasticModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="MicroElecInelasticModel")
 
virtual ~G4MicroElecInelasticModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
G4double TransferedEnergy (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
 
void SelectFasterComputation (G4bool input)
 
G4MicroElecInelasticModeloperator= (const G4MicroElecInelasticModel &right)=delete
 
 G4MicroElecInelasticModel (const G4MicroElecInelasticModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 58 of file G4MicroElecInelasticModel.hh.

Constructor & Destructor Documentation

◆ G4MicroElecInelasticModel() [1/2]

G4MicroElecInelasticModel::G4MicroElecInelasticModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "MicroElecInelasticModel" )

Definition at line 69 of file G4MicroElecInelasticModel.cc.

71:G4VEmModel(nam),isInitialised(false)
72{
74
75 verboseLevel= 0;
76 // Verbosity scale:
77 // 0 = nothing
78 // 1 = warning for energy non-conservation
79 // 2 = details of energy budget
80 // 3 = calculation of cross sections, file openings, sampling of atoms
81 // 4 = entering in methods
82
83 if( verboseLevel>0 )
84 {
85 G4cout << "MicroElec inelastic model is constructed " << G4endl;
86 }
87
88 //Mark this model as "applicable" for atomic deexcitation
90 fAtomDeexcitation = nullptr;
92
93 // default generator
95
96 // Selection of computation method
97 fasterCode = true; //false;
98}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

◆ ~G4MicroElecInelasticModel()

G4MicroElecInelasticModel::~G4MicroElecInelasticModel ( )
virtual

Definition at line 102 of file G4MicroElecInelasticModel.cc.

103{
104 // Cross section
105 for (auto & pos : tableData)
106 {
107 G4MicroElecCrossSectionDataSet* table = pos.second;
108 delete table;
109 }
110
111 // Final state
112 eVecm.clear();
113 pVecm.clear();
114
115}

◆ G4MicroElecInelasticModel() [2/2]

G4MicroElecInelasticModel::G4MicroElecInelasticModel ( const G4MicroElecInelasticModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4MicroElecInelasticModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 307 of file G4MicroElecInelasticModel.cc.

312{
313 if (verboseLevel > 3)
314 G4cout << "Calling CrossSectionPerVolume() of G4MicroElecInelasticModel" << G4endl;
315
316 G4double density = material->GetTotNbOfAtomsPerVolume();
317
318 // Calculate total cross section for model
319 G4double lowLim = 0;
320 G4double highLim = 0;
321 G4double sigma=0;
322
323 const G4String& particleName = particleDefinition->GetParticleName();
324 G4String nameLocal = particleName ;
325
326 G4double Zeff2 = 1.0;
327 G4double Mion_c2 = particleDefinition->GetPDGMass();
328
329 if (Mion_c2 > proton_mass_c2)
330 {
331 G4ionEffectiveCharge EffCharge ;
332 G4double Zeff = EffCharge.EffectiveCharge(particleDefinition, material,ekin);
333 Zeff2 = Zeff*Zeff;
334
335 if (verboseLevel > 3)
336 G4cout << "Before scaling : " << G4endl
337 << "Particle : " << nameLocal << ", mass : " << Mion_c2/proton_mass_c2 << "*mp, charge " << Zeff
338 << ", Ekin (eV) = " << ekin/eV << G4endl ;
339
340 ekin *= proton_mass_c2/Mion_c2 ;
341 nameLocal = "proton" ;
342
343 if (verboseLevel > 3)
344 G4cout << "After scaling : " << G4endl
345 << "Particle : " << nameLocal << ", Ekin (eV) = " << ekin/eV << G4endl ;
346 }
347
348 if (material == nistSi || material->GetBaseMaterial() == nistSi)
349 {
350 auto pos1 = lowEnergyLimit.find(nameLocal);
351 if (pos1 != lowEnergyLimit.end())
352 {
353 lowLim = pos1->second;
354 }
355
356 auto pos2 = highEnergyLimit.find(nameLocal);
357 if (pos2 != highEnergyLimit.end())
358 {
359 highLim = pos2->second;
360 }
361
362 if (ekin >= lowLim && ekin < highLim)
363 {
364 auto pos = tableData.find(nameLocal);
365 if (pos != tableData.end())
366 {
367 G4MicroElecCrossSectionDataSet* table = pos->second;
368 if (table != nullptr)
369 {
370 sigma = table->FindValue(ekin);
371 }
372 }
373 else
374 {
375 G4Exception("G4MicroElecInelasticModel::CrossSectionPerVolume","em0002",
376 FatalException,"Model not applicable to particle type.");
377 }
378 }
379 else
380 {
381 if (nameLocal!="e-")
382 {
383 // G4cout << "Particle : " << nameLocal << ", Ekin (eV) = " << ekin/eV << G4endl;
384 // G4cout << "### Warning: particle energy out of bounds! ###" << G4endl;
385 }
386 }
387
388 if (verboseLevel > 3)
389 {
390 G4cout << "---> Kinetic energy (eV)=" << ekin/eV << G4endl;
391 G4cout << " - Cross section per Si atom (cm^2)=" << sigma*Zeff2/cm2 << G4endl;
392 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density*Zeff2/(1./cm) << G4endl;
393 }
394
395 } // if (SiMaterial)
396 return sigma*density*Zeff2;
397}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
const G4Material * GetBaseMaterial() const
G4double GetTotNbOfAtomsPerVolume() const
G4double FindValue(G4double e, G4int componentId=0) const override
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)

◆ DifferentialCrossSection()

G4double G4MicroElecInelasticModel::DifferentialCrossSection ( G4ParticleDefinition * aParticleDefinition,
G4double k,
G4double energyTransfer,
G4int shell )

Definition at line 659 of file G4MicroElecInelasticModel.cc.

663{
664 G4double sigma = 0.;
665
666 if (energyTransfer >= SiStructure.Energy(LevelIndex))
667 {
668 G4double valueT1 = 0;
669 G4double valueT2 = 0;
670 G4double valueE21 = 0;
671 G4double valueE22 = 0;
672 G4double valueE12 = 0;
673 G4double valueE11 = 0;
674
675 G4double xs11 = 0;
676 G4double xs12 = 0;
677 G4double xs21 = 0;
678 G4double xs22 = 0;
679
680 if (particleDefinition == G4Electron::ElectronDefinition())
681 {
682 // k should be in eV and energy transfer eV also
683 auto t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
684 auto t1 = t2-1;
685 // The following condition avoids situations where energyTransfer >last vector element
686 if (energyTransfer <= eVecm[(*t1)].back() && energyTransfer <= eVecm[(*t2)].back() )
687 {
688 auto e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
689 auto e11 = e12-1;
690 auto e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
691 auto e21 = e22-1;
692
693 valueT1 =*t1;
694 valueT2 =*t2;
695 valueE21 =*e21;
696 valueE22 =*e22;
697 valueE12 =*e12;
698 valueE11 =*e11;
699
700 xs11 = eDiffCrossSectionData[LevelIndex][valueT1][valueE11];
701 xs12 = eDiffCrossSectionData[LevelIndex][valueT1][valueE12];
702 xs21 = eDiffCrossSectionData[LevelIndex][valueT2][valueE21];
703 xs22 = eDiffCrossSectionData[LevelIndex][valueT2][valueE22];
704 }
705 }
706
707 if (particleDefinition == G4Proton::ProtonDefinition())
708 {
709 // k should be in eV and energy transfer eV also
710 auto t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
711 auto t1 = t2-1;
712 if (energyTransfer <= pVecm[(*t1)].back() && energyTransfer <= pVecm[(*t2)].back() )
713 {
714 auto e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
715 auto e11 = e12-1;
716
717 auto e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer);
718 auto e21 = e22-1;
719
720 valueT1 =*t1;
721 valueT2 =*t2;
722 valueE21 =*e21;
723 valueE22 =*e22;
724 valueE12 =*e12;
725 valueE11 =*e11;
726
727 xs11 = pDiffCrossSectionData[LevelIndex][valueT1][valueE11];
728 xs12 = pDiffCrossSectionData[LevelIndex][valueT1][valueE12];
729 xs21 = pDiffCrossSectionData[LevelIndex][valueT2][valueE21];
730 xs22 = pDiffCrossSectionData[LevelIndex][valueT2][valueE22];
731 }
732 }
733
734 sigma = QuadInterpolator( valueE11, valueE12,
735 valueE21, valueE22,
736 xs11, xs12,
737 xs21, xs22,
738 valueT1, valueT2,
739 k, energyTransfer);
740 }
741 return sigma;
742}
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
G4double Energy(G4int level)
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85

◆ Initialise()

void G4MicroElecInelasticModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 119 of file G4MicroElecInelasticModel.cc.

121{
122
123 if (verboseLevel > 3)
124 G4cout << "Calling G4MicroElecInelasticModel::Initialise()" << G4endl;
125
126 // Energy limits
127
128 G4String fileElectron("microelec/sigma_inelastic_e_Si");
129 G4String fileProton("microelec/sigma_inelastic_p_Si");
130
133 G4String electron = electronDef->GetParticleName();
134 G4String proton = protonDef->GetParticleName();;
135
136 G4double scaleFactor = 1e-18 * cm *cm;
137
138 const char* path = G4FindDataDir("G4LEDATA");
139
140 // *** ELECTRON
141 tableFile[electron] = fileElectron;
142 lowEnergyLimit[electron] = 16.7 * eV;
143 highEnergyLimit[electron] = 100.0 * MeV;
144
145 // Cross section
147 tableE->LoadData(fileElectron);
148
149 tableData[electron] = tableE;
150
151 // Final state
152
153 std::ostringstream eFullFileName;
154
155 if (fasterCode) eFullFileName << path << "/microelec/sigmadiff_cumulated_inelastic_e_Si.dat";
156 else eFullFileName << path << "/microelec/sigmadiff_inelastic_e_Si.dat";
157
158 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
159
160 if (!eDiffCrossSection)
161 {
162 if (fasterCode) G4Exception("G4MicroElecInelasticModel::Initialise","em0003",
163 FatalException,"Missing data file:/microelec/sigmadiff_cumulated_inelastic_e_Si.dat");
164
165 else G4Exception("G4MicroElecInelasticModel::Initialise","em0003",
166 FatalException,"Missing data file:/microelec/sigmadiff_inelastic_e_Si.dat");
167 }
168
169 // Clear the arrays for re-initialization case (MT mode)
170 // Octobre 22nd, 2014 - Melanie Raine
171 eTdummyVec.clear();
172 pTdummyVec.clear();
173
174 eVecm.clear();
175 pVecm.clear();
176
177 for (int j=0; j<6; j++)
178 {
179 eProbaShellMap[j].clear();
180 pProbaShellMap[j].clear();
181
182 eDiffCrossSectionData[j].clear();
183 pDiffCrossSectionData[j].clear();
184
185 eNrjTransfData[j].clear();
186 pNrjTransfData[j].clear();
187 }
188
189 //
190 eTdummyVec.push_back(0.);
191 while(!eDiffCrossSection.eof())
192 {
193 G4double tDummy;
194 G4double eDummy;
195 eDiffCrossSection>>tDummy>>eDummy;
196 if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
197
198 G4double tmp;
199 for (int j=0; j<6; j++)
200 {
201 eDiffCrossSection>> tmp;
202
203 eDiffCrossSectionData[j][tDummy][eDummy] = tmp;
204
205 if (fasterCode)
206 {
207 eNrjTransfData[j][tDummy][eDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
208 eProbaShellMap[j][tDummy].push_back(eDiffCrossSectionData[j][tDummy][eDummy]);
209 }
210 else
211 {
212 // SI - only if eof is not reached !
213 if (!eDiffCrossSection.eof()) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
214 eVecm[tDummy].push_back(eDummy);
215 }
216
217 }
218 }
219 //
220
221 // *** PROTON
222 tableFile[proton] = fileProton;
223 lowEnergyLimit[proton] = 50. * keV;
224 highEnergyLimit[proton] = 10. * GeV;
225
226 // Cross section
228 tableP->LoadData(fileProton);
229 tableData[proton] = tableP;
230
231 // Final state
232 std::ostringstream pFullFileName;
233
234 if (fasterCode) pFullFileName << path << "/microelec/sigmadiff_cumulated_inelastic_p_Si.dat";
235 else pFullFileName << path << "/microelec/sigmadiff_inelastic_p_Si.dat";
236
237 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
238
239 if (!pDiffCrossSection)
240 {
241 if (fasterCode) G4Exception("G4MicroElecInelasticModel::Initialise","em0003",
242 FatalException,"Missing data file:/microelec/sigmadiff_cumulated_inelastic_p_Si.dat");
243
244 else G4Exception("G4MicroElecInelasticModel::Initialise","em0003",
245 FatalException,"Missing data file:/microelec/sigmadiff_inelastic_p_Si.dat");
246 }
247
248 pTdummyVec.push_back(0.);
249 while(!pDiffCrossSection.eof())
250 {
251 G4double tDummy;
252 G4double eDummy;
253 pDiffCrossSection>>tDummy>>eDummy;
254 if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
255 for (int j=0; j<6; j++)
256 {
257 pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
258
259 if (fasterCode)
260 {
261 pNrjTransfData[j][tDummy][pDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
262 pProbaShellMap[j][tDummy].push_back(pDiffCrossSectionData[j][tDummy][eDummy]);
263 }
264 else
265 {
266 // SI - only if eof is not reached !
267 if (!pDiffCrossSection.eof()) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
268 pVecm[tDummy].push_back(eDummy);
269 }
270 }
271 }
272
273 if (particle==electronDef)
274 {
275 SetLowEnergyLimit(lowEnergyLimit[electron]);
276 SetHighEnergyLimit(highEnergyLimit[electron]);
277 }
278
279 if (particle==protonDef)
280 {
281 SetLowEnergyLimit(lowEnergyLimit[proton]);
282 SetHighEnergyLimit(highEnergyLimit[proton]);
283 }
284
285 if( verboseLevel>0 )
286 {
287 G4cout << "MicroElec Inelastic model is initialized " << G4endl
288 << "Energy range: "
289 << LowEnergyLimit() / keV << " keV - "
290 << HighEnergyLimit() / MeV << " MeV for "
291 << particle->GetParticleName()
292 << " with mass (amu) " << particle->GetPDGMass()/proton_mass_c2
293 << " and charge " << particle->GetPDGCharge()
294 << G4endl << G4endl ;
295 }
296
297 //
298 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
299
300 if (isInitialised) { return; }
302 isInitialised = true;
303}
const char * G4FindDataDir(const char *)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4bool LoadData(const G4String &argFileName) override
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)

◆ operator=()

G4MicroElecInelasticModel & G4MicroElecInelasticModel::operator= ( const G4MicroElecInelasticModel & right)
delete

◆ SampleSecondaries()

void G4MicroElecInelasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * fvect,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * particle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 401 of file G4MicroElecInelasticModel.cc.

406{
407
408 if (verboseLevel > 3)
409 G4cout << "Calling SampleSecondaries() of G4MicroElecInelasticModel" << G4endl;
410
411 G4double lowLim = 0;
412 G4double highLim = 0;
413
414 G4double ekin = particle->GetKineticEnergy();
415 G4double k = ekin ;
416
417 G4ParticleDefinition* PartDef = particle->GetDefinition();
418 const G4String& particleName = PartDef->GetParticleName();
419 G4String nameLocal2 = particleName ;
420 G4double particleMass = particle->GetDefinition()->GetPDGMass();
421
422 if (particleMass > proton_mass_c2)
423 {
424 k *= proton_mass_c2/particleMass ;
425 PartDef = G4Proton::ProtonDefinition();
426 nameLocal2 = "proton" ;
427 }
428
429 auto pos1 = lowEnergyLimit.find(nameLocal2);
430 if (pos1 != lowEnergyLimit.end())
431 {
432 lowLim = pos1->second;
433 }
434
435 auto pos2 = highEnergyLimit.find(nameLocal2);
436 if (pos2 != highEnergyLimit.end())
437 {
438 highLim = pos2->second;
439 }
440
441 if (k >= lowLim && k < highLim)
442 {
443 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
444 G4double totalEnergy = ekin + particleMass;
445 G4double pSquare = ekin * (totalEnergy + particleMass);
446 G4double totalMomentum = std::sqrt(pSquare);
447
448 G4int Shell = 0;
449
450 /* if (!fasterCode)*/ Shell = RandomSelect(k,nameLocal2);
451
452 // SI: The following protection is necessary to avoid infinite loops :
453 // sigmadiff_ionisation_e_born.dat has non zero partial xs at 18 eV for shell 3 (ionizationShell ==2)
454 // sigmadiff_cumulated_ionisation_e_born.dat has zero cumulated partial xs at 18 eV for shell 3 (ionizationShell ==2)
455 // this is due to the fact that the max allowed transfered energy is (18+10.79)/2=17.025 eV and only transfered energies
456 // strictly above this value have non zero partial xs in sigmadiff_ionisation_e_born.dat (starting at trans = 17.12 eV)
457
458 G4double bindingEnergy = SiStructure.Energy(Shell);
459
460 if (verboseLevel > 3)
461 {
462 G4cout << "---> Kinetic energy (eV)=" << k/eV << G4endl ;
463 G4cout << "Shell: " << Shell << ", energy: " << bindingEnergy/eV << G4endl;
464 }
465
466 // sample deexcitation
467 std::size_t secNumberInit = 0; // need to know at a certain point the energy of secondaries
468 std::size_t secNumberFinal = 0; // So I'll make the difference and then sum the energies
469
470 //SI: additional protection if tcs interpolation method is modified
471 if (k<bindingEnergy) return;
472
473 G4int Z = 14;
474
475 if(fAtomDeexcitation && Shell > 2) {
476
478
479 if (Shell == 4)
480 {
482 }
483 else if (Shell == 3)
484 {
486 }
487
488 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
489 secNumberInit = fvect->size();
490 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
491 secNumberFinal = fvect->size();
492 }
493
494 G4double secondaryKinetic=-1000*eV;
495
496 if (!fasterCode)
497 {
498 secondaryKinetic = RandomizeEjectedElectronEnergy(PartDef,k,Shell);
499 }
500 else
501 {
502 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(PartDef,k,Shell);
503 }
504
505 if (verboseLevel > 3)
506 {
507 G4cout << "Ionisation process" << G4endl;
508 G4cout << "Shell: " << Shell << " Kin. energy (eV)=" << k/eV
509 << " Sec. energy (eV)=" << secondaryKinetic/eV << G4endl;
510 }
511
512 G4ThreeVector deltaDirection =
513 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
514 Z, Shell,
515 couple->GetMaterial());
516
518 {
519 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
520 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
521 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
522 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
523 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
524 finalPx /= finalMomentum;
525 finalPy /= finalMomentum;
526 finalPz /= finalMomentum;
527
528 G4ThreeVector direction;
529 direction.set(finalPx,finalPy,finalPz);
530
532 }
533 else
535
536 // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
537 G4double deexSecEnergy = 0;
538 for (std::size_t j=secNumberInit; j < secNumberFinal; ++j) {
539 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();}
540
541 fParticleChangeForGamma->SetProposedKineticEnergy(ekin-bindingEnergy-secondaryKinetic);
542 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy-deexSecEnergy);
543
544 if (secondaryKinetic>0)
545 {
546 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
547 fvect->push_back(dp);
548 }
549 }
550}
int G4int
Definition G4Types.hh:85
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4VEmAngularDistribution * GetAngularDistribution()
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectFasterComputation()

void G4MicroElecInelasticModel::SelectFasterComputation ( G4bool input)
inline

Definition at line 152 of file G4MicroElecInelasticModel.hh.

153{
154 fasterCode = input;
155}

◆ TransferedEnergy()

G4double G4MicroElecInelasticModel::TransferedEnergy ( G4ParticleDefinition * aParticleDefinition,
G4double incomingParticleEnergy,
G4int shell,
G4double random )

Definition at line 876 of file G4MicroElecInelasticModel.cc.

880{
881 G4double nrj = 0.;
882
883 G4double valueK1 = 0;
884 G4double valueK2 = 0;
885 G4double valuePROB21 = 0;
886 G4double valuePROB22 = 0;
887 G4double valuePROB12 = 0;
888 G4double valuePROB11 = 0;
889
890 G4double nrjTransf11 = 0;
891 G4double nrjTransf12 = 0;
892 G4double nrjTransf21 = 0;
893 G4double nrjTransf22 = 0;
894
895 G4double maximumEnergyTransfer1 = 0;
896 G4double maximumEnergyTransfer2 = 0;
897 G4double maximumEnergyTransferP = 4.* (electron_mass_c2 / proton_mass_c2) * k;
898 G4double bindingEnergy = SiStructure.Energy(ionizationLevelIndex)*1e6;
899
900 if (particleDefinition == G4Electron::ElectronDefinition())
901 {
902 // k should be in eV
903 auto k2 = std::upper_bound(eTdummyVec.begin(),
904 eTdummyVec.end(),
905 k);
906 auto k1 = k2 - 1;
907
908 /*
909 G4cout << "----> k=" << k
910 << " " << *k1
911 << " " << *k2
912 << " " << random
913 << " " << ionizationLevelIndex
914 << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
915 << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
916 << G4endl;
917 */
918
919 // SI : the following condition avoids situations where random >last vector element
920 if (random <= eProbaShellMap[ionizationLevelIndex][(*k1)].back()
921 && random <= eProbaShellMap[ionizationLevelIndex][(*k2)].back())
922 {
923 auto prob12 =
924 std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
925 eProbaShellMap[ionizationLevelIndex][(*k1)].end(),
926 random);
927 auto prob11 = prob12 - 1;
928 auto prob22 =
929 std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
930 eProbaShellMap[ionizationLevelIndex][(*k2)].end(),
931 random);
932 auto prob21 = prob22 - 1;
933
934 valueK1 = *k1;
935 valueK2 = *k2;
936 valuePROB21 = *prob21;
937 valuePROB22 = *prob22;
938 valuePROB12 = *prob12;
939 valuePROB11 = *prob11;
940
941 // The following condition avoid getting transfered energy < binding energy and forces cumxs = 1 for maximum energy transfer.
942 if(valuePROB11 == 0) nrjTransf11 = bindingEnergy;
943 else nrjTransf11 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
944 if(valuePROB12 == 1)
945 {
946 if ((valueK1+bindingEnergy)/2. > valueK1) maximumEnergyTransfer1=valueK1;
947 else maximumEnergyTransfer1 = (valueK1+bindingEnergy)/2.;
948
949 nrjTransf12 = maximumEnergyTransfer1;
950 }
951 else nrjTransf12 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
952
953 if(valuePROB21 == 0) nrjTransf21 = bindingEnergy;
954 else nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
955 if(valuePROB22 == 1)
956 {
957 if ((valueK2+bindingEnergy)/2. > valueK2) maximumEnergyTransfer2=valueK2;
958 else maximumEnergyTransfer2 = (valueK2+bindingEnergy)/2.;
959
960 nrjTransf22 = maximumEnergyTransfer2;
961 }
962 else
963 nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
964 }
965 // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
966 if (random > eProbaShellMap[ionizationLevelIndex][(*k1)].back())
967 {
968 auto prob22 =
969 std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
970 eProbaShellMap[ionizationLevelIndex][(*k2)].end(),
971 random);
972
973 auto prob21 = prob22 - 1;
974
975 valueK1 = *k1;
976 valueK2 = *k2;
977 valuePROB21 = *prob21;
978 valuePROB22 = *prob22;
979
980 nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
981 nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
982
983 G4double interpolatedvalue2 = Interpolate(valuePROB21,
984 valuePROB22,
985 random,
986 nrjTransf21,
987 nrjTransf22);
988
989 // zeros are explicitly set
990 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
991 return value;
992 }
993 }
994 //
995 else if (particleDefinition == G4Proton::ProtonDefinition())
996 {
997 // k should be in eV
998 auto k2 = std::upper_bound(pTdummyVec.begin(),
999 pTdummyVec.end(),
1000 k);
1001
1002 auto k1 = k2 - 1;
1003
1004 /*
1005 G4cout << "----> k=" << k
1006 << " " << *k1
1007 << " " << *k2
1008 << " " << random
1009 << " " << ionizationLevelIndex
1010 << " " << pProbaShellMap[ionizationLevelIndex][(*k1)].back()
1011 << " " << pProbaShellMap[ionizationLevelIndex][(*k2)].back()
1012 << G4endl;
1013 */
1014
1015 // SI : the following condition avoids situations where random > last vector element,
1016 // for eg. when the last element is zero
1017 if (random <= pProbaShellMap[ionizationLevelIndex][(*k1)].back()
1018 && random <= pProbaShellMap[ionizationLevelIndex][(*k2)].back())
1019 {
1020 auto prob12 =
1021 std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
1022 pProbaShellMap[ionizationLevelIndex][(*k1)].end(),
1023 random);
1024
1025 auto prob11 = prob12 - 1;
1026
1027 auto prob22 =
1028 std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1029 pProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1030 random);
1031
1032 auto prob21 = prob22 - 1;
1033
1034 valueK1 = *k1;
1035 valueK2 = *k2;
1036 valuePROB21 = *prob21;
1037 valuePROB22 = *prob22;
1038 valuePROB12 = *prob12;
1039 valuePROB11 = *prob11;
1040
1041 // The following condition avoid getting transfered energy < binding energy and forces cumxs = 1 for maximum energy transfer.
1042 if(valuePROB11 == 0) nrjTransf11 = bindingEnergy;
1043 else nrjTransf11 = pNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
1044 if(valuePROB12 == 1) nrjTransf12 = maximumEnergyTransferP;
1045 else nrjTransf12 = pNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
1046 if(valuePROB21 == 0) nrjTransf21 = bindingEnergy;
1047 else nrjTransf21 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1048 if(valuePROB22 == 1) nrjTransf22 = maximumEnergyTransferP;
1049 else nrjTransf22 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1050
1051 }
1052
1053 // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
1054 if (random > pProbaShellMap[ionizationLevelIndex][(*k1)].back())
1055 {
1056 auto prob22 =
1057 std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1058 pProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1059 random);
1060
1061 auto prob21 = prob22 - 1;
1062
1063 valueK1 = *k1;
1064 valueK2 = *k2;
1065 valuePROB21 = *prob21;
1066 valuePROB22 = *prob22;
1067
1068 nrjTransf21 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1069 nrjTransf22 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1070
1071 G4double interpolatedvalue2 = Interpolate(valuePROB21,
1072 valuePROB22,
1073 random,
1074 nrjTransf21,
1075 nrjTransf22);
1076
1077 // zeros are explicitly set
1078 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
1079
1080 return value;
1081 }
1082 }
1083 // End electron and proton cases
1084
1085 G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21
1086 * nrjTransf22;
1087
1088 if (nrjTransfProduct != 0.)
1089 {
1090 nrj = QuadInterpolator(valuePROB11,
1091 valuePROB12,
1092 valuePROB21,
1093 valuePROB22,
1094 nrjTransf11,
1095 nrjTransf12,
1096 nrjTransf21,
1097 nrjTransf22,
1098 valueK1,
1099 valueK2,
1100 k,
1101 random);
1102 }
1103
1104 return nrj;
1105}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MicroElecInelasticModel::fParticleChangeForGamma
protected

The documentation for this class was generated from the following files: