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

#include <G4MicroElecElasticModel.hh>

+ Inheritance diagram for G4MicroElecElasticModel:

Public Member Functions

 G4MicroElecElasticModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="MicroElecElasticModel")
 
virtual ~G4MicroElecElasticModel ()
 
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
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
G4MicroElecElasticModeloperator= (const G4MicroElecElasticModel &right)=delete
 
 G4MicroElecElasticModel (const G4MicroElecElasticModel &)=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
 
- 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 51 of file G4MicroElecElasticModel.hh.

Constructor & Destructor Documentation

◆ G4MicroElecElasticModel() [1/2]

G4MicroElecElasticModel::G4MicroElecElasticModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "MicroElecElasticModel" 
)

Definition at line 49 of file G4MicroElecElasticModel.cc.

51 :G4VEmModel(nam),isInitialised(false)
52{
54
55 killBelowEnergy = 16.7 * eV; // Minimum e- energy for energy loss by excitation
56 lowEnergyLimit = 0 * eV;
57 lowEnergyLimitOfModel = 5 * eV; // The model lower energy is 5 eV
58 highEnergyLimit = 100. * MeV;
59 SetLowEnergyLimit(lowEnergyLimit);
60 SetHighEnergyLimit(highEnergyLimit);
61
62 verboseLevel= 0;
63 // Verbosity scale:
64 // 0 = nothing
65 // 1 = warning for energy non-conservation
66 // 2 = details of energy budget
67 // 3 = calculation of cross sections, file openings, sampling of atoms
68 // 4 = entering in methods
69
70 if( verboseLevel>0 )
71 {
72 G4cout << "MicroElec Elastic model is constructed " << G4endl
73 << "Energy range: "
74 << lowEnergyLimit / eV << " eV - "
75 << highEnergyLimit / MeV << " MeV"
76 << G4endl;
77 }
79}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753

◆ ~G4MicroElecElasticModel()

G4MicroElecElasticModel::~G4MicroElecElasticModel ( )
virtual

Definition at line 83 of file G4MicroElecElasticModel.cc.

84{
85 // For total cross section
86 for (auto & pos : tableData)
87 {
88 G4MicroElecCrossSectionDataSet* table = pos.second;
89 delete table;
90 }
91
92 // For final state
93 eVecm.clear();
94}

◆ G4MicroElecElasticModel() [2/2]

G4MicroElecElasticModel::G4MicroElecElasticModel ( const G4MicroElecElasticModel )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 197 of file G4MicroElecElasticModel.cc.

202{
203 if (verboseLevel > 3)
204 G4cout << "Calling CrossSectionPerVolume() of G4MicroElecElasticModel" << G4endl;
205
206 // Calculate total cross section for model
207 G4double sigma=0;
208 G4double density = material->GetTotNbOfAtomsPerVolume();
209
210 if (material == nistSi || material->GetBaseMaterial() == nistSi)
211 {
212 const G4String& particleName = p->GetParticleName();
213
214 if (ekin < highEnergyLimit)
215 {
216 //SI : XS must not be zero otherwise sampling of secondaries method ignored
217 if (ekin < killBelowEnergy) return DBL_MAX;
218 //
219
220 auto pos = tableData.find(particleName);
221 if (pos != tableData.end())
222 {
223 G4MicroElecCrossSectionDataSet* table = pos->second;
224 if (table != nullptr)
225 {
226 sigma = table->FindValue(ekin);
227 }
228 }
229 else
230 {
231 G4Exception("G4MicroElecElasticModel::ComputeCrossSectionPerVolume","em0002",
232 FatalException,"Model not applicable to particle type.");
233 }
234 }
235
236 if (verboseLevel > 3)
237 {
238 G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
239 G4cout << " - Cross section per Si atom (cm^2)=" << sigma/cm/cm << G4endl;
240 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density/(1./cm) << G4endl;
241 }
242 }
243 return sigma*density;
244}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:228
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
G4double FindValue(G4double e, G4int componentId=0) const override
const G4String & GetParticleName() const
#define DBL_MAX
Definition: templates.hh:62

◆ GetKillBelowThreshold()

G4double G4MicroElecElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 74 of file G4MicroElecElasticModel.hh.

74{ return killBelowEnergy; }

