Geant4 11.4.0
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 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

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
121 SetAngularDistribution(new G4DeltaAngle());
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 *)

Referenced by G4MicroElecInelasticModel_new(), and operator=().

◆ ~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 1389 of file G4MicroElecInelasticModel_new.cc.

1390{
1391 // need atomic unit conversion
1392 G4double hbar = hbar_Planck, hbar2 = hbar*hbar, me = electron_mass_c2/c_squared, Ry = me*elm_coupling*elm_coupling/(2*hbar2);
1393 G4double hartree = 2*Ry, a0 = Bohr_radius, velocity = a0*hartree/hbar;
1395
1396 vp /= velocity;
1397
1398 G4double wp = Eplasmon/hartree;
1399 G4double a = std::pow(4./9./CLHEP::pi, 1./3.);
1400 G4double vF = std::pow(wp*wp/(3.*a*a*a), 1./3.);
1401 G4double c = 0.9;
1402 G4double vr = vrkreussler(vp /*in u.a*/, vF /*in u.a*/);
1403 G4double yr = vr/std::pow(Zp, 2./3.);
1404 G4double q = 0.;
1405 if(Zp==2) q = 1-exp(-c*vr/(Zp-5./16.));
1406 else q = 1.-exp(-c*(yr-0.07));
1407 G4double Neq = Zp*(1.-q);
1408 G4double l0 = 0.;
1409 if(Neq<=2) l0 = 3./(Zp-0.3*(Neq-1))/2.;
1410 else l0 = 0.48*std::pow(Neq, 2./3.)/(Zp-Neq/7.);
1411 if(Zp==2) c = 1.0;
1412 else c = 3./2.;
1413 return Zp*(q + c*(1.-q)/vF/vF/2.0 * log(1.+std::pow(2.*l0*vF,2.)));
1414}
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 1356 of file G4MicroElecInelasticModel_new.cc.

1356 {
1357 G4double v1i = ComputeRelativistVelocity(T1i, M1);
1358 G4double v2i = ComputeRelativistVelocity(T2i, M2);
1359
1360 G4double v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*-1*v2i;
1361 G4double vtransfer2a = v2f*v2f-v2i*v2i;
1362
1363 v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*v2i;
1364 G4double vtransfer2b = v2f*v2f-v2i*v2i;
1365
1366 G4double vtransfer2 = std::max(vtransfer2a, vtransfer2b);
1367 return 0.5*M2*vtransfer2;
1368}

◆ ComputeRelativistVelocity()

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

Definition at line 1349 of file G4MicroElecInelasticModel_new.cc.

1349 {
1350 G4double x = E/mass;
1351 return c_light*std::sqrt(x*(x + 2.0))/(x + 1.0);
1352}

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 480 of file G4MicroElecInelasticModel_new.cc.

