Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ScreeningMottCrossSection.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 // G4ScreeningMottCrossSection.hh
27 //-------------------------------------------------------------------
28 //
29 // GEANT4 Class header file
30 //
31 // File name: G4ScreeningMottCrossSection
32 //
33 // Author: Cristina Consolandi
34 //
35 // Creation date: 20.10.2011
36 //
37 // Modifications:
38 // 27-05-2012 Added Analytic Fitting to the Mott Cross Section by means of G4MottCoefficients class.
39 //
40 //
41 // Class Description:
42 // Computation of electron Coulomb Scattering Cross Section.
43 // Suitable for high energy electrons and light target materials.
44 //
45 // Reference:
46 // M.J. Boschini et al.
47 // "Non Ionizing Energy Loss induced by Electrons in the Space Environment"
48 // Proc. of the 13th International Conference on Particle Physics and Advanced Technology
49 // (13th ICPPAT, Como 3-7/10/2011), World Scientific (Singapore).
50 // Available at: http://arxiv.org/abs/1111.4042v4
51 //
52 // 1) Mott Differential Cross Section Approximation:
53 // For Target material up to Z=92 (U):
54 // As described in http://arxiv.org/abs/1111.4042v4
55 // par. 2.1 , eq. (16)-(17)
56 // Else (Z>92):
57 // W. A. McKinley and H. Fashbach, Phys. Rev. 74, (1948) 1759.
58 // 2) Screening coefficient:
59 // vomn G. Moliere, Z. Naturforsh A2 (1947), 133-145; A3 (1948), 78.
60 // 3) Nuclear Form Factor:
61 // A.V. Butkevich et al. Nucl. Instr. and Meth. in Phys. Res. A 488 (2002), 282-294.
62 //
63 // ----------------------------------------------------------------------------------------
64 
65 //
66 #ifndef G4ScreeningMottCrossSection_h
67 #define G4ScreeningMottCrossSection_h 1
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
71 
72 #include "G4MottCoefficients.hh"
73 #include "globals.hh"
74 #include "G4Material.hh"
75 #include "G4Element.hh"
76 #include "G4ElementVector.hh"
77 #include "G4NistManager.hh"
78 #include "G4ThreeVector.hh"
79 #include "G4Pow.hh"
80 #include "G4LossTableManager.hh"
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
84 static const G4int DIM = 750;
85 
87 
89 {
90 
91 public:
92 
93  explicit G4ScreeningMottCrossSection();
94 
96 
97  void Initialise(const G4ParticleDefinition*, G4double cosThetaLim);
98 
101 
102  inline void SetupParticle(const G4ParticleDefinition*);
103  void SetupKinematic(G4double kinEnergy ,G4double Z);
104 
107 
108  inline G4double GetMom2CM()const;
109  inline G4double GetMom2Lab()const;
110  inline G4double GetTrec() const;
111  inline G4double GetScreeningCoefficient() const;
112  inline G4double GetTotalCross() const;
113 
122 
123 private:
124 
127 
131 
133  //cost - min - max
135  G4double cosThetaMax;// def -1.0
136 
139 
140  //energy cut
143 
144  // projectile........................
146 
149 
150  //lab of incedent particle
154 
155  //relative system with nucleus
162 
163  // target nucleus
170 
171  //constants
175 
176  //angle
178 };
179 
180 
181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182 
184 {
185  particle = p;
186  mass = particle->GetPDGMass();
187  spin = particle->GetPDGSpin();
188  if(0.0 != spin) { spin = 0.5; }
189  tkin = 0.0;
190 }
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192 
193 
195 {
196  return mom2;
197 }
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199 
200 
202 {
203  return momLab2;
204 }
205 
206 
207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208 
210 {
211  return Trec;
212 }
213 
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215 
217 {
218  return As;
219 }
220 
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
223 
225 {
226  return TotalCross;
227 }
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
230 
231 #endif
const char * p
Definition: xmltok.h:285
G4ScreeningMottCrossSection & operator=(const G4ScreeningMottCrossSection &right)=delete
G4double GetPDGMass() const
void SetupParticle(const G4ParticleDefinition *)
Float_t Z
G4double GetPDGSpin() const
double G4double
Definition: G4Types.hh:76
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
G4ThreeVector GetNewDirection()
G4double SetDifferentialXSection(G4double, G4double, G4int)
static const G4int DIM
int G4int
Definition: G4Types.hh:78
Definition: G4Pow.hh:56
const G4ParticleDefinition * particle
void SetupKinematic(G4double kinEnergy, G4double Z)