Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ICRU73QOModel.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4ICRU73QOModel
34//
35// Author: Alexander Bagulya
36//
37// Creation date: 21.05.2010
38//
39// Modifications:
40//
41//
42// Class Description:
43//
44// Quantum Harmonic Oscillator Model for energy loss using atomic shell
45// structure of atoms taking into account Q^2 (main for projectile charge Q),
46// Q^3 and Q^4 terms for computation of energy loss due to binary collisions.
47// Can be applied on heavy negatively charged particles for the energy interval
48// 10 keV - 10 MeV scaled to the proton mass.
49//
50// Used data and formula of
51// 1. G4QAOLowEnergyLoss class, S.Chauvie, P.Nieminen, M.G.Pia. IEEE Trans.
52// Nucl. Sci. 54 (2007) 578.
53// 2. ShellStrength and ShellEnergy from ICRU'73 Report 2005,
54// 3. Data for Ta (Z=73) from P.Sigmund, A.Shinner. Eur. Phys. J. D15 (2001)
55// 165-172
56//
57// -------------------------------------------------------------------
58//
59
60#ifndef G4ICRU73QOModel_h
61#define G4ICRU73QOModel_h 1
62
64
65#include "G4VEmModel.hh"
66#include "G4AtomicShells.hh"
68
70
72{
73
74public:
75
77 const G4String& nam = "ICRU73QO");
78
79 virtual ~G4ICRU73QOModel();
80
81 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
82
85 G4double kineticEnergy,
86 G4double cutEnergy,
87 G4double maxEnergy);
88
91 G4double kineticEnergy,
92 G4double Z, G4double A,
93 G4double cutEnergy,
94 G4double maxEnergy);
95
98 G4double kineticEnergy,
99 G4double cutEnergy,
100 G4double maxEnergy);
101
104 G4double kineticEnergy,
105 G4double);
106
107 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109 const G4DynamicParticle*,
110 G4double tmin,
111 G4double maxEnergy);
112
113 // add correction to energy loss and compute non-ionizing energy loss
114 virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
115 const G4DynamicParticle*,
116 G4double& eloss,
117 G4double& niel,
118 G4double length);
119
120protected:
121
123 G4double kinEnergy);
124
125private:
126
127 inline void SetParticle(const G4ParticleDefinition* p);
128 inline void SetLowestKinEnergy(const G4double val);
129
130 G4double DEDX(const G4Material* material, G4double kineticEnergy);
131
132 G4double DEDXPerElement(G4int Z, G4double kineticEnergy);
133
134 // hide assignment operator
135 G4ICRU73QOModel & operator=(const G4ICRU73QOModel &right);
137
138 const G4ParticleDefinition* particle;
139 G4ParticleDefinition* theElectron;
140 G4ParticleChangeForLoss* fParticleChange;
141 G4DensityEffectData* denEffData;
142
143 G4double mass;
144 G4double charge;
145 G4double chargeSquare;
146 G4double massRate;
147 G4double ratio;
148 G4double lowestKinEnergy;
149
150 G4bool isInitialised;
151
152 // get number of shell, energy and oscillator strenghts for material
153 G4int GetNumberOfShells(G4int Z) const;
154
155 G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const;
156 G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const;
157 G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const;
158
159 // calculate stopping number for L's term
160 G4double GetL0(G4double normEnergy) const;
161 // terms in Z^2
162 G4double GetL1(G4double normEnergy) const;
163 // terms in Z^3
164 G4double GetL2(G4double normEnergy) const;
165 // terms in Z^4
166
167
168 // Z of element at now avaliable for the model
169 static const G4int NQOELEM = 26;
170 static const G4int NQODATA = 130;
171 static const G4int ZElementAvailable[NQOELEM];
172
173 // number, energy and oscillator strenghts
174 // for an harmonic oscillator model of material
175 static const G4int startElemIndex[NQOELEM];
176 static const G4int nbofShellsForElement[NQOELEM];
177 static const G4double ShellEnergy[NQODATA];
178 static const G4double SubShellOccupation[NQODATA]; // Z * ShellStrength
179
180 G4int indexZ[100];
181
182 // variable for calculation of stopping number of L's term
183 static const G4double L0[67][2];
184 static const G4double L1[22][2];
185 static const G4double L2[14][2];
186
187 G4int sizeL0;
188 G4int sizeL1;
189 G4int sizeL2;
190
191 static const G4double factorBethe[99];
192
193};
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196
197inline void G4ICRU73QOModel::SetParticle(const G4ParticleDefinition* p)
198{
199 particle = p;
200 mass = particle->GetPDGMass();
201 charge = particle->GetPDGCharge()/CLHEP::eplus;
202 chargeSquare = charge*charge;
203 massRate = mass/CLHEP::proton_mass_c2;
204 ratio = CLHEP::electron_mass_c2/mass;
205}
206
207//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
208
209inline G4int G4ICRU73QOModel::GetNumberOfShells(G4int Z) const
210{
211 G4int nShell = 0;
212
213 if(indexZ[Z] >= 0) { nShell = nbofShellsForElement[indexZ[Z]];
214 } else { nShell = G4AtomicShells::GetNumberOfShells(Z); }
215
216 return nShell;
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220
221inline G4double
222G4ICRU73QOModel::GetShellEnergy(G4int Z, G4int nbOfTheShell) const
223{
224 G4double shellEnergy = 0.;
225
226 G4int idx = indexZ[Z];
227
228 if(idx >= 0) { shellEnergy = ShellEnergy[startElemIndex[idx] + nbOfTheShell]*CLHEP::eV;
229 } else { shellEnergy = GetOscillatorEnergy(Z, nbOfTheShell); }
230
231 return shellEnergy;
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235
236inline G4double
237G4ICRU73QOModel::GetShellStrength(G4int Z, G4int nbOfTheShell) const
238{
239 G4double shellStrength = 0.;
240
241 G4int idx = indexZ[Z];
242
243 if(idx >= 0) { shellStrength = SubShellOccupation[startElemIndex[idx] + nbOfTheShell] / Z;
244 } else { shellStrength = G4double(G4AtomicShells::GetNumberOfElectrons(Z,nbOfTheShell))/Z; }
245
246 return shellStrength;
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
250
251inline void G4ICRU73QOModel::SetLowestKinEnergy(const G4double val)
252{
253 lowestKinEnergy = val;
254}
255
256#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
static G4int GetNumberOfShells(G4int Z)
virtual ~G4ICRU73QOModel()
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
G4double GetPDGCharge() const