485{
486 if (verboseLevel > 3) G4cout << "Calling CrossSectionPerVolume() of G4MicroElecInelasticModel" << G4endl;
487
488 G4double density = material->GetTotNbOfAtomsPerVolume();
489 currentMaterial = material->GetName().substr(3, material->GetName().size());
490
491 MapStructure::iterator structPos;
492 structPos = tableMaterialsStructures.find(currentMaterial);
493
494 // Calculate total cross section for model
495 TCSMap::iterator tablepos;
496 tablepos = tableTCS.find(currentMaterial);
497
498 if (tablepos == tableTCS.end() )
499 {
500 G4String str = "Material ";
501 str += currentMaterial + " TCS Table not found!";
502 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
503 return 0;
504 }
505 else if(structPos == tableMaterialsStructures.end())
506 {
507 G4String str = "Material ";
508 str += currentMaterial + " Structure not found!";
509 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
510 return 0;
511 }
512 else {
513 MapData* tableData = tablepos->second;
514 currentMaterialStructure = structPos->second;
515
516 G4double sigma = 0;
517
518 const G4String& particleName = particleDefinition->GetParticleName();
519 G4String nameLocal = particleName;
520 G4int pdg = particleDefinition->GetPDGEncoding();
521 G4int Z = particleDefinition->GetAtomicNumber();
522
523 G4double Zeff = 1.0, Zeff2 = Zeff*Zeff;
524 G4double Mion_c2 = particleDefinition->GetPDGMass();
525
526 if (Mion_c2 > proton_mass_c2)
527 {
528 ekin *= proton_mass_c2 / Mion_c2;
529 nameLocal = "proton";
530 }
531
532 G4double lowLim = currentMaterialStructure->GetInelasticModelLowLimit(pdg);
533 G4double highLim = currentMaterialStructure->GetInelasticModelHighLimit(pdg);
534
535 if (ekin >= lowLim && ekin < highLim)
536 {
537 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
538 pos = tableData->find(nameLocal); //find particle type
539
540 if (pos != tableData->end())
541 {
542 G4MicroElecCrossSectionDataSet_new* table = pos->second;
543 if (table != 0)
544 {
545 sigma = table->FindValue(ekin);
546
547 if (Mion_c2 > proton_mass_c2) {
548 sigma = 0.;
549 for (G4int i = 0; i < currentMaterialStructure->NumberOfLevels(); i++) {
550 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
551 Zeff2 = Zeff*Zeff;
552 sigma += Zeff2*table->FindShellValue(ekin, i);
553 // il faut utiliser le ekin mis à l'echelle pour chercher la bonne
554 // valeur dans les tables proton
555 }
556 }
557 else {
558 sigma = table->FindValue(ekin);
559 }
560 }
561 }
562 else
563 {
564 G4Exception("G4MicroElecInelasticModel_new::CrossSectionPerVolume",
565 "em0002", FatalException,
566 "Model not applicable to particle type.");
567 }
568 }
569 else
570 {
571 return 1 / DBL_MAX;
572 }
573
574 if (verboseLevel > 3)
575 {
576 G4cout << "---> Kinetic energy (eV)=" << ekin / eV << G4endl;
577 G4cout << " - Cross section per Si atom (cm^2)=" << sigma / cm2 << G4endl;
578 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) << G4endl;
579 }
580
581 return (sigma)*density;
582 }
583}
@ 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 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 1084 of file G4MicroElecInelasticModel_new.cc.

