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

#include <G4MicroElecInelasticModel_new.hh>

+ Inheritance diagram for G4MicroElecInelasticModel_new:

Public Member Functions

 G4MicroElecInelasticModel_new (const G4ParticleDefinition *p=nullptr, const G4String &nam="MicroElecInelasticModel")
 
 ~G4MicroElecInelasticModel_new () override
 
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 (const G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
G4double ComputeRelativistVelocity (G4double E, G4double mass)
 
G4double ComputeElasticQmax (G4double T1i, G4double T2i, G4double m1, G4double m2)
 
G4double BKZ (G4double Ep, G4double mp, G4int Zp, G4double EF)
 
G4double stepFunc (G4double x)
 
G4double vrkreussler (G4double v, G4double vF)
 
G4MicroElecInelasticModel_newoperator= (const G4MicroElecInelasticModel_new &right)=delete
 
 G4MicroElecInelasticModel_new (const G4MicroElecInelasticModel_new &)=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
 

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 *)
 
- 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
 

Detailed Description

Definition at line 93 of file G4MicroElecInelasticModel_new.hh.

Constructor & Destructor Documentation

◆ G4MicroElecInelasticModel_new() [1/2]

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

Definition at line 97 of file G4MicroElecInelasticModel_new.cc.

99 :G4VEmModel(nam),isInitialised(false)
100{
101
102 verboseLevel= 0;
103 // Verbosity scale:
104 // 0 = nothing
105 // 1 = warning for energy non-conservation
106 // 2 = details of energy budget
107 // 3 = calculation of cross sections, file openings, sampling of atoms
108 // 4 = entering in methods
109
110 if( verboseLevel>0 )
111 {
112 G4cout << "MicroElec inelastic model is constructed " << G4endl;
113 }
114
115 //Mark this model as "applicable" for atomic deexcitation
117 fAtomDeexcitation = nullptr;
118 fParticleChangeForGamma = nullptr;
119
120 // default generator
122
123 // Selection of computation method
124 fasterCode = true;
125 SEFromFermiLevel = false;
126}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

◆ ~G4MicroElecInelasticModel_new()

G4MicroElecInelasticModel_new::~G4MicroElecInelasticModel_new ( )
override

Definition at line 130 of file G4MicroElecInelasticModel_new.cc.

131{
132 // Cross section
133 // (0)
134 for (auto const& p : tableTCS) {
135 MapData* tableData = p.second;
136 for (auto const& pos : *tableData) { delete pos.second; }
137 delete tableData;
138 }
139 tableTCS.clear();
140
141 // (1)
142 for (auto const & obj : eDiffDatatable) {
143 auto ptr = obj.second;
144 ptr->clear();
145 delete ptr;
146 }
147
148 for (auto const & obj : pDiffDatatable) {
149 auto ptr = obj.second;
150 ptr->clear();
151 delete ptr;
152 }
153
154 // (2)
155 for (auto const & obj : eNrjTransStorage) {
156 delete obj.second;
157 }
158 for (auto const & obj : pNrjTransStorage) {
159 delete obj.second;
160 }
161
162 // (3)
163 for (auto const& p : eProbaShellStorage) {
164 delete p.second;
165 }
166
167 for (auto const& p : pProbaShellStorage) {
168 delete p.second;
169 }
170
171 // (4)
172 for (auto const& p : eIncidentEnergyStorage) {
173 delete p.second;
174 }
175
176 for (auto const& p : pIncidentEnergyStorage) {
177 delete p.second;
178 }
179
180 // (5)
181 for (auto const& p : eVecmStorage) {
182 delete p.second;
183 }
184
185 for (auto const& p : pVecmStorage) {
186 delete p.second;
187 }
188
189 // (6)
190 for (auto const& p : tableMaterialsStructures) {
191 delete p.second;
192 }
193}

◆ G4MicroElecInelasticModel_new() [2/2]

G4MicroElecInelasticModel_new::G4MicroElecInelasticModel_new ( const G4MicroElecInelasticModel_new & )
delete

Member Function Documentation

◆ BKZ()

G4double G4MicroElecInelasticModel_new::BKZ ( G4double Ep,
G4double mp,
G4int Zp,
G4double EF )

Definition at line 1410 of file G4MicroElecInelasticModel_new.cc.

