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

#include <G4MicroElecElasticModel_new.hh>

+ Inheritance diagram for G4MicroElecElasticModel_new:

Public Member Functions

 G4MicroElecElasticModel_new (const G4ParticleDefinition *p=0, const G4String &nam="MicroElecElasticModel")
 
 ~G4MicroElecElasticModel_new () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
G4double AcousticCrossSectionPerVolume (G4double ekin, G4double kbz, G4double rho, G4double cs, G4double Aac, G4double Eac, G4double prefactor)
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
G4double DamageEnergy (G4double T, G4double A, G4double Z)
 
- 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 90 of file G4MicroElecElasticModel_new.hh.

Constructor & Destructor Documentation

◆ G4MicroElecElasticModel_new()

G4MicroElecElasticModel_new::G4MicroElecElasticModel_new ( const G4ParticleDefinition * p = 0,
const G4String & nam = "MicroElecElasticModel" )

Definition at line 89 of file G4MicroElecElasticModel_new.cc.

91 :G4VEmModel(nam), isInitialised(false)
92{
93 killBelowEnergy = 0.1*eV; // Minimum e- energy for energy loss by excitation
94 lowEnergyLimit = 0.1 * eV;
95 lowEnergyLimitOfModel = 10 * eV; // The model lower energy is 10 eV
96 highEnergyLimit = 500. * keV;
97 SetLowEnergyLimit(lowEnergyLimit);
98 SetHighEnergyLimit(highEnergyLimit);
99
100 verboseLevel= 0;
101 // Verbosity scale:
102 // 0 = nothing
103 // 1 = warning for energy non-conservation
104 // 2 = details of energy budget
105 // 3 = calculation of cross sections, file openings, sampling of atoms
106 // 4 = entering in methods
107
108 if( verboseLevel>0 )
109 {
110 G4cout << "MicroElec Elastic model is constructed " << G4endl
111 << "Energy range: "
112 << lowEnergyLimit / eV << " eV - "
113 << highEnergyLimit / MeV << " MeV"
114 << G4endl;
115 }
117
118 killElectron = false;
119 acousticModelEnabled = false;
120 currentMaterialName = "";
121 isOkToBeInitialised = false;
122}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4MicroElecElasticModel_new()

G4MicroElecElasticModel_new::~G4MicroElecElasticModel_new ( )
override

Definition at line 126 of file G4MicroElecElasticModel_new.cc.

127{
128 // For total cross section
129 TCSMap::iterator pos2;
130 for (pos2 = tableTCS.begin(); pos2 != tableTCS.end(); ++pos2) {
131 MapData* tableData = pos2->second;
132 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
133 for (pos = tableData->begin(); pos != tableData->end(); ++pos)
134 {
136 delete table;
137 }
138 delete tableData;
139 }
140
141 //Clearing DCS maps
142
143 ThetaMap::iterator iterator_angle;
144 for (iterator_angle = thetaDataStorage.begin(); iterator_angle != thetaDataStorage.end(); ++iterator_angle) {
145 TriDimensionMap* eDiffCrossSectionData = iterator_angle->second;
146 eDiffCrossSectionData->clear();
147 delete eDiffCrossSectionData;
148 }
149
150 energyMap::iterator iterator_energy;
151 for (iterator_energy = eIncidentEnergyStorage.begin(); iterator_energy != eIncidentEnergyStorage.end(); ++iterator_energy) {
152 std::vector<G4double>* eTdummyVec = iterator_energy->second;
153 eTdummyVec->clear();
154 delete eTdummyVec;
155 }
156
157 ProbaMap::iterator iterator_proba;
158 for (iterator_proba = eProbaStorage.begin(); iterator_proba != eProbaStorage.end(); ++iterator_proba) {
159 VecMap* eVecm = iterator_proba->second;
160 eVecm->clear();
161 delete eVecm;
162 }
163
164}

Member Function Documentation

◆ AcousticCrossSectionPerVolume()

