Geant4 11.1.1
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 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 *, 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma = nullptr
 
- 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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 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:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607

◆ ~G4MicroElecInelasticModel_new()

G4MicroElecInelasticModel_new::~G4MicroElecInelasticModel_new ( )
override

Definition at line 130 of file G4MicroElecInelasticModel_new.cc.

131{
132 // Cross section
133 // (0)
134 TCSMap::iterator pos2;
135 for (pos2 = tableTCS.begin(); pos2 != tableTCS.end(); ++pos2) {
136 MapData* tableData = pos2->second;
137 for (auto pos = tableData->begin(); pos != tableData->end(); ++pos)
138 {
140 delete table;
141 }
142 delete tableData;
143 }
144 tableTCS.clear();
145
146 dataDiffCSMap::iterator iterator_proba;
147 // (1)
148 for (iterator_proba = eNrjTransStorage.begin(); iterator_proba != eNrjTransStorage.end(); ++iterator_proba) {
149 vector<TriDimensionMap>* eNrjTransfData = iterator_proba->second;
150 eNrjTransfData->clear();
151 delete eNrjTransfData;
152 }
153 eNrjTransStorage.clear();
154
155 for (iterator_proba = pNrjTransStorage.begin(); iterator_proba != pNrjTransStorage.end(); ++iterator_proba) {
156 vector<TriDimensionMap>* pNrjTransfData = iterator_proba->second;
157 pNrjTransfData->clear();
158 delete pNrjTransfData;
159 }
160 pNrjTransStorage.clear();
161
162 // (2)
163 for (iterator_proba = eDiffDatatable.begin(); iterator_proba != eDiffDatatable.end(); ++iterator_proba) {
164 vector<TriDimensionMap>* eDiffCrossSectionData = iterator_proba->second;
165 eDiffCrossSectionData->clear();
166 delete eDiffCrossSectionData;
167 }
168 eDiffDatatable.clear();
169
170 for (iterator_proba = pDiffDatatable.begin(); iterator_proba != pDiffDatatable.end(); ++iterator_proba) {
171 vector<TriDimensionMap>* pDiffCrossSectionData = iterator_proba->second;
172 pDiffCrossSectionData->clear();
173 delete pDiffCrossSectionData;
174 }
175 pDiffDatatable.clear();
176
177 // (3)
178 dataProbaShellMap::iterator iterator_probaShell;
179
180 for (iterator_probaShell = eProbaShellStorage.begin(); iterator_probaShell != eProbaShellStorage.end(); ++iterator_probaShell) {
181 vector<VecMap>* eProbaShellMap = iterator_probaShell->second;
182 eProbaShellMap->clear();
183 delete eProbaShellMap;
184 }
185 eProbaShellStorage.clear();
186
187 for (iterator_probaShell = pProbaShellStorage.begin(); iterator_probaShell != pProbaShellStorage.end(); ++iterator_probaShell) {
188 vector<VecMap>* pProbaShellMap = iterator_probaShell->second;
189 pProbaShellMap->clear();
190 delete pProbaShellMap;
191 }
192 pProbaShellStorage.clear();
193
194 // (4)
195 TranfEnergyMap::iterator iterator_nrjtransf;
196 for (iterator_nrjtransf = eVecmStorage.begin(); iterator_nrjtransf != eVecmStorage.end(); ++iterator_nrjtransf) {
197 VecMap* eVecm = iterator_nrjtransf->second;
198 eVecm->clear();
199 delete eVecm;
200 }
201 eVecmStorage.clear();
202 for (iterator_nrjtransf = pVecmStorage.begin(); iterator_nrjtransf != pVecmStorage.end(); ++iterator_nrjtransf) {
203 VecMap* pVecm = iterator_nrjtransf->second;
204 pVecm->clear();
205 delete pVecm;
206 }
207 pVecmStorage.clear();
208
209 // (5)
210 incidentEnergyMap::iterator iterator_energy;
211 for (iterator_energy = eIncidentEnergyStorage.begin(); iterator_energy != eIncidentEnergyStorage.end(); ++iterator_energy) {
212 std::vector<G4double>* eTdummyVec = iterator_energy->second;
213 eTdummyVec->clear();
214 delete eTdummyVec;
215 }
216 eIncidentEnergyStorage.clear();
217
218 for (iterator_energy = pIncidentEnergyStorage.begin(); iterator_energy != pIncidentEnergyStorage.end(); ++iterator_energy) {
219 std::vector<G4double>* pTdummyVec = iterator_energy->second;
220 pTdummyVec->clear();
221 delete pTdummyVec;
222 }
223 pIncidentEnergyStorage.clear();
224
225 // (6)
226 MapStructure::iterator iterator_matStructure;
227 for (iterator_matStructure = tableMaterialsStructures.begin();
228 iterator_matStructure != tableMaterialsStructures.end(); ++iterator_matStructure) {
229 currentMaterialStructure = iterator_matStructure->second;
230 delete currentMaterialStructure;
231 }
232 tableMaterialsStructures.clear();
233 currentMaterialStructure = nullptr;
234}

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

