Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AtimaEnergyLossModel.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// GEANT4 Class header file
28//
29// File name: G4AtimaEnergyLossModel
30//
31// Author: Jose Luis Rodriguez Sanchez on base of ATIMA code
32//
33// Creation date: 16.01.2018
34//
35// Modifications:
36//
37//
38// Class Description:
39//
40// Implementation of ATIMA model of energy loss
41// by heavy charged particles
42
43// -------------------------------------------------------------------
44//
45
46#ifndef G4AtimaEnergyLossModel_h
47#define G4AtimaEnergyLossModel_h 1
48
50
51#include "G4VEmModel.hh"
52#include "G4NistManager.hh"
53
54class G4EmCorrections;
56
58{
59
60public:
61
62 explicit G4AtimaEnergyLossModel(const G4ParticleDefinition* p = nullptr,
63 const G4String& nam = "Atima");
64
66
67 virtual void Initialise(const G4ParticleDefinition*,
68 const G4DataVector&) override;
69
71 const G4MaterialCutsCouple* couple) override;
72
75 G4double kineticEnergy,
76 G4double cutEnergy,
77 G4double maxEnergy);
78
81 G4double kineticEnergy,
83 G4double cutEnergy,
84 G4double maxEnergy) override;
85
88 G4double kineticEnergy,
89 G4double cutEnergy,
90 G4double maxEnergy) override;
91
94 G4double kineticEnergy,
95 G4double) override;
96
98 const G4Material* mat,
99 G4double kineticEnergy) override;
100
102 const G4Material* mat,
103 G4double kineticEnergy) override;
104
105 virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
106 const G4DynamicParticle*,
107 G4double&,
108 G4double&,
109 G4double) override;
110
111 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
113 const G4DynamicParticle*,
114 G4double tmin,
115 G4double maxEnergy) override;
116
117protected:
118
120 G4double kinEnergy) override;
121
122 inline G4double GetChargeSquareRatio() const;
123
124 inline void SetChargeSquareRatio(G4double val);
125
126private:
127
128 void SetupParameters();
129
130 inline void SetParticle(const G4ParticleDefinition* p);
131
132 inline void SetGenericIon(const G4ParticleDefinition* p);
133
134 G4double StoppingPower(G4double ap, G4double zp, G4double ep, G4double at, G4double zt);
135 G4double Bethek_dedx_e(G4double ap,G4double zp,G4double ep,G4double at,G4double zt);
136 G4double dedx_n(const G4double ap, const G4double zp, const G4double ep, const G4double at, const G4double zt);
137 G4double sezi_dedx_e(const G4double zp, const G4double ep, const G4double at, const G4double zt);
138 G4double sezi_p_se(const G4double energy, const G4double at, const G4double zt);
139 G4double EnergyTable_interpolate(G4double xval, const G4double* y);
140
141 // hide assignment operator
142 G4AtimaEnergyLossModel & operator=(const G4AtimaEnergyLossModel &right) = delete;
144
145 const G4ParticleDefinition* particle;
146 G4ParticleDefinition* theElectron;
147 G4EmCorrections* corr;
148 G4ParticleChangeForLoss* fParticleChange;
149 G4NistManager* nist;
150 G4Pow* g4calc;
151
152 G4double mass;
153 G4double tlimit;
154 G4double spin;
155 G4double magMoment2;
156 G4double chargeSquare;
157 G4double ratio;
158 G4double formfact;
159 G4double corrFactor;
160 G4bool isIon;
161 G4double MLN10;
162 G4double atomic_mass_unit;
163 G4double dedx_constant;
164 G4double electron_mass;
165 G4double fine_structure;
166 G4double domega2dx_constant;
167
168 static G4double stepE;
169 static G4double tableE[200];
170 static const G4double element_atomic_weights[110];
171 static const G4double ls_coefficients_a[110][200];
172 static const G4double ls_coefficients_ahi[110][200];
173 static const G4double proton_stopping_coef[92][8];
174 static const G4double ionisation_potentials_z[121];
175
176 static const G4double atima_vfermi[92];
177 static const G4double atima_lambda_screening[92];
178 static const G4double x0[92];
179 static const G4double x1[92];
180 static const G4double afermi[92];
181 static const G4double c[92];
182 static const G4double m0[92];
183 static const G4double del_0[92];
184
185};
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
189inline void G4AtimaEnergyLossModel::SetParticle(const G4ParticleDefinition* p)
190{
191 if(particle != p) {
192 particle = p;
193 if(p->GetBaryonNumber() > 3 || p->GetPDGCharge() > CLHEP::eplus)
194 { isIon = true; }
195 SetupParameters();
196 }
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200
201inline void G4AtimaEnergyLossModel::SetGenericIon(const G4ParticleDefinition* p)
202{
203 if(p && p->GetParticleName() == "GenericIon") { isIon = true; }
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207
209{
210 return chargeSquare;
211}
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214
216{
217 chargeSquare = val;
218}
219
220//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
221
222#endif
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &, G4double &, G4double) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
void SetChargeSquareRatio(G4double val)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double GetChargeSquareRatio() const
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double GetPDGCharge() const
const G4String & GetParticleName() const
Definition: G4Pow.hh:49