G4double G4MicroElecElasticModel_new::AcousticCrossSectionPerVolume ( G4double ekin,
G4double kbz,
G4double rho,
G4double cs,
G4double Aac,
G4double Eac,
G4double prefactor )

Definition at line 422 of file G4MicroElecElasticModel_new.cc.

429{
430
431 G4double e = 1.6e-19,
432 m0 = 9.10938356e-31,
433 h = 1.0546e-34,
434 kb = 1.38e-23;
435
436 G4double E = (ekin / eV) * e;
437 G4double D = (2 / (std::sqrt(2) * std::pow(pi, 2) * std::pow(h, 3))) * (1 + 2 * E) * std::pow(m0, 1.5) * std::sqrt(E);
438
439 // Parametres SiO2
440 G4double T = 300,
441 Ebz = (std::pow(h, 2) * std::pow(kbz, 2)) / (2 * m0),
442 hwbz = cs * kbz * h,
443 nbz = 1.0 / (exp(hwbz / (kb * T)) - 1),
444 Pac;
445
446 if (E < Ebz / 4.0)
447 {
448 Pac = ((pi * kb * T) / (h * std::pow(cs, 2) * rho)) * (std::pow(Eac, 2) * D) / (1 + (E / Aac));
449 }
450
451 else if (E > Ebz) //Screened relationship
452 {
453 Pac = ((2 * pi * m0 * (2 * nbz + 1)) / (h * rho * hwbz)) * std::pow(Eac, 2) * D * E * 2 * std::pow((Aac / E), 2) * (((-E / Aac) / (1 + (E / Aac))) + log(1 + (E / Aac)));
454 }
455 else //Linear interpolation
456 {
457 G4double fEbz = ((2 * pi * m0 * (2 * nbz + 1)) / (h * rho * hwbz)) * std::pow(Eac, 2) * D * Ebz * 2 * std::pow((Aac / Ebz), 2) * (((-Ebz / Aac) / (1 + (Ebz / Aac))) + log(1 + (Ebz / Aac)));
458 G4double fEbz4 = ((pi * kb * T) / (h * std::pow(cs, 2) * rho)) * (std::pow(Eac, 2) * D) / (1 + ((Ebz / 4) / Aac));
459 G4double alpha = ((fEbz - fEbz4) / (Ebz - (Ebz / 4)));
460 Pac = alpha * E + (fEbz - alpha * Ebz);
461 }
462
463 G4double MFP = (std::sqrt(2 * E / m0) / (prefactor * Pac)) * m;
464
465 return 1 / MFP;
466}
G4double D(G4double temp)
double G4double
Definition G4Types.hh:83

Referenced by CrossSectionPerVolume().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 292 of file G4MicroElecElasticModel_new.cc.