1089{
1090 G4double sigma = 0.;
1091
1092 if (energyTransfer >= currentMaterialStructure->GetLimitEnergy(LevelIndex))
1093 {
1094 G4double valueT1 = 0;
1095 G4double valueT2 = 0;
1096 G4double valueE21 = 0;
1097 G4double valueE22 = 0;
1098 G4double valueE12 = 0;
1099 G4double valueE11 = 0;
1100 G4double xs11 = 0;
1101 G4double xs12 = 0;
1102 G4double xs21 = 0;
1103 G4double xs22 = 0;
1104
1105 if (particleDefinition == G4Electron::ElectronDefinition())
1106 {
1107
1108 dataDiffCSMap::iterator iterator_Proba;
1109 iterator_Proba = eDiffDatatable.find(currentMaterial);
1110
1111 incidentEnergyMap::iterator iterator_Nrj;
1112 iterator_Nrj = eIncidentEnergyStorage.find(currentMaterial);
1113
1114 TranfEnergyMap::iterator iterator_TransfNrj;
1115 iterator_TransfNrj = eVecmStorage.find(currentMaterial);
1116
1117 if (iterator_Proba != eDiffDatatable.end() && iterator_Nrj != eIncidentEnergyStorage.end()
1118 && iterator_TransfNrj!= eVecmStorage.end())
1119 {
1120 vector<TriDimensionMap>* eDiffCrossSectionData = (iterator_Proba->second);
1121 vector<G4double>* eTdummyVec = iterator_Nrj->second; //Incident energies for interpolation
1122 VecMap* eVecm = iterator_TransfNrj->second;
1123
1124 // k should be in eV and energy transfer eV also
1125 auto t2 = std::upper_bound(eTdummyVec->begin(), eTdummyVec->end(), k);
1126 auto t1 = t2 - 1;
1127 // SI : the following condition avoids situations where energyTransfer >last vector element
1128 if (energyTransfer <= (*eVecm)[(*t1)].back() && energyTransfer <= (*eVecm)[(*t2)].back())
1129 {
1130 auto e12 = std::upper_bound((*eVecm)[(*t1)].begin(), (*eVecm)[(*t1)].end(), energyTransfer);
1131 auto e11 = e12 - 1;
1132 auto e22 = std::upper_bound((*eVecm)[(*t2)].begin(), (*eVecm)[(*t2)].end(), energyTransfer);
1133 auto e21 = e22 - 1;
1134
1135 valueT1 = *t1;
1136 valueT2 = *t2;
1137 valueE21 = *e21;
1138 valueE22 = *e22;
1139 valueE12 = *e12;
1140 valueE11 = *e11;
1141
1142 xs11 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1143 xs12 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1144 xs21 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1145 xs22 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1146 }
1147 }
1148 else {
1149 G4String str = "Material ";
1150 str += currentMaterial + " not found!";
1151 G4Exception("G4MicroElecDielectricModels::DifferentialCrossSection", "em0002", FatalException, str);
1152 }
1153 }
1154
1155 if (particleDefinition == G4Proton::ProtonDefinition())
1156 {
1157 dataDiffCSMap::iterator iterator_Proba;
1158 iterator_Proba = pDiffDatatable.find(currentMaterial);
1159
1160 incidentEnergyMap::iterator iterator_Nrj;
1161 iterator_Nrj = pIncidentEnergyStorage.find(currentMaterial);
1162
1163 TranfEnergyMap::iterator iterator_TransfNrj;
1164 iterator_TransfNrj = pVecmStorage.find(currentMaterial);
1165
1166 if (iterator_Proba != pDiffDatatable.end() && iterator_Nrj != pIncidentEnergyStorage.end()
1167 && iterator_TransfNrj != pVecmStorage.end())
1168 {
1169 vector<TriDimensionMap>* pDiffCrossSectionData = (iterator_Proba->second);
1170 vector<G4double>* pTdummyVec = iterator_Nrj->second; //Incident energies for interpolation
1171 VecMap* pVecm = iterator_TransfNrj->second;
1172
1173 // k should be in eV and energy transfer eV also
1174 auto t2 =
1175 std::upper_bound(pTdummyVec->begin(), pTdummyVec->end(), k);
1176 auto t1 = t2 - 1;
1177 if (energyTransfer <= (*pVecm)[(*t1)].back() && energyTransfer <= (*pVecm)[(*t2)].back())
1178 {
1179 auto e12 = std::upper_bound((*pVecm)[(*t1)].begin(), (*pVecm)[(*t1)].end(), energyTransfer);
1180 auto e11 = e12 - 1;
1181 auto e22 = std::upper_bound((*pVecm)[(*t2)].begin(), (*pVecm)[(*t2)].end(), energyTransfer);
1182 auto e21 = e22 - 1;
1183
1184 valueT1 = *t1;
1185 valueT2 = *t2;
1186 valueE21 = *e21;
1187 valueE22 = *e22;
1188 valueE12 = *e12;
1189 valueE11 = *e11;
1190
1191 xs11 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1192 xs12 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1193 xs21 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1194 xs22 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1195 }
1196 }
1197 else {
1198 G4String str = "Material ";
1199 str += currentMaterial + " not found!";
1200 G4Exception("G4MicroElecDielectricModels::DifferentialCrossSection", "em0002", FatalException, str);
1201 }
1202 }
1203
1204 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
1205 if (xsProduct != 0.)
1206 {
1207 sigma = QuadInterpolator( valueE11, valueE12,
1208 valueE21, valueE22,
1209 xs11, xs12,
1210 xs21, xs22,
1211 valueT1, valueT2,
1212 k, energyTransfer);
1213 }
1214
1215 }
1216
1217 return sigma;
1218}
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
201 if (verboseLevel > 3)
202 G4cout << "Calling G4MicroElecInelasticModel_new::Initialise()" << G4endl;
203
204 const char* path = G4FindDataDir("G4LEDATA");
205 if (!path)
206 {
207 G4Exception("G4MicroElecElasticModel_new::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
208 return;
209 }
210
211 G4String modelName = "mermin";
212 G4cout << "****************************" << G4endl;
213 G4cout << modelName << " model loaded !" << G4endl;
214
215 // Energy limits
216 G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
217 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
218 G4String electron = electronDef->GetParticleName();
219 G4String proton = protonDef->GetParticleName();
220
221 G4double scaleFactor = 1.0;
222
223 // *** ELECTRON
224 lowEnergyLimit[electron] = 2 * eV;
225 highEnergyLimit[electron] = 10.0 * MeV;
226
227 // Cross section
228 G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
229 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
230
231 for (G4int i = 0; i < numOfCouples; ++i) {
232 const G4Material* material = theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
233 G4cout << "Material " << i + 1 << " / " << numOfCouples << " : " << material->GetName() << G4endl;
234 if (material->GetName() == "Vacuum") continue;
235 G4String mat = material->GetName().substr(3, material->GetName().size());
236 MapData* tableData = new MapData;
237 currentMaterialStructure = new G4MicroElecMaterialStructure(mat);
238
239 tableMaterialsStructures[mat] = currentMaterialStructure;
240 if (particle == electronDef) {
241 //TCS
242 G4String fileElectron("Inelastic/" + modelName + "_sigma_inelastic_e-_" + mat);
243 G4cout << fileElectron << G4endl;
244 G4MicroElecCrossSectionDataSet_new* tableE = new G4MicroElecCrossSectionDataSet_new(new G4LogLogInterpolation, MeV, scaleFactor);
245 tableE->LoadData(fileElectron);
246 tableData->insert(make_pair(electron, tableE));
247
248 // DCS
249 std::ostringstream eFullFileName;
250 if (fasterCode) {
251 eFullFileName << path << "/microelec/Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat";
252 G4cout << "Faster code = true" << G4endl;
253 G4cout << "Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
254 }
255 else {
256 eFullFileName << path << "/microelec/Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat";
257 G4cout << "Faster code = false" << G4endl;
258 G4cout << "Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
259 }
260
261 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
262 if (!eDiffCrossSection)
263 {
264 std::stringstream ss;
265 ss << "Missing data " << eFullFileName.str().c_str();
266 std::string sortieString = ss.str();
267
268 if (fasterCode) G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
269 FatalException, sortieString.c_str());
270 else {
271 G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
272 FatalException, "Missing data file:/microelec/sigmadiff_inelastic_e_Si.dat");
273 }
274 }
275
276 // Clear the arrays for re-initialization case (MT mode)
277 // Octobre 22nd, 2014 - Melanie Raine
278 //Creating vectors of maps for DCS and Cumulated DCS for the current material.
279 //Each vector is storing one map for each shell.
280 G4bool isUsed1 = false;
281 vector<TriDimensionMap>* eDiffCrossSectionData = new vector<TriDimensionMap>; //Storage of [IncidentEnergy, TransfEnergy, DCS values], used in slower code
282 vector<TriDimensionMap>* eNrjTransfData = new vector<TriDimensionMap>; //Storage of possible transfer energies by shell
283 vector<VecMap>* eProbaShellMap = new vector<VecMap>; //Storage of the vectors containing all cumulated DCS values for an initial energy, by shell
284 vector<G4double>* eTdummyVec = new vector<G4double>; //Storage of incident energies for interpolation
285 VecMap* eVecm = new VecMap; //Transfered energy map for slower code
286
287 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); ++j) //Filling the map vectors with an empty map for each shell
288 {
289 eDiffCrossSectionData->push_back(TriDimensionMap());
290 eNrjTransfData->push_back(TriDimensionMap());
291 eProbaShellMap->push_back(VecMap());
292 }
293
294 eTdummyVec->push_back(0.);
295 while (!eDiffCrossSection.eof())
296 {
297 G4double tDummy; //incident energy
298 G4double eDummy; //transfered energy
299 eDiffCrossSection >> tDummy >> eDummy;
300 if (tDummy != eTdummyVec->back()) eTdummyVec->push_back(tDummy);
301
302 G4double tmp; //probability
303 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); ++j)
304 {
305 eDiffCrossSection >> tmp;
306 (*eDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
307
308 if (fasterCode)
309 {
310 (*eNrjTransfData)[j][tDummy][(*eDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
311 (*eProbaShellMap)[j][tDummy].push_back((*eDiffCrossSectionData)[j][tDummy][eDummy]);
312 }
313 else { // SI - only if eof is not reached !
314 (*eDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
315 (*eVecm)[tDummy].push_back(eDummy);
316 }
317 }
318 }
319 //G4cout << "add to material vector" << G4endl;
320 //Filing maps for the current material into the master maps
321 if (fasterCode) {
322 isUsed1 = true;
323 eNrjTransStorage[mat] = eNrjTransfData;
324 eProbaShellStorage[mat] = eProbaShellMap;
325 }
326 else {
327 eDiffDatatable[mat] = eDiffCrossSectionData;
328 eVecmStorage[mat] = eVecm;
329 }
330 eIncidentEnergyStorage[mat] = eTdummyVec;
331
332 //Cleanup support vectors
333 if (!isUsed1) {
334 delete eProbaShellMap;
335 delete eNrjTransfData;
336 }
337 else {
338 delete eDiffCrossSectionData;
339 delete eVecm;
340 }
341 }
342
343 // *** PROTON
344 if (particle == protonDef)
345 {
346 // Cross section
347 G4String fileProton("Inelastic/" + modelName + "_sigma_inelastic_p_" + mat); G4cout << fileProton << G4endl;
348 G4MicroElecCrossSectionDataSet_new* tableP = new G4MicroElecCrossSectionDataSet_new(new G4LogLogInterpolation, MeV, scaleFactor);
349 tableP->LoadData(fileProton);
350 tableData->insert(make_pair(proton, tableP));
351
352 // DCS
353 std::ostringstream pFullFileName;
354 if (fasterCode) {
355 pFullFileName << path << "/microelec/Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat";
356 G4cout << "Faster code = true" << G4endl;
357 G4cout << "Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat" << G4endl;
358 }
359 else {
360 pFullFileName << path << "/microelec/Inelastic/" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat";
361 G4cout << "Faster code = false" << G4endl;
362 G4cout << "Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
363 }
364
365 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
366 if (!pDiffCrossSection)
367 {
368 if (fasterCode) G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
369 FatalException, "Missing data file:/microelec/sigmadiff_cumulated_inelastic_p_Si.dat");
370 else {
371 G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
372 FatalException, "Missing data file:/microelec/sigmadiff_inelastic_p_Si.dat");
373 }
374 }
375 // Clear the arrays for re-initialization case (MT mode)
376 // Octobre 22nd, 2014 - Melanie Raine
377 // Creating vectors of maps for DCS and Cumulated DCS for the current material.
378 // Each vector is storing one map for each shell.
379
380 G4bool isUsed1 = false;
381 vector<TriDimensionMap>* pDiffCrossSectionData =
382 new vector<TriDimensionMap>; //Storage of [IncidentEnergy, TransfEnergy, DCS values], used in slower code
383 vector<TriDimensionMap>* pNrjTransfData =
384 new vector<TriDimensionMap>; //Storage of possible transfer energies by shell
385 vector<VecMap>* pProbaShellMap =
386 new vector<VecMap>; //Storage of the vectors containing all cumulated DCS values for an initial energy, by shell
387 vector<G4double>* pTdummyVec =
388 new vector<G4double>; //Storage of incident energies for interpolation
389 VecMap* pVecm = new VecMap; //Transfered energy map for slower code
390
391 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); ++j)
392 //Filling the map vectors with an empty map for each shell
393 {
394 pDiffCrossSectionData->push_back(TriDimensionMap());
395 pNrjTransfData->push_back(TriDimensionMap());
396 pProbaShellMap->push_back(VecMap());
397 }
398
399 pTdummyVec->push_back(0.);
400 while (!pDiffCrossSection.eof())
401 {
402 G4double tDummy; //incident energy
403 G4double eDummy; //transfered energy
404 pDiffCrossSection >> tDummy >> eDummy;
405 if (tDummy != pTdummyVec->back()) pTdummyVec->push_back(tDummy);
406
407 G4double tmp; //probability
408 for (G4int j = 0; j < currentMaterialStructure->NumberOfLevels(); j++)
409 {
410 pDiffCrossSection >> tmp;
411 (*pDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
412 // ArrayofMaps[j] -> fill with 3DMap(incidentEnergy,
413 // 2Dmap (transferedEnergy,proba=tmp) ) -> fill map for shell j
414 // with proba for transfered energy eDummy
415 if (fasterCode)
416 {
417 (*pNrjTransfData)[j][tDummy][(*pDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
418 (*pProbaShellMap)[j][tDummy].push_back((*pDiffCrossSectionData)[j][tDummy][eDummy]);
419 }
420 else { // SI - only if eof is not reached !
421 (*pDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
422 (*pVecm)[tDummy].push_back(eDummy);
423 }
424 }
425 }
426 //Filing maps for the current material into the master maps
427 if (fasterCode) {
428 isUsed1 = true;
429 pNrjTransStorage[mat] = pNrjTransfData;
430 pProbaShellStorage[mat] = pProbaShellMap;
431 }
432 else {
433 pDiffDatatable[mat] = pDiffCrossSectionData;
434 pVecmStorage[mat] = pVecm;
435 }
436 pIncidentEnergyStorage[mat] = pTdummyVec;
437
438 //Cleanup support vectors
439 if (!isUsed1) {
440 delete pProbaShellMap;
441 delete pNrjTransfData;
442 } else {
443 delete pDiffCrossSectionData;
444 delete pVecm;
445 }
446 }
447 tableTCS[mat] = tableData;
448 }
449 if (particle==electronDef)
450 {
451 SetLowEnergyLimit(lowEnergyLimit[electron]);
452 SetHighEnergyLimit(highEnergyLimit[electron]);
453 }
454
455 if (particle==protonDef)
456 {
457 SetLowEnergyLimit(100*eV);
458 SetHighEnergyLimit(10*MeV);
459 }
460
461 if( verboseLevel>1 )
462 {
463 G4cout << "MicroElec Inelastic model is initialized " << G4endl
464 << "Energy range: "
465 << LowEnergyLimit() / keV << " keV - "
466 << HighEnergyLimit() / MeV << " MeV for "
467 << particle->GetParticleName()
468 << " with mass (amu) " << particle->GetPDGMass()/proton_mass_c2
469 << " and charge " << particle->GetPDGCharge()
470 << G4endl << G4endl ;
471 }
472
473 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
474 fParticleChangeForGamma = GetParticleChangeForGamma();
475 isInitialised = true;
476}
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 587 of file G4MicroElecInelasticModel_new.cc.

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

◆ stepFunc()

G4double G4MicroElecInelasticModel_new::stepFunc ( G4double x)

Definition at line 1372 of file G4MicroElecInelasticModel_new.cc.

1372 {
1373 return (x < 0.) ? 1.0 : 0.0;
1374}

Referenced by vrkreussler().

◆ vrkreussler()

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

Definition at line 1378 of file G4MicroElecInelasticModel_new.cc.

1379{
1380 G4double r = vF*( std::pow(v/vF+1., 3.) - fabs(std::pow(v/vF-1., 3.))
1381 + 4.*(v/vF)*(v/vF) ) + stepFunc(v/vF-1.) * (3./2.*v/vF -
1382 4.*(v/vF)*(v/vF) + 3.*std::pow(v/vF, 3.)
1383 - 0.5*std::pow(v/vF, 5.));
1384 return r/(10.*v/vF);
1385}

Referenced by BKZ().


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