1418{
1419 // need atomic unit conversion
1420 G4double hbar = hbar_Planck, hbar2 = hbar*hbar, me = electron_mass_c2/c_squared, Ry = me*elm_coupling*elm_coupling/(2*hbar2);
1421 G4double hartree = 2*Ry, a0 = Bohr_radius, velocity = a0*hartree/hbar;
1423
1424 vp /= velocity;
1425
1426 G4double wp = Eplasmon/hartree;
1427 G4double a = std::pow(4./9./CLHEP::pi, 1./3.);
1428 G4double vF = std::pow(wp*wp/(3.*a*a*a), 1./3.);
1429 G4double c = 0.9;
1430 G4double vr = vrkreussler(vp /*in u.a*/, vF /*in u.a*/);
1431 G4double yr = vr/std::pow(Zp, 2./3.);
1432 G4double q = 0.;
1433 if(Zp==2) q = 1-exp(-c*vr/(Zp-5./16.));
1434 else q = 1.-exp(-c*(yr-0.07));
1435 G4double Neq = Zp*(1.-q);
1436 G4double l0 = 0.;
1437 if(Neq<=2) l0 = 3./(Zp-0.3*(Neq-1))/2.;
1438 else l0 = 0.48*std::pow(Neq, 2./3.)/(Zp-Neq/7.);
1439 if(Zp==2) c = 1.0;
1440 else c = 3./2.;
1441 return Zp*(q + c*(1.-q)/vF/vF/2.0 * log(1.+std::pow(2.*l0*vF,2.)));
1442}
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 1384 of file G4MicroElecInelasticModel_new.cc.

1384 {
1385 G4double v1i = ComputeRelativistVelocity(T1i, M1);
1386 G4double v2i = ComputeRelativistVelocity(T2i, M2);
1387
1388 G4double v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*-1*v2i;
1389 G4double vtransfer2a = v2f*v2f-v2i*v2i;
1390
1391 v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*v2i;
1392 G4double vtransfer2b = v2f*v2f-v2i*v2i;
1393
1394 G4double vtransfer2 = std::max(vtransfer2a, vtransfer2b);
1395 return 0.5*M2*vtransfer2;
1396}

◆ ComputeRelativistVelocity()

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

Definition at line 1377 of file G4MicroElecInelasticModel_new.cc.

1377 {
1378 G4double x = E/mass;
1379 return c_light*std::sqrt(x*(x + 2.0))/(x + 1.0);
1380}

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