297{
298 if (verboseLevel > 3)
299 G4cout << "Calling CrossSectionPerVolume() of G4MicroElecElasticModel" << G4endl;
300
301 isOkToBeInitialised = true;
302 currentMaterialName = material->GetName().substr(3, material->GetName().size());
303 const G4DataVector cuts;
304 Initialise(p, cuts);
305 // Calculate total cross section for model
306 MapEnergy::iterator lowEPos;
307 lowEPos = lowEnergyLimitTable.find(currentMaterialName);
308
309 MapEnergy::iterator highEPos;
310 highEPos = highEnergyLimitTable.find(currentMaterialName);
311
312 MapEnergy::iterator killEPos;
313 killEPos = workFunctionTable.find(currentMaterialName);
314
315 if (lowEPos == lowEnergyLimitTable.end() || highEPos == highEnergyLimitTable.end() || killEPos == workFunctionTable.end())
316 {
317 G4String str = "Material ";
318 str += currentMaterialName + " not found!";
319 G4Exception("G4MicroElecElasticModel_new::EnergyLimits", "em0002", FatalException, str);
320 return 0;
321 }
322 else {
323 // G4cout << "normal elastic " << G4endl;
324 lowEnergyLimit = lowEPos->second;
325 highEnergyLimit = highEPos->second;
326 killBelowEnergy = killEPos->second;
327
328 }
329
330 if (ekin < killBelowEnergy) {
331
332 return DBL_MAX; }
333
334 G4double sigma=0;
335
336 //Phonon for SiO2
337 if (currentMaterialName == "SILICON_DIOXIDE" && ekin < 100 * eV) {
338 acousticModelEnabled = true;
339
340 //Values for SiO2
341 G4double kbz = 11.54e9,
342 rho = 2.2 * 1000, // [g/cm3] * 1000
343 cs = 3560, //Sound speed
344 Ebz = 5.1 * 1.6e-19,
345 Aac = 17 * Ebz, //A screening parameter
346 Eac = 3.5 * 1.6e-19, //C deformation potential
347 prefactor = 2.2;// Facteur pour modifier les MFP
348
349 return AcousticCrossSectionPerVolume(ekin, kbz, rho, cs, Aac, Eac, prefactor);
350 }
351
352else if (currentMaterialName == "ALUMINUM_OXIDE" && ekin < 20 * eV) {
353 acousticModelEnabled = true;
354
355 //Values for Al2O3
356 G4double kbz = 8871930614.247564,
357 rho = 3.97 * 1000, // [g/cm3] * 1000
358 cs = 233329.07733059773, //Sound speed
359 Aac = 2.9912494342262614e-19, //A screening parameter
360 Eac = 2.1622471654789847e-18, //C deformation potential
361 prefactor = 1;
362 return AcousticCrossSectionPerVolume(ekin, kbz, rho, cs, Aac, Eac, prefactor);
363 }
364 //Elastic
365 else {
366 acousticModelEnabled = false;
367
368 G4double density = material->GetTotNbOfAtomsPerVolume();
369 const G4String& particleName = p->GetParticleName();
370
371 TCSMap::iterator tablepos;
372 tablepos = tableTCS.find(currentMaterialName);
373
374 if (tablepos != tableTCS.end())
375 {
376 MapData* tableData = tablepos->second;
377
378 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
379 {
380 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
381 pos = tableData->find(particleName);
382
383 if (pos != tableData->end())
384 {
386 if (table != 0)
387 {
388 sigma = table->FindValue(ekin);
389 }
390 }
391 else
392 {
393 G4Exception("G4MicroElecElasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, "Model not applicable to particle type.");
394 }
395 }
396 else return 1 / DBL_MAX;
397 }
398 else
399 {
400 G4String str = "Material ";
401 str += currentMaterialName + " TCS table not found!";
402 G4Exception("G4MicroElecElasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
403 }
404
405 if (verboseLevel > 3)
406 {
407 G4cout << "---> Kinetic energy(eV)=" << ekin / eV << G4endl;
408 G4cout << " - Cross section per Si atom (cm^2)=" << sigma / cm / cm << G4endl;
409 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) << G4endl;
410 }
411
412 // Hsing-YinChangaAndrewAlvaradoaTreyWeberaJaimeMarianab Monte Carlo modeling of low-energy electron-induced secondary electron emission yields in micro-architected boron nitride surfaces - ScienceDirect, (n.d.). https://www.sciencedirect.com/science/article/pii/S0168583X19304069 (accessed April 1, 2022).
413 if (currentMaterialName == "BORON_NITRIDE") {
414 sigma = sigma * tanh(0.5 * pow(ekin / 5.2e-6, 2));
415 }
416 return sigma*density;
417 }
418}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double GetTotNbOfAtomsPerVolume() const
const G4String & GetName() const
G4double FindValue(G4double e, G4int componentId=0) const override
G4double AcousticCrossSectionPerVolume(G4double ekin, G4double kbz, G4double rho, G4double cs, G4double Aac, G4double Eac, G4double prefactor)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4String & GetParticleName() const
#define DBL_MAX
Definition templates.hh:62

◆ DamageEnergy()

G4double G4MicroElecElasticModel_new::DamageEnergy ( G4double T,
G4double A,
G4double Z )

Definition at line 523 of file G4MicroElecElasticModel_new.cc.

524{
525 //.................. T in eV!!!!!!!!!!!!!
526 G4double Z2= Z;
527 G4double M2= A;
528 G4double k_d;
529 G4double epsilon_d;
530 G4double g_epsilon_d;
531 G4double E_nu;
532
533 k_d=0.1334*std::pow(Z2,(2./3.))*std::pow(M2,(-1./2.));
534 epsilon_d=0.01014*std::pow(Z2,(-7./3.))*(T/eV);
535 g_epsilon_d= epsilon_d+0.40244*std::pow(epsilon_d,(3./4.))+3.4008*std::pow(epsilon_d,(1./6.));
536
537 E_nu=1./(1.+ k_d*g_epsilon_d);
538
539 return E_nu;
540}
const G4double A[17]

◆ GetKillBelowThreshold()

G4double G4MicroElecElasticModel_new::GetKillBelowThreshold ( )
inline

Definition at line 118 of file G4MicroElecElasticModel_new.hh.

118{ return killBelowEnergy; }

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 168 of file G4MicroElecElasticModel_new.cc.

170{
171 if (isOkToBeInitialised == true && isInitialised == false) {
172
173 if (verboseLevel > -1)
174 G4cout << "Calling G4MicroElecElasticModel_new::Initialise()" << G4endl;
175 // Energy limits
176 // Reading of data files
177
178 G4double scaleFactor = 1e-18 * cm * cm;
179
180 G4ProductionCutsTable* theCoupleTable =
182 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
183
184 for (G4int i = 0; i < numOfCouples; ++i) {
185 const G4Material* material =
186 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
187
188 //theCoupleTable->GetMaterialCutsCouple(i)->;
189
190 G4cout << "MicroElasticModel, Material " << i + 1 << " / " << numOfCouples << " : " << material->GetName() << G4endl;
191 if (material->GetName() == "Vacuum") continue;
192
193 G4String matName = material->GetName().substr(3, material->GetName().size());
194 G4cout<< matName<< G4endl;
195
196 currentMaterialStructure = new G4MicroElecMaterialStructure(matName);
197 lowEnergyLimitTable[matName]=currentMaterialStructure->GetElasticModelLowLimit();
198 highEnergyLimitTable[matName]=currentMaterialStructure->GetElasticModelHighLimit();
199 workFunctionTable[matName] = currentMaterialStructure->GetWorkFunction();
200
201 delete currentMaterialStructure;
202
203 G4cout << "Reading TCS file" << G4endl;
204 G4String fileElectron = "Elastic/elsepa_elastic_cross_e_" + matName;
205 G4cout << "Elastic Total Cross file : " << fileElectron << G4endl;
206
208 G4String electron = electronDef->GetParticleName();
209
210 // For total cross section
211 MapData* tableData = new MapData();
212
214 tableE->LoadData(fileElectron);
215 tableData->insert(make_pair(electron, tableE));
216 tableTCS[matName] = tableData; //Storage of TCS
217
218 // For final state
219 const char* path = G4FindDataDir("G4LEDATA");
220 if (!path)
221 {
222 G4Exception("G4MicroElecElasticModel_new::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
223 return;
224 }
225
226 //Reading DCS file
227 std::ostringstream eFullFileName;
228 eFullFileName << path << "/microelec/Elastic/elsepa_elastic_cumulated_diffcross_e_" + matName + ".dat";
229 G4cout << "Elastic Cumulated Diff Cross : " << eFullFileName.str().c_str() << G4endl;
230 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
231
232 if (!eDiffCrossSection)
233 G4Exception("G4MicroElecElasticModel_new::Initialise", "em0003", FatalException, "Missing data file: /microelec/sigmadiff_cumulated_elastic_e_Si.dat");
234
235 // October 21th, 2014 - Melanie Raine
236 // Added clear for MT
237 // Diff Cross Sections in cumulated mode
238 TriDimensionMap* eDiffCrossSectionData = new TriDimensionMap(); //Angles
239 std::vector<G4double>* eTdummyVec = new std::vector<G4double>; //Incident energy vector
240 VecMap* eProbVec = new VecMap; //Probabilities
241
242 eTdummyVec->push_back(0.);
243
244 while (!eDiffCrossSection.eof())
245 {
246 G4double tDummy; //incident energy
247 G4double eProb; //Proba
248 eDiffCrossSection >> tDummy >> eProb;
249
250 // SI : mandatory eVecm initialization
251 if (tDummy != eTdummyVec->back())
252 {
253 eTdummyVec->push_back(tDummy); //adding values for incident energy points
254 (*eProbVec)[tDummy].push_back(0.); //adding probability for the first angle, equal to 0
255 }
256
257 eDiffCrossSection >> (*eDiffCrossSectionData)[tDummy][eProb]; //adding Angle Value to map
258
259 if (eProb != (*eProbVec)[tDummy].back()) {
260 (*eProbVec)[tDummy].push_back(eProb); //Adding cumulated proba to map
261 }
262
263 }
264
265 //Filling maps for the material
266 thetaDataStorage[matName] = eDiffCrossSectionData;
267 eIncidentEnergyStorage[matName] = eTdummyVec;
268 eProbaStorage[matName] = eProbVec;
269 }
270 // End final state
271
272 if (verboseLevel > 2)
273 G4cout << "Loaded cross section files for MicroElec Elastic model" << G4endl;
274
275 if (verboseLevel > 0)
276 {
277 G4cout << "MicroElec Elastic model is initialized " << G4endl
278 << "Energy range: "
279 << LowEnergyLimit() / eV << " eV - "
280 << HighEnergyLimit() / MeV << " MeV"
281 << G4endl; // system("pause"); linux doesn't like
282 }
283
284 if (isInitialised) { return; }
286 isInitialised = true;
287 }
288}
const char * G4FindDataDir(const char *)
int G4int
Definition G4Types.hh:85
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
const G4Material * GetMaterial() const
G4bool LoadData(const G4String &argFileName) override
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const

Referenced by CrossSectionPerVolume().

◆ SampleSecondaries()

void G4MicroElecElasticModel_new::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * aDynamicElectron,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 471 of file G4MicroElecElasticModel_new.cc.

476{
477
478 if (verboseLevel > 3)
479 G4cout << "Calling SampleSecondaries() of G4MicroElecElasticModel" << G4endl;
480
481 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
482
483 if (electronEnergy0 < killBelowEnergy)
484 {
488 return;
489 }
490
491 if (electronEnergy0 < highEnergyLimit)
492 {
493 G4double cosTheta = 0;
494 if (acousticModelEnabled)
495 {
496 cosTheta = 1 - 2 * G4UniformRand(); //Isotrope
497 }
498 else if (electronEnergy0 >= lowEnergyLimit)
499 {
500 cosTheta = RandomizeCosTheta(electronEnergy0);
501 }
502
503 G4double phi = 2. * pi * G4UniformRand();
504
505 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
506 G4ThreeVector xVers = zVers.orthogonal();
507 G4ThreeVector yVers = zVers.cross(xVers);
508
509 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
510 G4double yDir = xDir;
511 xDir *= std::cos(phi);
512 yDir *= std::sin(phi);
513
514 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
515
518 }
519}
@ fStopAndKill
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetKillBelowThreshold()

void G4MicroElecElasticModel_new::SetKillBelowThreshold ( G4double threshold)

Definition at line 699 of file G4MicroElecElasticModel_new.cc.

700{
701 killBelowEnergy = threshold;
702
703 if (threshold < 5*CLHEP::eV)
704 {
705 G4Exception ("*** WARNING : the G4MicroElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
706 threshold = 5*CLHEP::eV;
707 }
708}
@ JustWarning

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MicroElecElasticModel_new::fParticleChangeForGamma
protected

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