Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ICRU73QOModel.hh
이 파일의 문서화 페이지로 가기
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: G4ICRU73QOModel.hh 108737 2018-03-02 13:49:56Z gcosmo $
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"
67 #include "G4DensityEffectData.hh"
68 
70 
72 {
73 
74 public:
75 
76  explicit G4ICRU73QOModel(const G4ParticleDefinition* p = 0,
77  const G4String& nam = "ICRU73QO");
78 
79  ~G4ICRU73QOModel() = default;
80 
81  virtual void Initialise(const G4ParticleDefinition*,
82  const G4DataVector&) override;
83 
85  const G4ParticleDefinition*,
87  G4double cutEnergy,
88  G4double maxEnergy) final;
89 
91  const G4ParticleDefinition*,
92  G4double kineticEnergy,
94  G4double cutEnergy,
95  G4double maxEnergy) override;
96 
98  const G4ParticleDefinition*,
99  G4double kineticEnergy,
100  G4double cutEnergy,
101  G4double maxEnergy) override;
102 
104  const G4ParticleDefinition*,
105  G4double kineticEnergy,
106  G4double) override;
107 
108  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109  const G4MaterialCutsCouple*,
110  const G4DynamicParticle*,
111  G4double tmin,
112  G4double maxEnergy) override;
113 
114  // add correction to energy loss and compute non-ionizing energy loss
115  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
116  const G4DynamicParticle*,
117  G4double& eloss,
118  G4double& niel,
119  G4double length) override;
120 
121 protected:
122 
124  G4double kinEnergy) final;
125 
126 private:
127 
128  inline void SetParticle(const G4ParticleDefinition* p);
129  inline void SetLowestKinEnergy(G4double val);
130 
131  G4double DEDX(const G4Material* material, G4double kineticEnergy);
132 
133  G4double DEDXPerElement(G4int Z, G4double kineticEnergy);
134 
135  // get number of shell, energy and oscillator strenghts for material
136  G4int GetNumberOfShells(G4int Z) const;
137 
138  G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const;
139  G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const;
140  G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const;
141 
142  // calculate stopping number for L's term
143  G4double GetL0(G4double normEnergy) const;
144  // terms in Z^2
145  G4double GetL1(G4double normEnergy) const;
146  // terms in Z^3
147  G4double GetL2(G4double normEnergy) const;
148  // terms in Z^4
149 
150  // hide assignment operator
151  G4ICRU73QOModel & operator=(const G4ICRU73QOModel &right) = delete;
152  G4ICRU73QOModel(const G4ICRU73QOModel&) = delete;
153 
158 
165 
167 
168  // Z of element at now avaliable for the model
169  static const G4int NQOELEM = 26;
170  static const G4int NQODATA = 130;
172 
173  // number, energy and oscillator strenghts
174  // for an harmonic oscillator model of material
175  static const G4int startElemIndex[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 
190 
191  static const G4double factorBethe[99];
192 
193 };
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
198 {
199  particle = p;
200  mass = particle->GetPDGMass();
205 }
206 
207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
208 
210 {
211  lowestKinEnergy = val;
212 }
213 
214 #endif
G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const
G4double DEDXPerElement(G4int Z, G4double kineticEnergy)
G4double lowestKinEnergy
G4ICRU73QOModel(const G4ParticleDefinition *p=0, const G4String &nam="ICRU73QO")
G4int GetNumberOfShells(G4int Z) const
G4double GetL2(G4double normEnergy) const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
static const G4double ShellEnergy[NQODATA]
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static const G4double L1[22][2]
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
const char * p
Definition: xmltok.h:285
void SetLowestKinEnergy(G4double val)
G4double GetPDGCharge() const
~G4ICRU73QOModel()=default
static const G4int ZElementAvailable[NQOELEM]
static const G4int NQOELEM
G4double GetPDGMass() const
static const G4int startElemIndex[NQOELEM]
G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const
G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const
static constexpr double proton_mass_c2
Float_t Z
static const G4double SubShellOccupation[NQODATA]
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ICRU73QOModel & operator=(const G4ICRU73QOModel &right)=delete
static const G4double L2[14][2]
static constexpr double electron_mass_c2
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
double A(double temperature)
G4double DEDX(const G4Material *material, G4double kineticEnergy)
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length) override
G4DensityEffectData * denEffData
static const G4double factorBethe[99]
const G4ParticleDefinition * particle
static const G4double L0[67][2]
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
int G4int
Definition: G4Types.hh:78
void SetParticle(const G4ParticleDefinition *p)
static const G4int nbofShellsForElement[NQOELEM]
static const G4int NQODATA
G4double GetL1(G4double normEnergy) const
G4double GetL0(G4double normEnergy) const
G4ParticleChangeForLoss * fParticleChange
G4ParticleDefinition * theElectron
static constexpr double eplus
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override