Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4WentzelVIModel.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: G4WentzelVIModel.hh 104307 2017-05-24 09:01:45Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4WentzelVIModel
35 //
36 // Author: V.Ivanchenko
37 //
38 // Creation date: 09.04.2008 from G4MuMscModel
39 //
40 // Modifications:
41 // 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
42 // compute cross sections and sample scattering angle
43 //
44 // Class Description:
45 //
46 // Implementation of the model of multiple scattering based on
47 // G.Wentzel, Z. Phys. 40 (1927) 590.
48 // H.W.Lewis, Phys Rev 78 (1950) 526.
49 // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
50 // L.Urban, CERN-OPEN-2006-077.
51 
52 // -------------------------------------------------------------------
53 //
54 
55 #ifndef G4WentzelVIModel_h
56 #define G4WentzelVIModel_h 1
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
60 #include "G4VMscModel.hh"
61 #include "G4MaterialCutsCouple.hh"
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67 {
68 
69 public:
70 
71  explicit G4WentzelVIModel(G4bool comb=true, const G4String& nam = "WentzelVIUni");
72 
73  virtual ~G4WentzelVIModel();
74 
75  virtual void Initialise(const G4ParticleDefinition*,
76  const G4DataVector&) override;
77 
78  virtual void InitialiseLocal(const G4ParticleDefinition*,
79  G4VEmModel* masterModel) override;
80 
81  virtual void StartTracking(G4Track*) override;
82 
84  G4double KineticEnergy,
85  G4double AtomicNumber,
86  G4double AtomicWeight=0.,
87  G4double cut = DBL_MAX,
88  G4double emax= DBL_MAX) override;
89 
91  G4double safety) override;
92 
93  virtual G4double
95  G4double& currentMinimalStep) override;
96 
97  virtual G4double ComputeGeomPathLength(G4double truePathLength) override;
98 
99  virtual G4double ComputeTrueStepLength(G4double geomStepLength) override;
100 
101  // defines low energy limit on energy transfer to atomic electron
102  inline void SetFixedCut(G4double);
103 
104  // low energy limit on energy transfer to atomic electron
105  inline G4double GetFixedCut() const;
106 
107  // access to cross section class
109 
111 
112  inline void SetUseSecondMoment(G4bool);
113 
114  inline G4bool UseSecondMoment() const;
115 
117 
119  const G4MaterialCutsCouple*,
121 
123 
125 
126 protected:
127 
129 
130  inline void SetupParticle(const G4ParticleDefinition*);
131 
132 private:
133 
135  G4double kineticEnergy);
136 
137  // hide assignment operator
139  G4WentzelVIModel(const G4WentzelVIModel&) = delete;
140 
141 protected:
142 
144 
148 
149  // cache kinematics
156 
157  // cache material
161 
165 
168 
169  // cache kinematics
171 
172  // single scattering parameters
175 
177  size_t idx2;
178 
179  // data for single scattering mode
182  std::vector<G4double> xsecn;
183  std::vector<G4double> prob;
184 
187 
188  // projectile
190 
191  // flags
195 };
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199 
201 {
202  // Initialise mass and charge
203  if(p != particle) {
204  particle = p;
205  wokvi->SetupParticle(p);
206  }
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
210 
212 {
213  fixedCut = val;
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217 
219 {
220  return fixedCut;
221 }
222 
223 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224 
226 {
227  if(ptr != wokvi) {
228  delete wokvi;
229  wokvi = ptr;
230  }
231 }
232 
233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234 
236 {
237  return wokvi;
238 }
239 
240 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
241 
243 {
244  useSecondMoment = val;
245 }
246 
247 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
248 
250 {
251  return useSecondMoment;
252 }
253 
254 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
255 
257 {
258  return fSecondMoments;
259 }
260 
261 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
262 
263 inline G4double
265  const G4MaterialCutsCouple* couple,
266  G4double ekin)
267 {
268  G4double x = 0.0;
269  if(useSecondMoment) {
270  DefineMaterial(couple);
271  x = (fSecondMoments) ?
272  (*fSecondMoments)[(*theDensityIdx)[currentMaterialIndex]]->Value(ekin, idx2)
273  *(*theDensityFactor)[currentMaterialIndex]/(ekin*ekin)
274  : ComputeSecondMoment(part, ekin);
275  }
276  return x;
277 }
278 
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280 
281 #endif
282 
Float_t x
Definition: compare.C:6
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
G4WentzelOKandVIxSection * GetWVICrossSection()
void DefineMaterial(const G4MaterialCutsCouple *)
virtual ~G4WentzelVIModel()
G4WentzelVIModel & operator=(const G4WentzelVIModel &right)=delete
G4WentzelVIModel(G4bool comb=true, const G4String &nam="WentzelVIUni")
G4PhysicsTable * GetSecondMomentTable()
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
const char * p
Definition: xmltok.h:285
G4ParticleChangeForMSC * fParticleChange
std::vector< G4double > prob
void SetSingleScatteringFactor(G4double)
static const G4double emax
G4double SecondMoment(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
void SetWVICrossSection(G4WentzelOKandVIxSection *)
G4double ComputeTransportXSectionPerVolume(G4double cosTheta)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
virtual G4double ComputeTrueStepLength(G4double geomStepLength) override
virtual G4double ComputeGeomPathLength(G4double truePathLength) override
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
const G4MaterialCutsCouple * currentCouple
void SetUseSecondMoment(G4bool)
G4double ComputeSecondMoment(const G4ParticleDefinition *, G4double kineticEnergy)
const G4ParticleDefinition * particle
int G4int
Definition: G4Types.hh:78
virtual void StartTracking(G4Track *) override
G4bool UseSecondMoment() const
void SetFixedCut(G4double)
void SetupParticle(const G4ParticleDefinition *)
G4PhysicsTable * fSecondMoments
const G4DataVector * currentCuts
G4double GetFixedCut() const
const G4Material * currentMaterial
#define DBL_MAX
Definition: templates.hh:83
std::vector< G4double > xsecn