1411{
1412 // need atomic unit conversion
1413 G4double hbar = hbar_Planck, hbar2 = hbar*hbar, me = electron_mass_c2/c_squared, Ry = me*elm_coupling*elm_coupling/(2*hbar2);
1414 G4double hartree = 2*Ry, a0 = Bohr_radius, velocity = a0*hartree/hbar;
1416
1417 vp /= velocity;
1418
1419 G4double wp = Eplasmon/hartree;
1420 G4double a = std::pow(4./9./CLHEP::pi, 1./3.);
1421 G4double vF = std::pow(wp*wp/(3.*a*a*a), 1./3.);
1422 G4double c = 0.9;
1423 G4double vr = vrkreussler(vp /*in u.a*/, vF /*in u.a*/);
1424 G4double yr = vr/std::pow(Zp, 2./3.);
1425 G4double q = 0.;
1426 if(Zp==2) q = 1-exp(-c*vr/(Zp-5./16.));
1427 else q = 1.-exp(-c*(yr-0.07));
1428 G4double Neq = Zp*(1.-q);
1429 G4double l0 = 0.;
1430 if(Neq<=2) l0 = 3./(Zp-0.3*(Neq-1))/2.;
1431 else l0 = 0.48*std::pow(Neq, 2./3.)/(Zp-Neq/7.);
1432 if(Zp==2) c = 1.0;
1433 else c = 3./2.;
1434 return Zp*(q + c*(1.-q)/vF/vF/2.0 * log(1.+std::pow(2.*l0*vF,2.)));
1435}
const G4double a0
double G4double
Definition G4Types.hh:83
G4double ComputeRelativistVelocity(G4double E, G4double mass)
G4double vrkreussler(G4double v, G4double vF)

Referenced by CrossSectionPerVolume().

◆ ComputeElasticQmax()

G4double G4MicroElecInelasticModel_new::ComputeElasticQmax ( G4double T1i,
G4double T2i,
G4double m1,
G4double m2 )

Definition at line 1377 of file G4MicroElecInelasticModel_new.cc.

1377 {
1378 G4double v1i = ComputeRelativistVelocity(T1i, M1);
1379 G4double v2i = ComputeRelativistVelocity(T2i, M2);
1380
1381 G4double v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*-1*v2i;
1382 G4double vtransfer2a = v2f*v2f-v2i*v2i;
1383
1384 v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*v2i;
1385 G4double vtransfer2b = v2f*v2f-v2i*v2i;
1386
1387 G4double vtransfer2 = std::max(vtransfer2a, vtransfer2b);
1388 return 0.5*M2*vtransfer2;
1389}

◆ ComputeRelativistVelocity()

G4double G4MicroElecInelasticModel_new::ComputeRelativistVelocity ( G4double E,
G4double mass )

Definition at line 1370 of file G4MicroElecInelasticModel_new.cc.

1370 {
1371 G4double x = E/mass;
1372 return c_light*std::sqrt(x*(x + 2.0))/(x + 1.0);
1373}

Referenced by BKZ(), and ComputeElasticQmax().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 489 of file G4MicroElecInelasticModel_new.cc.