◆ Initialise()

void G4MicroElecElasticModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 98 of file G4MicroElecElasticModel.cc.

100{
101 if (verboseLevel > 3)
102 G4cout << "Calling G4MicroElecElasticModel::Initialise()" << G4endl;
103
104 // Energy limits
105 if (LowEnergyLimit() < lowEnergyLimit)
106 {
107 G4cout << "G4MicroElecElasticModel: low energy limit increased from " <<
108 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
109 SetLowEnergyLimit(lowEnergyLimit);
110 }
111
112 if (HighEnergyLimit() > highEnergyLimit)
113 {
114 G4cout << "G4MicroElecElasticModel: high energy limit decreased from " <<
115 HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
116 SetHighEnergyLimit(highEnergyLimit);
117 }
118
119 // Reading of data files
120
121 G4double scaleFactor = 1e-18 * cm * cm;
122 G4String fileElectron("microelec/sigma_elastic_e_Si");
123
126
127 // For total cross section
128 electron = electronDef->GetParticleName();
129 tableFile[electron] = fileElectron;
130
132 tableE->LoadData(fileElectron);
133 tableData[electron] = tableE;
134
135 // For final state
136 const char* path = G4FindDataDir("G4LEDATA");
137
138 if (!path)
139 {
140 G4Exception("G4MicroElecElasticModel::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
141 return;
142 }
143
144 std::ostringstream eFullFileName;
145 eFullFileName << path << "/microelec/sigmadiff_cumulated_elastic_e_Si.dat";
146 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
147
148 if (!eDiffCrossSection)
149 G4Exception("G4MicroElecElasticModel::Initialise","em0003",FatalException,
150 "Missing data file: /microelec/sigmadiff_cumulated_elastic_e_Si.dat");
151
152 // Added clear for MT
153 eTdummyVec.clear();
154 eVecm.clear();
155 eDiffCrossSectionData.clear();
156
157 //
158 eTdummyVec.push_back(0.);
159
160 while(!eDiffCrossSection.eof())
161 {
162 double tDummy;
163 double eDummy;
164 eDiffCrossSection>>tDummy>>eDummy;
165
166 if (tDummy != eTdummyVec.back())
167 {
168 eTdummyVec.push_back(tDummy);
169 eVecm[tDummy].push_back(0.);
170 }
171
172 eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
173
174 if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
175 }
176 // End final state
177
178 if (verboseLevel > 2)
179 G4cout << "Loaded cross section files for MicroElec Elastic model" << G4endl;
180
181 if( verboseLevel>0 )
182 {
183 G4cout << "MicroElec Elastic model is initialized " << G4endl
184 << "Energy range: "
185 << LowEnergyLimit() / eV << " eV - "
186 << HighEnergyLimit() / MeV << " MeV"
187 << G4endl;
188 }
189
190 if (isInitialised) { return; }
192 isInitialised = true;
193}
const char * G4FindDataDir(const char *)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
G4bool LoadData(const G4String &argFileName) override
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 248 of file G4MicroElecElasticModel.cc.

253{
254
255 if (verboseLevel > 3)
256 G4cout << "Calling SampleSecondaries() of G4MicroElecElasticModel" << G4endl;
257
258 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
259
260 if (electronEnergy0 < killBelowEnergy)
261 {
265 return ;
266 }
267
268 if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
269 {
270 G4double cosTheta = RandomizeCosTheta(electronEnergy0);
271 G4double phi = 2. * pi * G4UniformRand();
272 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
273 G4ThreeVector xVers = zVers.orthogonal();
274 G4ThreeVector yVers = zVers.cross(xVers);
275
276 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
277 G4double yDir = xDir;
278 xDir *= std::cos(phi);
279 yDir *= std::sin(phi);
280
281 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
282
285 }
286}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4double pi

◆ SetKillBelowThreshold()

void G4MicroElecElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 125 of file G4MicroElecElasticModel.hh.

126{
127 killBelowEnergy = threshold;
128
129 if (threshold < 5*CLHEP::eV)
130 {
131 G4Exception ("*** WARNING : the G4MicroElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
132 threshold = 5*CLHEP::eV;
133 }
134
135}
@ JustWarning

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MicroElecElasticModel::fParticleChangeForGamma
protected

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