Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4AtimaEnergyLossModel.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 //
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 
54 class G4EmCorrections;
56 
58 {
59 
60 public:
61 
62  explicit G4AtimaEnergyLossModel(const G4ParticleDefinition* p = nullptr,
63  const G4String& nam = "Atima");
64 
65  virtual ~G4AtimaEnergyLossModel();
66 
67  virtual void Initialise(const G4ParticleDefinition*,
68  const G4DataVector&) override;
69 
71  const G4MaterialCutsCouple* couple) override;
72 
74  const G4ParticleDefinition*,
76  G4double cutEnergy,
77  G4double maxEnergy);
78 
80  const G4ParticleDefinition*,
81  G4double kineticEnergy,
83  G4double cutEnergy,
84  G4double maxEnergy) override;
85 
87  const G4ParticleDefinition*,
88  G4double kineticEnergy,
89  G4double cutEnergy,
90  G4double maxEnergy) override;
91 
93  const G4ParticleDefinition*,
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*>*,
112  const G4MaterialCutsCouple*,
113  const G4DynamicParticle*,
114  G4double tmin,
115  G4double maxEnergy) override;
116 
117 protected:
118 
120  G4double kinEnergy) override;
121 
122  inline G4double GetChargeSquareRatio() const;
123 
124  inline void SetChargeSquareRatio(G4double val);
125 
126 private:
127 
128  void SetupParameters();
129 
130  inline void SetParticle(const G4ParticleDefinition* p);
131 
132  inline void SetGenericIon(const G4ParticleDefinition* p);
133 
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);
140 
141  // hide assignment operator
144 
151 
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];
175 
176  static const G4double atima_vfermi[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 
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 
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
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
static const G4double x0[92]
G4double sezi_dedx_e(const G4double zp, const G4double ep, const G4double at, const G4double zt)
static const G4double del_0[92]
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
G4double GetPDGCharge() const
G4double StoppingPower(G4double ap, G4double zp, G4double ep, G4double at, G4double zt)
static const G4double atima_lambda_screening[92]
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
G4AtimaEnergyLossModel & operator=(const G4AtimaEnergyLossModel &right)=delete
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double sezi_p_se(const G4double energy, const G4double at, const G4double zt)
Float_t Z
G4ParticleChangeForLoss * fParticleChange
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static const G4double element_atomic_weights[110]
static const G4double x1[92]
static const G4double proton_stopping_coef[92][8]
G4double dedx_n(const G4double ap, const G4double zp, const G4double ep, const G4double at, const G4double zt)
double energy
Definition: plottest35.C:25
void SetChargeSquareRatio(G4double val)
double A(double temperature)
static const G4double atima_vfermi[92]
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double GetChargeSquareRatio() const
Float_t mat
static const G4double ionisation_potentials_z[121]
static const G4double m0[92]
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double EnergyTable_interpolate(G4double xval, const G4double *y)
static const G4double ls_coefficients_a[110][200]
void SetGenericIon(const G4ParticleDefinition *p)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &, G4double &, G4double) override
Definition: G4Pow.hh:56
G4ParticleDefinition * theElectron
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4AtimaEnergyLossModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Atima")
const G4ParticleDefinition * particle
G4double Bethek_dedx_e(G4double ap, G4double zp, G4double ep, G4double at, G4double zt)
void SetParticle(const G4ParticleDefinition *p)
static const G4double afermi[92]
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
static const G4double ls_coefficients_ahi[110][200]
static constexpr double eplus
static const G4double c[92]