494{
495 if (verboseLevel > 3) G4cout << "Calling CrossSectionPerVolume() of G4MicroElecInelasticModel" << G4endl;
496
497 G4double density = material->GetTotNbOfAtomsPerVolume();
498 currentMaterial = material->GetName().substr(3, material->GetName().size());
499
500 MapStructure::iterator structPos;
501 structPos = tableMaterialsStructures.find(currentMaterial);
502
503 // Calculate total cross section for model
504 TCSMap::iterator tablepos;
505 tablepos = tableTCS.find(currentMaterial);
506
507 if (tablepos == tableTCS.end() )
508 {
509 G4String str = "Material ";
510 str += currentMaterial + " TCS Table not found!";
511 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
512 return 0;
513 }
514 else if(structPos == tableMaterialsStructures.end())
515 {
516 G4String str = "Material ";
517 str += currentMaterial + " Structure not found!";
518 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
519 return 0;
520 }
521 else {
522 MapData* tableData = tablepos->second;
523 currentMaterialStructure = structPos->second;
524
525 G4double sigma = 0;
526
527 const G4String& particleName = particleDefinition->GetParticleName();
528 G4String nameLocal = particleName;
529 G4int pdg = particleDefinition->GetPDGEncoding();
530 G4int Z = particleDefinition->GetAtomicNumber();
531
532 G4double Zeff = 1.0, Zeff2 = Zeff*Zeff;
533 G4double Mion_c2 = particleDefinition->GetPDGMass();
534
535 if (Mion_c2 > proton_mass_c2)
536 {
537 ekin *= proton_mass_c2 / Mion_c2;
538 nameLocal = "proton";
539 }
540
541 G4double lowLim = currentMaterialStructure->GetInelasticModelLowLimit(pdg);
542 G4double highLim = currentMaterialStructure->GetInelasticModelHighLimit(pdg);
543
544 if (ekin >= lowLim && ekin < highLim)
545 {
546 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
547 pos = tableData->find(nameLocal); //find particle type
548
549 if (pos != tableData->end())
550 {
552 if (table != 0)
553 {
554 sigma = table->FindValue(ekin);
555
556 if (Mion_c2 > proton_mass_c2) {
557 sigma = 0.;
558 for (G4int i = 0; i < currentMaterialStructure->NumberOfLevels(); i++) {
559 Zeff = BKZ(ekin / (proton_mass_c2 / Mion_c2), Mion_c2 / c_squared, Z, currentMaterialStructure->Energy(i)); // il faut garder le vrai ekin car le calcul à l'interieur de la methode convertie l'énergie en vitesse
560 Zeff2 = Zeff*Zeff;
561 sigma += Zeff2*table->FindShellValue(ekin, i);
562// il faut utiliser le ekin mis à l'echelle pour chercher la bonne
563// valeur dans les tables proton
564
565 }
566 }
567 else {
568 sigma = table->FindValue(ekin);
569 }
570 }
571 }
572 else
573 {
574 G4Exception("G4MicroElecInelasticModel_new::CrossSectionPerVolume",
575 "em0002", FatalException,
576 "Model not applicable to particle type.");
577 }
578 }
579 else
580 {
581 return 1 / DBL_MAX;
582 }
583
584 if (verboseLevel > 3)
585 {
586 G4cout << "---> Kinetic energy (eV)=" << ekin / eV << G4endl;
587 G4cout << " - Cross section per Si atom (cm^2)=" << sigma / cm2 << G4endl;
588 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) << G4endl;
589 }
590
591 return (sigma)*density;
592 }
593}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
int G4int
Definition G4Types.hh:85
G4double GetTotNbOfAtomsPerVolume() const
const G4String & GetName() const
G4double FindShellValue(G4double argEnergy, G4int shell) const
G4double FindValue(G4double e, G4int componentId=0) const override
G4double BKZ(G4double Ep, G4double mp, G4int Zp, G4double EF)
#define DBL_MAX
Definition templates.hh:62

◆ DifferentialCrossSection()

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

Definition at line 1109 of file G4MicroElecInelasticModel_new.cc.

