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

#include <G4MuElecElasticModel.hh>

+ Inheritance diagram for G4MuElecElasticModel:

Public Member Functions

 G4MuElecElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="MuElecElasticModel")
 
virtual ~G4MuElecElasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 G4MuElecElasticModel.hh.

Constructor & Destructor Documentation

◆ G4MuElecElasticModel()

G4MuElecElasticModel::G4MuElecElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MuElecElasticModel" 
)

Definition at line 48 of file G4MuElecElasticModel.cc.

50:G4VEmModel(nam),isInitialised(false)
51{
53
54 killBelowEnergy = 16.7 * eV; // Minimum e- energy for energy loss by excitation
55 lowEnergyLimit = 0 * eV;
56 lowEnergyLimitOfModel = 5 * eV; // The model lower energy is 5 eV
57 highEnergyLimit = 100. * MeV;
58 SetLowEnergyLimit(lowEnergyLimit);
59 SetHighEnergyLimit(highEnergyLimit);
60
61 verboseLevel= 0;
62 // Verbosity scale:
63 // 0 = nothing
64 // 1 = warning for energy non-conservation
65 // 2 = details of energy budget
66 // 3 = calculation of cross sections, file openings, sampling of atoms
67 // 4 = entering in methods
68
69 if( verboseLevel>0 )
70 {
71 G4cout << "MuElec Elastic model is constructed " << G4endl
72 << "Energy range: "
73 << lowEnergyLimit / eV << " eV - "
74 << highEnergyLimit / keV << " keV"
75 << G4endl;
76 }
78}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4MuElecElasticModel()

G4MuElecElasticModel::~G4MuElecElasticModel ( )
virtual

Definition at line 82 of file G4MuElecElasticModel.cc.

83{
84 // For total cross section
85
86 std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
87 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
88 {
89 G4MuElecCrossSectionDataSet* table = pos->second;
90 delete table;
91 }
92
93 // For final state
94
95 eVecm.clear();
96
97}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4MuElecElasticModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 207 of file G4MuElecElasticModel.cc.

212{
213 if (verboseLevel > 3)
214 G4cout << "Calling CrossSectionPerVolume() of G4MuElecElasticModel" << G4endl;
215
216 // Calculate total cross section for model
217
218 G4double sigma=0;
219
220 G4double density = material->GetTotNbOfAtomsPerVolume();
221
222 if (material == nistSi || material->GetBaseMaterial() == nistSi)
223 {
224 const G4String& particleName = p->GetParticleName();
225
226 if (ekin < highEnergyLimit)
227 {
228 //SI : XS must not be zero otherwise sampling of secondaries method ignored
229 if (ekin < lowEnergyLimitOfModel) ekin = lowEnergyLimitOfModel;
230 //
231
232 std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
233 pos = tableData.find(particleName);
234
235 if (pos != tableData.end())
236 {
237 G4MuElecCrossSectionDataSet* table = pos->second;
238 if (table != 0)
239 {
240 sigma = table->FindValue(ekin);
241 }
242 }
243 else
244 {
245 G4Exception("G4MuElecElasticModel::ComputeCrossSectionPerVolume","em0002",FatalException,"Model not applicable to particle type.");
246 }
247 }
248
249 if (verboseLevel > 3)
250 {
251 G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
252 G4cout << " - Cross section per Si atom (cm^2)=" << sigma/cm/cm << G4endl;
253 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density/(1./cm) << G4endl;
254 }
255
256 }
257
258 return sigma*density;
259}
@ FatalException
double G4double
Definition: G4Types.hh:64
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:232
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:208
virtual G4double FindValue(G4double e, G4int componentId=0) const
const G4String & GetParticleName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ GetKillBelowThreshold()

G4double G4MuElecElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 76 of file G4MuElecElasticModel.hh.

76{ return killBelowEnergy; }

◆ Initialise()

void G4MuElecElasticModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 101 of file G4MuElecElasticModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 263 of file G4MuElecElasticModel.cc.

268{
269
270 if (verboseLevel > 3)
271 G4cout << "Calling SampleSecondaries() of G4MuElecElasticModel" << G4endl;
272
273 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
274
275 if (electronEnergy0 < killBelowEnergy)
276 {
279 return ;
280 }
281
282 if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
283 {
284 G4double cosTheta = RandomizeCosTheta(electronEnergy0);
285
286 G4double phi = 2. * pi * G4UniformRand();
287
288 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
289 G4ThreeVector xVers = zVers.orthogonal();
290 G4ThreeVector yVers = zVers.cross(xVers);
291
292 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
293 G4double yDir = xDir;
294 xDir *= std::cos(phi);
295 yDir *= std::sin(phi);
296
297 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
298
300
302 }
303
304}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4double pi

◆ SetKillBelowThreshold()

void G4MuElecElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 138 of file G4MuElecElasticModel.hh.

139{
140 killBelowEnergy = threshold;
141
142 if (threshold < 5*CLHEP::eV)
143 {
144 G4Exception ("*** WARNING : the G4MuElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
145 threshold = 0.025*CLHEP::eV;
146 }
147
148}
@ JustWarning

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MuElecElasticModel::fParticleChangeForGamma
protected

Definition at line 80 of file G4MuElecElasticModel.hh.

Referenced by G4MuElecElasticModel(), Initialise(), and SampleSecondaries().


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