525{
526 if (verboseLevel > 3) G4cout << "Calling CrossSectionPerVolume() of G4MicroElecInelasticModel" << G4endl;
527
528 G4double density = material->GetTotNbOfAtomsPerVolume();
529 currentMaterial = material->GetName().substr(3, material->GetName().size());
530
531 MapStructure::iterator structPos;
532 structPos = tableMaterialsStructures.find(currentMaterial);
533
534 // Calculate total cross section for model
535 TCSMap::iterator tablepos;
536 tablepos = tableTCS.find(currentMaterial);
537
538 if (tablepos == tableTCS.end() )
539 {
540 G4String str = "Material ";
541 str += currentMaterial + " TCS Table not found!";
542 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
543 return 0;
544 }
545 else if(structPos == tableMaterialsStructures.end())
546 {
547 G4String str = "Material ";
548 str += currentMaterial + " Structure not found!";
549 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
550 return 0;
551 }
552 else {
553 MapData* tableData = tablepos->second;
554 currentMaterialStructure = structPos->second;
555
556 G4double sigma = 0;
557
558 const G4String& particleName = particleDefinition->GetParticleName();
559 G4String nameLocal = particleName;
560 G4int pdg = particleDefinition->GetPDGEncoding();
561 G4int Z = particleDefinition->GetAtomicNumber();
562
563 G4double Zeff = 1.0, Zeff2 = Zeff*Zeff;
564 G4double Mion_c2 = particleDefinition->GetPDGMass();
565
566 if (Mion_c2 > proton_mass_c2)
567 {
568 ekin *= proton_mass_c2 / Mion_c2;
569 nameLocal = "proton";
570 }
571
572 G4double lowLim = currentMaterialStructure->GetInelasticModelLowLimit(pdg);
573 G4double highLim = currentMaterialStructure->GetInelasticModelHighLimit(pdg);
574
575 if (ekin >= lowLim && ekin < highLim)
576 {
577 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
578 pos = tableData->find(nameLocal); //find particle type
579
580 if (pos != tableData->end())
581 {
583 if (table != 0)
584 {
585 sigma = table->FindValue(ekin);
586
587 if (Mion_c2 > proton_mass_c2) {
588 sigma = 0.;
589 for (G4int i = 0; i < currentMaterialStructure->NumberOfLevels(); i++) {
590 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
591 Zeff2 = Zeff*Zeff;
592 sigma += Zeff2*table->FindShellValue(ekin, i);
593// il faut utiliser le ekin mis à l'echelle pour chercher la bonne
594// valeur dans les tables proton
595
596 }
597 }
598 else {
599 sigma = table->FindValue(ekin);
600 }
601 }
602 }
603 else
604 {
605 G4Exception("G4MicroElecInelasticModel_new::CrossSectionPerVolume",
606 "em0002", FatalException,
607 "Model not applicable to particle type.");
608 }
609 }
610 else
611 {
612 return 1 / DBL_MAX;
613 }
614
615 if (verboseLevel > 3)
616 {
617 G4cout << "---> Kinetic energy (eV)=" << ekin / eV << G4endl;
618 G4cout << " - Cross section per Si atom (cm^2)=" << sigma / cm2 << G4endl;
619 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) << G4endl;
620 }
621
622 return (sigma)*density;}
623
624}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4String & GetName() const
Definition: G4Material.hh:172
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 1118 of file G4MicroElecInelasticModel_new.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 238 of file G4MicroElecInelasticModel_new.cc.

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

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