1114{
1115 G4double sigma = 0.;
1116
1117 if (energyTransfer >= currentMaterialStructure->GetLimitEnergy(LevelIndex))
1118 {
1119 G4double valueT1 = 0;
1120 G4double valueT2 = 0;
1121 G4double valueE21 = 0;
1122 G4double valueE22 = 0;
1123 G4double valueE12 = 0;
1124 G4double valueE11 = 0;
1125
1126 G4double xs11 = 0;
1127 G4double xs12 = 0;
1128 G4double xs21 = 0;
1129 G4double xs22 = 0;
1130
1131 if (particleDefinition == G4Electron::ElectronDefinition())
1132 {
1133
1134 dataDiffCSMap::iterator iterator_Proba;
1135 iterator_Proba = eDiffDatatable.find(currentMaterial);
1136
1137 incidentEnergyMap::iterator iterator_Nrj;
1138 iterator_Nrj = eIncidentEnergyStorage.find(currentMaterial);
1139
1140 TranfEnergyMap::iterator iterator_TransfNrj;
1141 iterator_TransfNrj = eVecmStorage.find(currentMaterial);
1142
1143 if (iterator_Proba != eDiffDatatable.end() && iterator_Nrj != eIncidentEnergyStorage.end()
1144 && iterator_TransfNrj!= eVecmStorage.end())
1145 {
1146 vector<TriDimensionMap>* eDiffCrossSectionData = (iterator_Proba->second);
1147 vector<G4double>* eTdummyVec = iterator_Nrj->second; //Incident energies for interpolation
1148 VecMap* eVecm = iterator_TransfNrj->second;
1149
1150 // k should be in eV and energy transfer eV also
1151 auto t2 = std::upper_bound(eTdummyVec->begin(), eTdummyVec->end(), k);
1152 auto t1 = t2 - 1;
1153 // SI : the following condition avoids situations where energyTransfer >last vector element
1154 if (energyTransfer <= (*eVecm)[(*t1)].back() && energyTransfer <= (*eVecm)[(*t2)].back())
1155 {
1156 auto e12 = std::upper_bound((*eVecm)[(*t1)].begin(), (*eVecm)[(*t1)].end(), energyTransfer);
1157 auto e11 = e12 - 1;
1158 auto e22 = std::upper_bound((*eVecm)[(*t2)].begin(), (*eVecm)[(*t2)].end(), energyTransfer);
1159 auto e21 = e22 - 1;
1160
1161 valueT1 = *t1;
1162 valueT2 = *t2;
1163 valueE21 = *e21;
1164 valueE22 = *e22;
1165 valueE12 = *e12;
1166 valueE11 = *e11;
1167
1168 xs11 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1169 xs12 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1170 xs21 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1171 xs22 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1172 }
1173 }
1174 else {
1175 G4String str = "Material ";
1176 str += currentMaterial + " not found!";
1177 G4Exception("G4MicroElecDielectricModels::DifferentialCrossSection", "em0002", FatalException, str);
1178 }
1179 }
1180
1181 if (particleDefinition == G4Proton::ProtonDefinition())
1182 {
1183 dataDiffCSMap::iterator iterator_Proba;
1184 iterator_Proba = pDiffDatatable.find(currentMaterial);
1185
1186 incidentEnergyMap::iterator iterator_Nrj;
1187 iterator_Nrj = pIncidentEnergyStorage.find(currentMaterial);
1188
1189 TranfEnergyMap::iterator iterator_TransfNrj;
1190 iterator_TransfNrj = pVecmStorage.find(currentMaterial);
1191
1192 if (iterator_Proba != pDiffDatatable.end() && iterator_Nrj != pIncidentEnergyStorage.end()
1193 && iterator_TransfNrj != pVecmStorage.end())
1194 {
1195 vector<TriDimensionMap>* pDiffCrossSectionData = (iterator_Proba->second);
1196 vector<G4double>* pTdummyVec = iterator_Nrj->second; //Incident energies for interpolation
1197 VecMap* pVecm = iterator_TransfNrj->second;
1198
1199 // k should be in eV and energy transfer eV also
1200 auto t2 =
1201 std::upper_bound(pTdummyVec->begin(), pTdummyVec->end(), k);
1202 auto t1 = t2 - 1;
1203 if (energyTransfer <= (*pVecm)[(*t1)].back() && energyTransfer <= (*pVecm)[(*t2)].back())
1204 {
1205 auto e12 = std::upper_bound((*pVecm)[(*t1)].begin(), (*pVecm)[(*t1)].end(), energyTransfer);
1206 auto e11 = e12 - 1;
1207 auto e22 = std::upper_bound((*pVecm)[(*t2)].begin(), (*pVecm)[(*t2)].end(), energyTransfer);
1208 auto e21 = e22 - 1;
1209
1210 valueT1 = *t1;
1211 valueT2 = *t2;
1212 valueE21 = *e21;
1213 valueE22 = *e22;
1214 valueE12 = *e12;
1215 valueE11 = *e11;
1216
1217 xs11 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1218 xs12 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1219 xs21 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1220 xs22 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1221 }
1222 }
1223 else {
1224 G4String str = "Material ";
1225 str += currentMaterial + " not found!";
1226 G4Exception("G4MicroElecDielectricModels::DifferentialCrossSection", "em0002", FatalException, str);
1227 }
1228 }
1229
1230 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
1231 if (xsProduct != 0.)
1232 {
1233 sigma = QuadInterpolator( valueE11, valueE12,
1234 valueE21, valueE22,
1235 xs11, xs12,
1236 xs21, xs22,
1237 valueT1, valueT2,
1238 k, energyTransfer);
1239 }
1240
1241 }
1242
1243 return sigma;
1244}
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 197 of file G4MicroElecInelasticModel_new.cc.

