Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmCorrections.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: G4EmCorrections.hh 108386 2018-02-09 15:38:32Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4EmCorrections
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 13.01.2005
38 //
39 // Modifications:
40 // 28.04.2006 General cleanup, add finite size corrections (V.Ivanchenko)
41 // 13.05.2006 Add corrections for ion stopping (V.Ivanhcenko)
42 // 20.05.2008 Removed Finite Size correction (V.Ivanchenko)
43 // 12.09.2008 Added inlined interfaces to effective charge (V.Ivanchenko)
44 // 19.04.2012 Fix reproducibility problem (A.Ribon)
45 //
46 // Class Description:
47 //
48 // This class provides calculation of EM corrections to ionisation
49 //
50 
51 // -------------------------------------------------------------------
52 //
53 
54 #ifndef G4EmCorrections_h
55 #define G4EmCorrections_h 1
56 
58 
59 #include "globals.hh"
60 #include "G4ionEffectiveCharge.hh"
61 #include "G4Material.hh"
62 #include "G4ParticleDefinition.hh"
63 
64 class G4VEmModel;
65 class G4PhysicsVector;
66 class G4IonTable;
69 class G4Pow;
70 
72 {
73 
74 public:
75 
76  explicit G4EmCorrections(G4int verb);
77 
78  virtual ~G4EmCorrections();
79 
81  const G4Material*,
83  G4double cutEnergy);
84 
86  const G4MaterialCutsCouple*,
87  G4double kineticEnergy);
88 
90  const G4Material*,
91  G4double kineticEnergy);
92 
94  const G4Material*,
95  G4double kineticEnergy);
96 
98  const G4Material*,
99  G4double kineticEnergy);
100 
102  const G4Material*,
103  G4double kineticEnergy);
104 
106  const G4Material*,
107  G4double kineticEnergy);
108 
110  const G4Material*,
111  G4double kineticEnergy);
112 
114  const G4Material*,
115  G4double kineticEnergy);
116 
118  const G4Material*,
119  G4double kineticEnergy);
120 
122  const G4Material*,
123  G4double kineticEnergy);
124 
126  const G4Material*,
127  G4double kineticEnergy);
128 
130  const G4Material*,
131  G4double kineticEnergy);
132 
134  const G4Material*,
135  G4double kineticEnergy);
136 
138  G4PhysicsVector* dVector);
139 
140  void InitialiseForNewRun();
141 
142  // effective charge correction using stopping power data
144  const G4Material*,
145  G4double kineticEnergy);
146 
147  // effective charge of an ion
149  const G4Material*,
150  G4double kineticEnergy);
151 
152  inline
154  const G4Material*,
155  G4double kineticEnergy);
156 
157  // ionisation models for ions
158  inline void SetIonisationModels(G4VEmModel* m1 = nullptr,
159  G4VEmModel* m2 = nullptr);
160 
161  inline G4int GetNumberOfStoppingVectors() const;
162 
163  inline void SetVerbose(G4int verb);
164 
165 private:
166 
167  void Initialise();
168 
169  void BuildCorrectionVector();
170 
172  const G4Material*,
173  G4double kineticEnergy);
174 
175  G4double KShell(G4double theta, G4double eta);
176 
177  G4double LShell(G4double theta, G4double eta);
178 
179  G4int Index(G4double x, const G4double* y, G4int n) const;
180 
182  G4double y1, G4double y2) const;
183 
185  G4double y1, G4double y2, G4double z11, G4double z21,
186  G4double z12, G4double z22) const;
187 
188  // hide assignment operator
189  G4EmCorrections & operator=(const G4EmCorrections &right) = delete;
190  G4EmCorrections(const G4EmCorrections&) = delete;
191 
193 
194  static const G4double ZD[11];
195  static const G4double UK[20];
196  static const G4double VK[20];
197  static G4double ZK[20];
198  static const G4double Eta[29];
199  static G4double CK[20][29];
200  static G4double CL[26][28];
201  static const G4double UL[26];
202  static G4double VL[26];
203 
207 
209 
210  std::vector<const G4Material*> currmat;
211  std::map< G4int, std::vector<G4double> > thcorr;
212  size_t ncouples;
213 
220 
222 
226 
242 
244 
249 
252 
253  // Ion stopping data
257  std::vector<G4int> Zion;
258  std::vector<G4int> Aion;
259  std::vector<G4String> materialName;
260 
261  std::vector<const G4ParticleDefinition*> ionList;
262 
263  std::vector<const G4Material*> materialList;
264  std::vector<G4PhysicsVector*> stopData;
265 
268 };
269 
270 inline G4int
272 {
273  G4int iddd = n-1;
274  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
275  do {--iddd;} while (iddd>0 && x<y[iddd]);
276  return iddd;
277 }
278 
280  G4double y1, G4double y2) const
281 {
282  return y1 + (y2 - y1)*(xv - x1)/(x2 - x1);
283 }
284 
288  G4double z11, G4double z21,
289  G4double z12, G4double z22) const
290 {
291  return (z11*(x2-xv)*(y2-yv) + z22*(xv-x1)*(yv-y1) +
292  0.5*(z12*((x2-xv)*(yv-y1)+(xv-x1)*(y2-yv))+
293  z21*((xv-x1)*(y2-yv)+(yv-y1)*(x2-xv))))
294  / ((x2-x1)*(y2-y1));
295 }
296 
297 inline void
299 {
300  if(mod1) { ionLEModel = mod1; }
301  if(mod2) { ionHEModel = mod2; }
302 }
303 
305 {
306  return nIons;
307 }
308 
309 inline G4double
311  const G4Material* mat,
313 {
314  return effCharge.EffectiveCharge(p,mat,kineticEnergy);
315 }
316 
317 inline G4double
319  const G4Material* mat,
321 {
322  return effCharge.EffectiveChargeSquareRatio(p,mat,kineticEnergy);
323 }
324 
326 {
327  verbose = verb;
328 }
329 
330 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
331 
332 #endif
static const G4double UL[26]
Float_t x
Definition: compare.C:6
static G4double CL[26][28]
const G4ParticleDefinition * particle
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
static G4double CK[20][29]
std::map< G4int, std::vector< G4double > > thcorr
const G4ElementVector * theElementVector
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double MottCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Float_t y1[n_points_granero]
Definition: compare.C:5
std::vector< const G4Material * > materialList
static const G4double ZD[11]
static constexpr double m2
Definition: G4SIunits.hh:130
Float_t x1[n_points_granero]
Definition: compare.C:5
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Float_t y
Definition: compare.C:6
G4double BlochCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const char * p
Definition: xmltok.h:285
static G4LPhysicsFreeVector * BarkasCorr
G4double SpinCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4VEmModel * ionHEModel
std::vector< G4String > materialName
const G4Material * material
G4double Bethe(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4IonTable * ionTable
std::vector< G4PhysicsVector * > stopData
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double LShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double ComputeIonCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static const G4double UK[20]
G4double Value(G4double xv, G4double x1, G4double x2, G4double y1, G4double y2) const
Float_t y2[n_points_geant4]
Definition: compare.C:26
Float_t Z
std::vector< G4int > Zion
G4int GetNumberOfStoppingVectors() const
void SetIonisationModels(G4VEmModel *m1=nullptr, G4VEmModel *m2=nullptr)
G4double KShell(G4double theta, G4double eta)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void BuildCorrectionVector()
const G4double * atomDensity
static const G4double VK[20]
double A(double temperature)
G4double IonHighOrderCorrections(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
static G4LPhysicsFreeVector * ThetaK
static G4LPhysicsFreeVector * ThetaL
Float_t mat
G4double KShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
void AddStoppingData(G4int Z, G4int A, const G4String &materialName, G4PhysicsVector *dVector)
G4EmCorrections(G4int verb)
const G4Material * curMaterial
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
G4EmCorrections & operator=(const G4EmCorrections &right)=delete
G4VEmModel * ionLEModel
const G4ParticleDefinition * curParticle
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double LShell(G4double theta, G4double eta)
static G4double VL[26]
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
G4ionEffectiveCharge effCharge
void SetVerbose(G4int verb)
G4int Index(G4double x, const G4double *y, G4int n) const
Definition: G4Pow.hh:56
std::vector< G4int > Aion
G4double ShellCorrectionSTD(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static G4double ZK[20]
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
G4double Value2(G4double xv, G4double yv, G4double x1, G4double x2, G4double y1, G4double y2, G4double z11, G4double z21, G4double z12, G4double z22) const
G4double DensityCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4PhysicsVector * curVector
Char_t n[5]
static const G4double Eta[29]
std::vector< const G4Material * > currmat
void SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Float_t x2[n_points_geant4]
Definition: compare.C:26
std::vector< const G4ParticleDefinition * > ionList
virtual ~G4EmCorrections()