633{
634
635 if (verboseLevel > 3)
636 G4cout << "Calling SampleSecondaries() of G4MicroElecInelasticModel" << G4endl;
637
638 G4int pdg = particle->GetParticleDefinition()->GetPDGEncoding();
639 G4double lowLim = currentMaterialStructure->GetInelasticModelLowLimit(pdg);
640 G4double highLim = currentMaterialStructure->GetInelasticModelHighLimit(pdg);
641
642 G4double ekin = particle->GetKineticEnergy();
643 G4double k = ekin ;
644
645 G4ParticleDefinition* PartDef = particle->GetDefinition();
646 const G4String& particleName = PartDef->GetParticleName();
647 G4String nameLocal2 = particleName ;
648 G4double particleMass = particle->GetDefinition()->GetPDGMass();
649 G4double originalMass = particleMass; // a passer en argument dans samplesecondaryenergy pour évaluer correctement Qmax
650 G4int originalZ = particle->GetDefinition()->GetAtomicNumber();
651
652 if (particleMass > proton_mass_c2)
653 {
654 k *= proton_mass_c2/particleMass ;
655 PartDef = G4Proton::ProtonDefinition();
656 nameLocal2 = "proton" ;
657 }
658
659 if (k >= lowLim && k < highLim)
660 {
661 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
662 G4double totalEnergy = ekin + particleMass;
663 G4double pSquare = ekin * (totalEnergy + particleMass);
664 G4double totalMomentum = std::sqrt(pSquare);
665
666 G4int Shell = 1;
667
668 Shell = RandomSelect(k,nameLocal2,originalMass, originalZ);
669
670 G4double bindingEnergy = currentMaterialStructure->Energy(Shell);
671 G4double limitEnergy = currentMaterialStructure->GetLimitEnergy(Shell);
672
673 if (verboseLevel > 3)
674 {
675 G4cout << "---> Kinetic energy (eV)=" << k/eV << G4endl ;
676 G4cout << "Shell: " << Shell << ", energy: " << bindingEnergy/eV << G4endl;
677 }
678
679 // sample deexcitation
680
681 std::size_t secNumberInit = 0; // need to know at a certain point the energy of secondaries
682 std::size_t secNumberFinal = 0; // So I'll make the difference and then sum the energies
683
684 //SI: additional protection if tcs interpolation method is modified
685 //if (k<bindingEnergy) return;
686 if (k<limitEnergy) return;
687 // G4cout << currentMaterial << G4endl;
688 G4int Z = currentMaterialStructure->GetZ(Shell);
689 G4int shellEnum = currentMaterialStructure->GetEADL_Enumerator(Shell);
690 if (currentMaterialStructure->IsShellWeaklyBound(Shell)) { shellEnum = -1; }
691
692 if(fAtomDeexcitation && shellEnum >=0)
693 {
694 // G4cout << "enter if deex and shell 0" << G4endl;
696 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
697 secNumberInit = fvect->size();
698 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
699 secNumberFinal = fvect->size();
700 }
701
702 G4double secondaryKinetic=-1000*eV;
703 SEFromFermiLevel = false;
704 if (!fasterCode)
705 {
706 secondaryKinetic = RandomizeEjectedElectronEnergy(PartDef, k, Shell, originalMass, originalZ);
707 }
708 else
709 {
710 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(PartDef, k, Shell) ;
711 }
712
713 if (verboseLevel > 3)
714 {
715 G4cout << "Ionisation process" << G4endl;
716 G4cout << "Shell: " << Shell << " Kin. energy (eV)=" << k/eV
717 << " Sec. energy (eV)=" << secondaryKinetic/eV << G4endl;
718 }
719 G4ThreeVector deltaDirection =
720 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
721 Z, Shell,
722 couple->GetMaterial());
723
725 {
726 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
727
728 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
729 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
730 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
731 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
732 finalPx /= finalMomentum;
733 finalPy /= finalMomentum;
734 finalPz /= finalMomentum;
735
736 G4ThreeVector direction;
737 direction.set(finalPx,finalPy,finalPz);
738
740 }
742
743 // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
744 G4double deexSecEnergy = 0;
745 for (std::size_t j=secNumberInit; j < secNumberFinal; ++j) {
746 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
747 }
748 if (SEFromFermiLevel) limitEnergy = currentMaterialStructure->GetEnergyGap();
749 fParticleChangeForGamma->SetProposedKineticEnergy(ekin - secondaryKinetic - limitEnergy); //Ef = Ei-(Q-El)-El = Ei-Q
750 fParticleChangeForGamma->ProposeLocalEnergyDeposit(limitEnergy - deexSecEnergy);
751
752 if (secondaryKinetic>0)
753 {
754 G4DynamicParticle* dp = new G4DynamicParticle(G4Electron::Electron(), deltaDirection, secondaryKinetic); //Esec = Q-El
755 fvect->push_back(dp);
756 }
757 }
758}
G4AtomicShellEnumerator
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:93
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()
Definition: G4VEmModel.hh:600
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ stepFunc()

G4double G4MicroElecInelasticModel_new::stepFunc ( G4double  x)

Definition at line 1400 of file G4MicroElecInelasticModel_new.cc.

1400 {
1401 return (x < 0.) ? 1.0 : 0.0;
1402}

Referenced by vrkreussler().

◆ vrkreussler()

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

Definition at line 1406 of file G4MicroElecInelasticModel_new.cc.

1407{
1408 G4double r = vF*( std::pow(v/vF+1., 3.) - fabs(std::pow(v/vF-1., 3.))
1409 + 4.*(v/vF)*(v/vF) ) + stepFunc(v/vF-1.) * (3./2.*v/vF -
1410 4.*(v/vF)*(v/vF) + 3.*std::pow(v/vF, 3.)
1411 - 0.5*std::pow(v/vF, 5.));
1412 return r/(10.*v/vF);
1413}

Referenced by BKZ().

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MicroElecInelasticModel_new::fParticleChangeForGamma = nullptr
protected

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