199{
200 if (isInitialised) { return; }
201
202 if (verboseLevel > 3)
203 G4cout << "Calling G4MicroElecInelasticModel_new::Initialise()" << G4endl;
204
205 const char* path = G4FindDataDir("G4LEDATA");
206 if (!path)
207 {
208 G4Exception("G4MicroElecElasticModel_new::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
209 return;
210 }
211
212 G4String modelName = "mermin";
213 G4cout << "****************************" << G4endl;
214 G4cout << modelName << " model loaded !" << G4endl;
215
216 // Energy limits
219 G4String electron = electronDef->GetParticleName();
220 G4String proton = protonDef->GetParticleName();
221
222 G4double scaleFactor = 1.0;
223
224 // *** ELECTRON
225 lowEnergyLimit[electron] = 2 * eV;
226 highEnergyLimit[electron] = 10.0 * MeV;
227
228 // Cross section
230 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
231
232 for (G4int i = 0; i < numOfCouples; ++i) {
233 const G4Material* material = theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
234 G4cout << "Material " << i + 1 << " / " << numOfCouples << " : " << material->GetName() << G4endl;
235 if (material->GetName() == "Vacuum") continue;
236 G4String mat = material->GetName().substr(3, material->GetName().size());
237 MapData* tableData = new MapData;
238 currentMaterialStructure = new G4MicroElecMaterialStructure(mat);
239
240 tableMaterialsStructures[mat] = currentMaterialStructure;
241 if (particle == electronDef) {
242 //TCS
243 G4String fileElectron("Inelastic/" + modelName + "_sigma_inelastic_e-_" + mat);
244 G4cout << fileElectron << G4endl;
246 tableE->LoadData(fileElectron);
247 tableData->insert(make_pair(electron, tableE));
248
249 // DCS
250 std::ostringstream eFullFileName;
251 if (fasterCode) {
252 eFullFileName << path << "/microelec/Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat";
253 G4cout << "Faster code = true" << G4endl;
254 G4cout << "Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
255 }
256 else {
257 eFullFileName << path << "/microelec/Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat";
258 G4cout << "Faster code = false" << G4endl;
259 G4cout << "Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
260 }
261
262 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
263 if (!eDiffCrossSection)
264 {
265 std::stringstream ss;
266 ss << "Missing data " << eFullFileName.str().c_str();
267 std::string sortieString = ss.str();
268
269 if (fasterCode) G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
270 FatalException, sortieString.c_str());
271 else {
272 G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
273 FatalException, "Missing data file:/microelec/sigmadiff_inelastic_e_Si.dat");
274 }
275 }
276
277 // Clear the arrays for re-initialization case (MT mode)
278 // Octobre 22nd, 2014 - Melanie Raine
279 //Creating vectors of maps for DCS and Cumulated DCS for the current material.
280 //Each vector is storing one map for each shell.
281 G4bool isUsed1 = false;
282 vector<TriDimensionMap>* eDiffCrossSectionData =
283 new vector<TriDimensionMap>; //Storage of [IncidentEnergy, TransfEnergy, DCS values], used in slower code
284 vector<TriDimensionMap>* eNrjTransfData =
285 new vector<TriDimensionMap>; //Storage of possible transfer energies by shell
286 vector<VecMap>* eProbaShellMap = new vector<VecMap>; //Storage of the vectors containing all cumulated DCS values for an initial energy, by shell
287 vector<G4double>* eTdummyVec = new vector<G4double>; //Storage of incident energies for interpolation
288 VecMap* eVecm = new VecMap; //Transfered energy map for slower code
289
290 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); ++j) //Filling the map vectors with an empty map for each shell
291 {
292 eDiffCrossSectionData->push_back(TriDimensionMap());
293 eNrjTransfData->push_back(TriDimensionMap());
294 eProbaShellMap->push_back(VecMap());
295 }
296
297 eTdummyVec->push_back(0.);
298 while (!eDiffCrossSection.eof())
299 {
300 G4double tDummy; //incident energy
301 G4double eDummy; //transfered energy
302 eDiffCrossSection >> tDummy >> eDummy;
303 if (tDummy != eTdummyVec->back()) eTdummyVec->push_back(tDummy);
304
305 G4double tmp; //probability
306 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); ++j)
307 {
308 eDiffCrossSection >> tmp;
309 (*eDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
310
311 if (fasterCode)
312 {
313 (*eNrjTransfData)[j][tDummy][(*eDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
314 (*eProbaShellMap)[j][tDummy].push_back((*eDiffCrossSectionData)[j][tDummy][eDummy]);
315 }
316 else { // SI - only if eof is not reached !
317 (*eDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
318 (*eVecm)[tDummy].push_back(eDummy);
319 }
320 }
321 }
322 //
323 //G4cout << "add to material vector" << G4endl;
324
325 //Filing maps for the current material into the master maps
326 if (fasterCode) {
327 isUsed1 = true;
328 eNrjTransStorage[mat] = eNrjTransfData;
329 eProbaShellStorage[mat] = eProbaShellMap;
330 }
331 else {
332 eDiffDatatable[mat] = eDiffCrossSectionData;
333 eVecmStorage[mat] = eVecm;
334 }
335 eIncidentEnergyStorage[mat] = eTdummyVec;
336
337 //Cleanup support vectors
338 if (!isUsed1) {
339 delete eProbaShellMap;
340 delete eNrjTransfData;
341 } else {
342 delete eDiffCrossSectionData;
343 delete eVecm;
344 }
345 }
346
347 // *** PROTON
348 if (particle == protonDef)
349 {
350 // Cross section
351 G4String fileProton("Inelastic/" + modelName + "_sigma_inelastic_p_" + mat); G4cout << fileProton << G4endl;
353 tableP->LoadData(fileProton);
354 tableData->insert(make_pair(proton, tableP));
355
356 // DCS
357 std::ostringstream pFullFileName;
358 if (fasterCode) {
359 pFullFileName << path << "/microelec/Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat";
360 G4cout << "Faster code = true" << G4endl;
361 G4cout << "Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat" << G4endl;
362 }
363 else {
364 pFullFileName << path << "/microelec/Inelastic/" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat";
365 G4cout << "Faster code = false" << G4endl;
366 G4cout << "Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
367 }
368
369 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
370 if (!pDiffCrossSection)
371 {
372 if (fasterCode) G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
373 FatalException, "Missing data file:/microelec/sigmadiff_cumulated_inelastic_p_Si.dat");
374 else {
375 G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
376 FatalException, "Missing data file:/microelec/sigmadiff_inelastic_p_Si.dat");
377 }
378 }
379
380 //
381 // Clear the arrays for re-initialization case (MT mode)
382 // Octobre 22nd, 2014 - Melanie Raine
383 //Creating vectors of maps for DCS and Cumulated DCS for the current material.
384 //Each vector is storing one map for each shell.
385
386 G4bool isUsed1 = false;
387 vector<TriDimensionMap>* pDiffCrossSectionData =
388 new vector<TriDimensionMap>; //Storage of [IncidentEnergy, TransfEnergy, DCS values], used in slower code
389 vector<TriDimensionMap>* pNrjTransfData =
390 new vector<TriDimensionMap>; //Storage of possible transfer energies by shell
391 vector<VecMap>* pProbaShellMap =
392 new vector<VecMap>; //Storage of the vectors containing all cumulated DCS values for an initial energy, by shell
393 vector<G4double>* pTdummyVec =
394 new vector<G4double>; //Storage of incident energies for interpolation
395 VecMap* pVecm = new VecMap; //Transfered energy map for slower code
396
397 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); ++j)
398 //Filling the map vectors with an empty map for each shell
399 {
400 pDiffCrossSectionData->push_back(TriDimensionMap());
401 pNrjTransfData->push_back(TriDimensionMap());
402 pProbaShellMap->push_back(VecMap());
403 }
404
405 pTdummyVec->push_back(0.);
406 while (!pDiffCrossSection.eof())
407 {
408 G4double tDummy; //incident energy
409 G4double eDummy; //transfered energy
410 pDiffCrossSection >> tDummy >> eDummy;
411 if (tDummy != pTdummyVec->back()) pTdummyVec->push_back(tDummy);
412
413 G4double tmp; //probability
414 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); j++)
415 {
416 pDiffCrossSection >> tmp;
417 (*pDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
418 // ArrayofMaps[j] -> fill with 3DMap(incidentEnergy,
419 // 2Dmap (transferedEnergy,proba=tmp) ) -> fill map for shell j
420 // with proba for transfered energy eDummy
421
422 if (fasterCode)
423 {
424 (*pNrjTransfData)[j][tDummy][(*pDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
425 (*pProbaShellMap)[j][tDummy].push_back((*pDiffCrossSectionData)[j][tDummy][eDummy]);
426 }
427 else { // SI - only if eof is not reached !
428 (*pDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
429 (*pVecm)[tDummy].push_back(eDummy);
430 }
431 }
432 }
433
434 //Filing maps for the current material into the master maps
435 if (fasterCode) {
436 isUsed1 = true;
437 pNrjTransStorage[mat] = pNrjTransfData;
438 pProbaShellStorage[mat] = pProbaShellMap;
439 }
440 else {
441 pDiffDatatable[mat] = pDiffCrossSectionData;
442 pVecmStorage[mat] = pVecm;
443 }
444 pIncidentEnergyStorage[mat] = pTdummyVec;
445
446 //Cleanup support vectors
447 if (!isUsed1) {
448 delete pProbaShellMap;
449 delete pNrjTransfData;
450 } else {
451 delete pDiffCrossSectionData;
452 delete pVecm;
453 }
454 }
455 tableTCS[mat] = tableData;
456 }
457 if (particle==electronDef)
458 {
459 SetLowEnergyLimit(lowEnergyLimit[electron]);
460 SetHighEnergyLimit(highEnergyLimit[electron]);
461 }
462
463 if (particle==protonDef)
464 {
465 SetLowEnergyLimit(100*eV);
466 SetHighEnergyLimit(10*MeV);
467 }
468
469 if( verboseLevel>1 )
470 {
471 G4cout << "MicroElec Inelastic model is initialized " << G4endl
472 << "Energy range: "
473 << LowEnergyLimit() / keV << " keV - "
474 << HighEnergyLimit() / MeV << " MeV for "
475 << particle->GetParticleName()
476 << " with mass (amu) " << particle->GetPDGMass()/proton_mass_c2
477 << " and charge " << particle->GetPDGCharge()
478 << G4endl << G4endl ;
479 }
480
481 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
482
483 fParticleChangeForGamma = GetParticleChangeForGamma();
484 isInitialised = true;
485}
const char * G4FindDataDir(const char *)
bool G4bool
Definition G4Types.hh:86
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
G4bool LoadData(const G4String &argFileName) override
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 597 of file G4MicroElecInelasticModel_new.cc.

