Geant4 11.4.0
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 G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, 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)
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)
void SetMasterThread (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}
#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 125 of file G4MicroElecElasticModel_new.cc.

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

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 405 of file G4MicroElecElasticModel_new.cc.

412{
413
414 G4double e = 1.6e-19,
415 m0 = 9.10938356e-31,
416 h = 1.0546e-34,
417 kb = 1.38e-23;
418
419 G4double E = (ekin / eV) * e;
420 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);
421
422 // Parametres SiO2
423 G4double T = 300,
424 Ebz = (std::pow(h, 2) * std::pow(kbz, 2)) / (2 * m0),
425 hwbz = cs * kbz * h,
426 nbz = 1.0 / (exp(hwbz / (kb * T)) - 1),
427 Pac;
428
429 if (E < Ebz / 4.0)
430 {
431 Pac = ((pi * kb * T) / (h * std::pow(cs, 2) * rho)) * (std::pow(Eac, 2) * D) / (1 + (E / Aac));
432 }
433
434 else if (E > Ebz) //Screened relationship
435 {
436 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)));
437 }
438 else //Linear interpolation
439 {
440 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)));
441 G4double fEbz4 = ((pi * kb * T) / (h * std::pow(cs, 2) * rho)) * (std::pow(Eac, 2) * D) / (1 + ((Ebz / 4) / Aac));
442 G4double alpha = ((fEbz - fEbz4) / (Ebz - (Ebz / 4)));
443 Pac = alpha * E + (fEbz - alpha * Ebz);
444 }
445
446 G4double MFP = (std::sqrt(2 * E / m0) / (prefactor * Pac)) * m;
447
448 return 1 / MFP;
449}
G4double D(G4double temp)
double G4double
Definition G4Types.hh:83
const G4double pi

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 283 of file G4MicroElecElasticModel_new.cc.

