Geant4-11
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 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
81
84 G4double kineticEnergy,
85 G4double cutEnergy,
86 G4double maxEnergy);
87
90 G4double kineticEnergy,
92 G4double cutEnergy,
93 G4double maxEnergy) override;
94
97 G4double kineticEnergy,
98 G4double cutEnergy,
99 G4double maxEnergy) override;
100
103 G4double kineticEnergy,
104 G4double) override;
105
106 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
108 const G4DynamicParticle*,
109 G4double tmin,
110 G4double maxEnergy) override;
111
112 // add correction to energy loss and compute non-ionizing energy loss
114 const G4DynamicParticle*,
115 const G4double& length,
116 G4double& eloss) override;
117
118 // hide assignment operator
119 G4ICRU73QOModel & operator=(const G4ICRU73QOModel &right) = delete;
121
122protected:
123
125 G4double kinEnergy) final;
126
127private:
128
129 inline void SetParticle(const G4ParticleDefinition* p);
130 inline void SetLowestKinEnergy(G4double val);
131
132 G4double DEDX(const G4Material* material, G4double kineticEnergy);
133
134 G4double DEDXPerElement(G4int Z, G4double kineticEnergy);
135
136 // get number of shell, energy and oscillator strengths for material
138
139 G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const;
140 G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const;
141 G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const;
142
143 // calculate stopping number for L's term
144 G4double GetL0(G4double normEnergy) const;
145 // terms in Z^2
146 G4double GetL1(G4double normEnergy) const;
147 // terms in Z^3
148 G4double GetL2(G4double normEnergy) const;
149 // terms in Z^4
150
155
162
164
165 // Z of element at now avaliable for the model
166 static const G4int NQOELEM = 26;
167 static const G4int NQODATA = 130;
169
170 // number, energy and oscillator strengths
171 // for an harmonic oscillator model of material
175 static const G4double SubShellOccupation[NQODATA]; // Z * ShellStrength
176
178
179 // variable for calculation of stopping number of L's term
180 static const G4double L0[67][2];
181 static const G4double L1[22][2];
182 static const G4double L2[14][2];
183
187
188 static const G4double factorBethe[99];
189
190};
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193
195{
196 particle = p;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
207{
208 lowestKinEnergy = val;
209}
210
211#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
static const G4double L1[22][2]
G4double lowestKinEnergy
static const G4int ZElementAvailable[NQOELEM]
static const G4double L2[14][2]
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const
static const G4double ShellEnergy[NQODATA]
static const G4int nbofShellsForElement[NQOELEM]
G4ParticleDefinition * theElectron
static const G4double L0[67][2]
const G4ParticleDefinition * particle
G4double GetL1(G4double normEnergy) const
G4double DEDXPerElement(G4int Z, G4double kineticEnergy)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double GetL2(G4double normEnergy) const
G4ParticleChangeForLoss * fParticleChange
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4int GetNumberOfShells(G4int Z) const
static const G4double factorBethe[99]
G4ICRU73QOModel & operator=(const G4ICRU73QOModel &right)=delete
G4DensityEffectData * denEffData
void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss) override
G4ICRU73QOModel(const G4ICRU73QOModel &)=delete
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static const G4int startElemIndex[NQOELEM]
void SetLowestKinEnergy(G4double val)
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
static const G4int NQOELEM
G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const
G4double GetL0(G4double normEnergy) const
G4double DEDX(const G4Material *material, G4double kineticEnergy)
~G4ICRU73QOModel()=default
static const G4double SubShellOccupation[NQODATA]
void SetParticle(const G4ParticleDefinition *p)
G4ICRU73QOModel(const G4ParticleDefinition *p=0, const G4String &nam="ICRU73QO")
G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
static const G4int NQODATA
G4double GetPDGCharge() const
static constexpr double eplus
static constexpr double electron_mass_c2
static constexpr double proton_mass_c2
string material
Definition: eplot.py:19