Geant4 10.7.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4ICRU73QOModel
33//
34// Author: Alexander Bagulya
35//
36// Creation date: 21.05.2010
37//
38// Modifications:
39//
40//
41// Class Description:
42//
43// Quantum Harmonic Oscillator Model for energy loss using atomic shell
44// structure of atoms taking into account Q^2 (main for projectile charge Q),
45// Q^3 and Q^4 terms for computation of energy loss due to binary collisions.
46// Can be applied on heavy negatively charged particles for the energy interval
47// 10 keV - 10 MeV scaled to the proton mass.
48//
49// Used data and formula of
50// 1. G4QAOLowEnergyLoss class, S.Chauvie, P.Nieminen, M.G.Pia. IEEE Trans.
51// Nucl. Sci. 54 (2007) 578.
52// 2. ShellStrength and ShellEnergy from ICRU'73 Report 2005,
53// 3. Data for Ta (Z=73) from P.Sigmund, A.Shinner. Eur. Phys. J. D15 (2001)
54// 165-172
55//
56// -------------------------------------------------------------------
57//
58
59#ifndef G4ICRU73QOModel_h
60#define G4ICRU73QOModel_h 1
61
63
64#include "G4VEmModel.hh"
65#include "G4AtomicShells.hh"
67
69
71{
72
73public:
74
75 explicit G4ICRU73QOModel(const G4ParticleDefinition* p = 0,
76 const G4String& nam = "ICRU73QO");
77
78 ~G4ICRU73QOModel() = default;
79
80 virtual void Initialise(const G4ParticleDefinition*,
81 const G4DataVector&) override;
82
85 G4double kineticEnergy,
86 G4double cutEnergy,
87 G4double maxEnergy) final;
88
91 G4double kineticEnergy,
92 G4double Z, G4double A,
93 G4double cutEnergy,
94 G4double maxEnergy) override;
95
98 G4double kineticEnergy,
99 G4double cutEnergy,
100 G4double maxEnergy) override;
101
104 G4double kineticEnergy,
105 G4double) override;
106
107 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109 const G4DynamicParticle*,
110 G4double tmin,
111 G4double maxEnergy) override;
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) override;
119
120protected:
121
123 G4double kinEnergy) final;
124
125private:
126
127 inline void SetParticle(const G4ParticleDefinition* p);
128 inline void SetLowestKinEnergy(G4double val);
129
130 G4double DEDX(const G4Material* material, G4double kineticEnergy);
131
132 G4double DEDXPerElement(G4int Z, G4double kineticEnergy);
133
134 // get number of shell, energy and oscillator strengths for material
135 G4int GetNumberOfShells(G4int Z) const;
136
137 G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const;
138 G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const;
139 G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const;
140
141 // calculate stopping number for L's term
142 G4double GetL0(G4double normEnergy) const;
143 // terms in Z^2
144 G4double GetL1(G4double normEnergy) const;
145 // terms in Z^3
146 G4double GetL2(G4double normEnergy) const;
147 // terms in Z^4
148
149 // hide assignment operator
150 G4ICRU73QOModel & operator=(const G4ICRU73QOModel &right) = delete;
151 G4ICRU73QOModel(const G4ICRU73QOModel&) = delete;
152
153 const G4ParticleDefinition* particle;
154 G4ParticleDefinition* theElectron;
155 G4ParticleChangeForLoss* fParticleChange;
156 G4DensityEffectData* denEffData;
157
158 G4double mass;
159 G4double charge;
160 G4double chargeSquare;
161 G4double massRate;
162 G4double ratio;
163 G4double lowestKinEnergy;
164
165 G4bool isInitialised;
166
167 // Z of element at now avaliable for the model
168 static const G4int NQOELEM = 26;
169 static const G4int NQODATA = 130;
170 static const G4int ZElementAvailable[NQOELEM];
171
172 // number, energy and oscillator strengths
173 // for an harmonic oscillator model of material
174 static const G4int startElemIndex[NQOELEM];
175 static const G4int nbofShellsForElement[NQOELEM];
176 static const G4double ShellEnergy[NQODATA];
177 static const G4double SubShellOccupation[NQODATA]; // Z * ShellStrength
178
179 G4int indexZ[100];
180
181 // variable for calculation of stopping number of L's term
182 static const G4double L0[67][2];
183 static const G4double L1[22][2];
184 static const G4double L2[14][2];
185
186 G4int sizeL0;
187 G4int sizeL1;
188 G4int sizeL2;
189
190 static const G4double factorBethe[99];
191
192};
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195
196inline void G4ICRU73QOModel::SetParticle(const G4ParticleDefinition* p)
197{
198 particle = p;
199 mass = particle->GetPDGMass();
200 charge = particle->GetPDGCharge()/CLHEP::eplus;
201 chargeSquare = charge*charge;
202 massRate = mass/CLHEP::proton_mass_c2;
203 ratio = CLHEP::electron_mass_c2/mass;
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
207
208inline void G4ICRU73QOModel::SetLowestKinEnergy(G4double val)
209{
210 lowestKinEnergy = val;
211}
212
213#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length) override
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
~G4ICRU73QOModel()=default
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
G4double GetPDGCharge() const