288{
289 if (verboseLevel > 3)
290 G4cout << "Calling CrossSectionPerVolume() of G4MicroElecElasticModel" << G4endl;
291
292 currentMaterialName = material->GetName().substr(3, material->GetName().size());
293 const G4DataVector cuts;
294 // Calculate total cross section for model
295 MapEnergy::iterator lowEPos;
296 lowEPos = lowEnergyLimitTable.find(currentMaterialName);
297
298 MapEnergy::iterator highEPos;
299 highEPos = highEnergyLimitTable.find(currentMaterialName);
300
301 MapEnergy::iterator killEPos;
302 killEPos = workFunctionTable.find(currentMaterialName);
303
304 if (lowEPos == lowEnergyLimitTable.end() || highEPos == highEnergyLimitTable.end() || killEPos == workFunctionTable.end())
305 {
306 G4String str = "Material ";
307 str += currentMaterialName + " not found!";
308 G4Exception("G4MicroElecElasticModel_new::EnergyLimits", "em0002", FatalException, str);
309 return 0;
310 }
311 else {
312 // G4cout << "normal elastic " << G4endl;
313 lowEnergyLimit = lowEPos->second;
314 highEnergyLimit = highEPos->second;
315 killBelowEnergy = killEPos->second;
316 }
317
318 if (ekin < killBelowEnergy) { return DBL_MAX; }
319
320 G4double sigma=0;
321
322 //Phonon for SiO2
323 if (currentMaterialName == "SILICON_DIOXIDE" && ekin < 100 * eV) {
324 acousticModelEnabled = true;
325
326 //Values for SiO2
327 G4double kbz = 11.54e9,
328 rho = 2.2 * 1000, // [g/cm3] * 1000
329 cs = 3560, //Sound speed
330 Ebz = 5.1 * 1.6e-19,
331 Aac = 17 * Ebz, //A screening parameter
332 Eac = 3.5 * 1.6e-19, //C deformation potential
333 prefactor = 2.2;// Facteur pour modifier les MFP
334
335 return AcousticCrossSectionPerVolume(ekin, kbz, rho, cs, Aac, Eac, prefactor);
336 }
337
338 else if (currentMaterialName == "ALUMINUM_OXIDE" && ekin < 20 * eV) {
339 acousticModelEnabled = true;
340
341 //Values for Al2O3
342 G4double kbz = 8871930614.247564,
343 rho = 3.97 * 1000, // [g/cm3] * 1000
344 cs = 233329.07733059773, //Sound speed
345 Aac = 2.9912494342262614e-19, //A screening parameter
346 Eac = 2.1622471654789847e-18, //C deformation potential
347 prefactor = 1;
348 return AcousticCrossSectionPerVolume(ekin, kbz, rho, cs, Aac, Eac, prefactor);
349 }
350 //Elastic
351 else {
352 acousticModelEnabled = false;
353
354 G4double density = material->GetTotNbOfAtomsPerVolume();
355 const G4String& particleName = p->GetParticleName();
356
357 TCSMap::iterator tablepos;
358 tablepos = tableTCS.find(currentMaterialName);
359
360 if (tablepos != tableTCS.end())
361 {
362 MapData* tableData = tablepos->second;
363
364 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
365 {
366 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
367 pos = tableData->find(particleName);
368
369 if (pos != tableData->end()){
370 G4MicroElecCrossSectionDataSet_new* table = pos->second;
371 if (table != 0)
372 {
373 sigma = table->FindValue(ekin);
374 }
375 }
376 else
377 {
378 G4Exception("G4MicroElecElasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, "Model not applicable to particle type.");
379 }
380 }
381 else return 1 / DBL_MAX;
382 }
383 else {
384 G4String str = "Material ";
385 str += currentMaterialName + " TCS table not found!";
386 G4Exception("G4MicroElecElasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
387 }
388
389 if (verboseLevel > 3) {
390 G4cout << "---> Kinetic energy(eV)=" << ekin / eV << G4endl;
391 G4cout << " - Cross section per Si atom (cm^2)=" << sigma / cm / cm << G4endl;
392 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) << G4endl;
393 }
394
395 // 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).
396 if (currentMaterialName == "BORON_NITRIDE") {
397 sigma = sigma * tanh(0.5 * pow(ekin / 5.2e-6, 2));
398 }
399 return sigma*density;
400 }
401}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double GetTotNbOfAtomsPerVolume() const
const G4String & GetName() const
G4double AcousticCrossSectionPerVolume(G4double ekin, G4double kbz, G4double rho, G4double cs, G4double Aac, G4double Eac, G4double prefactor)
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 504 of file G4MicroElecElasticModel_new.cc.

505{
506 //.................. T in eV!!!!!!!!!!!!!
507 G4double Z2= Z;
508 G4double M2= A;
509 G4double k_d;
510 G4double epsilon_d;
511 G4double g_epsilon_d;
512 G4double E_nu;
513
514 k_d=0.1334*std::pow(Z2,(2./3.))*std::pow(M2,(-1./2.));
515 epsilon_d=0.01014*std::pow(Z2,(-7./3.))*(T/eV);
516 g_epsilon_d= epsilon_d+0.40244*std::pow(epsilon_d,(3./4.))+3.4008*std::pow(epsilon_d,(1./6.));
517
518 E_nu=1./(1.+ k_d*g_epsilon_d);
519
520 return E_nu;
521}
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 167 of file G4MicroElecElasticModel_new.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 454 of file G4MicroElecElasticModel_new.cc.

459{
460 if (verboseLevel > 3)
461 G4cout << "Calling SampleSecondaries() of G4MicroElecElasticModel" << G4endl;
462
463 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
464
465 if (electronEnergy0 < killBelowEnergy)
466 {
467 fParticleChangeForGamma->SetProposedKineticEnergy(0.);
468 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
469 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
470 return;
471 }
472
473 if (electronEnergy0 < highEnergyLimit)
474 {
475 G4double cosTheta = 0;
476 if (acousticModelEnabled)
477 {
478 cosTheta = 1 - 2 * G4UniformRand(); //Isotrope
479 }
480 else if (electronEnergy0 >= lowEnergyLimit)
481 {
482 cosTheta = RandomizeCosTheta(electronEnergy0);
483 }
484 G4double phi = 2. * pi * G4UniformRand();
485
486 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
487 G4ThreeVector xVers = zVers.orthogonal();
488 G4ThreeVector yVers = zVers.cross(xVers);
489
490 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
491 G4double yDir = xDir;
492 xDir *= std::cos(phi);
493 yDir *= std::sin(phi);
494
495 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
496
497 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit());
498 fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
499 }
500}
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const

◆ SetKillBelowThreshold()

void G4MicroElecElasticModel_new::SetKillBelowThreshold ( G4double threshold)

Definition at line 679 of file G4MicroElecElasticModel_new.cc.

680{
681 killBelowEnergy = threshold;
682
683 if (threshold < 5*CLHEP::eV)
684 {
685 G4Exception ("*** WARNING : the G4MicroElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
686 threshold = 5*CLHEP::eV;
687 }
688}
@ JustWarning

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MicroElecElasticModel_new::fParticleChangeForGamma
protected

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