Geant4 10.7.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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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

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

◆ ~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 {
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<double>* 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 404 of file G4MicroElecElasticModel_new.cc.

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

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

506{
507 //.................. T in eV!!!!!!!!!!!!!
508 G4double Z2= Z;
509 G4double M2= A;
510 G4double k_d;
511 G4double epsilon_d;
512 G4double g_epsilon_d;
513 G4double E_nu;
514
515 k_d=0.1334*std::pow(Z2,(2./3.))*std::pow(M2,(-1./2.));
516 epsilon_d=0.01014*std::pow(Z2,(-7./3.))*(T/eV);
517 g_epsilon_d= epsilon_d+0.40244*std::pow(epsilon_d,(3./4.))+3.4008*std::pow(epsilon_d,(1./6.));
518
519 E_nu=1./(1.+ k_d*g_epsilon_d);
520
521 return E_nu;
522}
double A(double temperature)

◆ GetKillBelowThreshold()

G4double G4MicroElecElasticModel_new::GetKillBelowThreshold ( )
inline

Definition at line 119 of file G4MicroElecElasticModel_new.hh.

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

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

458{
459
460 if (verboseLevel > 3)
461 G4cout << "Calling SampleSecondaries() of G4MicroElecElasticModel" << G4endl;
462
463 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
464
465 if (electronEnergy0 < killBelowEnergy)
466 {
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
485 G4double phi = 2. * pi * G4UniformRand();
486
487 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
488 G4ThreeVector xVers = zVers.orthogonal();
489 G4ThreeVector yVers = zVers.cross(xVers);
490
491 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
492 G4double yDir = xDir;
493 xDir *= std::cos(phi);
494 yDir *= std::sin(phi);
495
496 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
497
500 }
501}
@ 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(G4double Px, G4double Py, G4double Pz)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetKillBelowThreshold()

void G4MicroElecElasticModel_new::SetKillBelowThreshold ( G4double  threshold)

Definition at line 681 of file G4MicroElecElasticModel_new.cc.

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

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MicroElecElasticModel_new::fParticleChangeForGamma
protected

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