602{
603
604 if (verboseLevel > 3)
605 G4cout << "Calling SampleSecondaries() of G4MicroElecInelasticModel" << G4endl;
606
607 G4int pdg = particle->GetParticleDefinition()->GetPDGEncoding();
608 G4double lowLim = currentMaterialStructure->GetInelasticModelLowLimit(pdg);
609 G4double highLim = currentMaterialStructure->GetInelasticModelHighLimit(pdg);
610
611 G4double ekin = particle->GetKineticEnergy();
612 G4double k = ekin ;
613
614 G4ParticleDefinition* PartDef = particle->GetDefinition();
615 const G4String& particleName = PartDef->GetParticleName();
616 G4String nameLocal2 = particleName ;
617 G4double particleMass = particle->GetDefinition()->GetPDGMass();
618 G4double originalMass = particleMass; // a passer en argument dans samplesecondaryenergy pour évaluer correctement Qmax
619 G4int originalZ = particle->GetDefinition()->GetAtomicNumber();
620
621 if (particleMass > proton_mass_c2)
622 {
623 k *= proton_mass_c2/particleMass ;
624 PartDef = G4Proton::ProtonDefinition();
625 nameLocal2 = "proton" ;
626 }
627
628 if (k >= lowLim && k < highLim)
629 {
630 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
631 G4double totalEnergy = ekin + particleMass;
632 G4double pSquare = ekin * (totalEnergy + particleMass);
633 G4double totalMomentum = std::sqrt(pSquare);
634
635 G4int Shell = 1;
636
637 Shell = RandomSelect(k,nameLocal2,originalMass, originalZ);
638
639 G4double bindingEnergy = currentMaterialStructure->Energy(Shell);
640 G4double limitEnergy = currentMaterialStructure->GetLimitEnergy(Shell);
641
642 G4bool weaklyBound = currentMaterialStructure->IsShellWeaklyBound(Shell);
643
644 if (verboseLevel > 3)
645 {
646 G4cout << "---> Kinetic energy (eV)=" << k/eV << G4endl ;
647 G4cout << "Shell: " << Shell << ", energy: " << bindingEnergy/eV << G4endl;
648 }
649
650
651 if (k<limitEnergy) {
652 if (weaklyBound && k > currentMaterialStructure->GetEnergyGap()) {
653 limitEnergy = currentMaterialStructure->GetEnergyGap();
654 }
655 else return; }
656
657 // sample deexcitation
658
659 std::size_t secNumberInit = 0; // need to know at a certain point the energy of secondaries
660 std::size_t secNumberFinal = 0; // So I'll make the difference and then sum the energies
661
662 // G4cout << currentMaterial << G4endl;
663 G4int Z = currentMaterialStructure->GetZ(Shell);
664 G4int shellEnum = currentMaterialStructure->GetEADL_Enumerator(Shell);
665 if (currentMaterialStructure->IsShellWeaklyBound(Shell)) { shellEnum = -1; }
666
667 if(fAtomDeexcitation && shellEnum >=0)
668 {
669 // G4cout << "enter if deex and shell 0" << G4endl;
671 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
672 secNumberInit = fvect->size();
673 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
674 secNumberFinal = fvect->size();
675 }
676
677 G4double secondaryKinetic=-1000*eV;
678 SEFromFermiLevel = false;
679 if (!fasterCode)
680 {
681 secondaryKinetic = RandomizeEjectedElectronEnergy(PartDef, k, Shell, originalMass, originalZ);
682 }
683 else
684 {
685 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(PartDef, k, Shell) ;
686 }
687
688 if (verboseLevel > 3)
689 {
690 G4cout << "Ionisation process" << G4endl;
691 G4cout << "Shell: " << Shell << " Kin. energy (eV)=" << k/eV
692 << " Sec. energy (eV)=" << secondaryKinetic/eV << G4endl;
693 }
694 G4ThreeVector deltaDirection =
695 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
696 Z, Shell,
697 couple->GetMaterial());
698
700 {
701 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
702
703 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
704 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
705 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
706 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
707 finalPx /= finalMomentum;
708 finalPy /= finalMomentum;
709 finalPz /= finalMomentum;
710
711 G4ThreeVector direction;
712 direction.set(finalPx,finalPy,finalPz);
713
714 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit());
715 }
716 else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection);
717
718 // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
719 G4double deexSecEnergy = 0;
720 for (std::size_t j=secNumberInit; j < secNumberFinal; ++j) {
721 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
722 }
723 // correction CI 12/01/2023 limit energy = gap
724 //if (SEFromFermiLevel) limitEnergy = currentMaterialStructure->GetEnergyGap();
725 //fParticleChangeForGamma->SetProposedKineticEnergy(ekin - secondaryKinetic - limitEnergy); //Ef = Ei-(Q-El)-El = Ei-Q
726 //fParticleChangeForGamma->ProposeLocalEnergyDeposit(limitEnergy - deexSecEnergy);
727
728 // correction CI 09/03/2022 limit energy = gap
729 //if (!SEFromFermiLevel && weaklyBound) limitEnergy += currentMaterialStructure->GetInitialEnergy();
730 //if (!SEFromFermiLevel && weaklyBound) limitEnergy += currentMaterialStructure->GetEnergyGap();
731 fParticleChangeForGamma->SetProposedKineticEnergy(ekin - secondaryKinetic-limitEnergy); //Ef = Ei-(Q-El)-El = Ei-Q
732 fParticleChangeForGamma->ProposeLocalEnergyDeposit(limitEnergy-deexSecEnergy);
733
734 if (secondaryKinetic>0)
735 {
736 G4DynamicParticle* dp = new G4DynamicParticle(G4Electron::Electron(), deltaDirection, secondaryKinetic); //Esec = Q-El
737 fvect->push_back(dp);
738 }
739 }
740}
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4int GetAtomicNumber() const
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)

◆ stepFunc()

G4double G4MicroElecInelasticModel_new::stepFunc ( G4double x)

Definition at line 1393 of file G4MicroElecInelasticModel_new.cc.

1393 {
1394 return (x < 0.) ? 1.0 : 0.0;
1395}

Referenced by vrkreussler().

◆ vrkreussler()

G4double G4MicroElecInelasticModel_new::vrkreussler ( G4double v,
G4double vF )

Definition at line 1399 of file G4MicroElecInelasticModel_new.cc.

1400{
1401 G4double r = vF*( std::pow(v/vF+1., 3.) - fabs(std::pow(v/vF-1., 3.))
1402 + 4.*(v/vF)*(v/vF) ) + stepFunc(v/vF-1.) * (3./2.*v/vF -
1403 4.*(v/vF)*(v/vF) + 3.*std::pow(v/vF, 3.)
1404 - 0.5*std::pow(v/vF, 5.));
1405 return r/(10.*v/vF);
1406}

Referenced by BKZ().


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