Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4LindhardSorensenIonModel.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: G4LindhardSorensenIonModel.hh 108808 2018-03-08 09:45:40Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4LindhardSorensenIonModel
34 //
35 // Author: Alexander Bagulya & Vladimir Ivanchenko
36 //
37 // Creation date: 16.04.2018
38 //
39 //
40 // Class Description:
41 //
42 // Implementation of ion ionisation energy loss and delta-electron
43 // production by heavy charged particles according to
44 // J. Lindhard & A.H. Sorensen, Phys. Rev. A 53 (1996) 2443-2455
45 
46 // -------------------------------------------------------------------
47 //
48 
49 #ifndef G4LindhardSorensenIonModel_h
50 #define G4LindhardSorensenIonModel_h 1
51 
53 
54 #include "G4VEmModel.hh"
55 #include "G4NistManager.hh"
56 
57 class G4EmCorrections;
60 
62 {
63 
64 public:
65 
66  explicit G4LindhardSorensenIonModel(const G4ParticleDefinition* p = nullptr,
67  const G4String& nam = "LindhardSorensen");
68 
70 
71  virtual void Initialise(const G4ParticleDefinition*,
72  const G4DataVector&) override;
73 
75  const G4MaterialCutsCouple* couple) override;
76 
78  const G4ParticleDefinition*,
80  G4double cutEnergy,
81  G4double maxEnergy);
82 
84  const G4ParticleDefinition*,
85  G4double kineticEnergy,
87  G4double cutEnergy,
88  G4double maxEnergy) override;
89 
91  const G4ParticleDefinition*,
92  G4double kineticEnergy,
93  G4double cutEnergy,
94  G4double maxEnergy) override;
95 
97  const G4ParticleDefinition*,
98  G4double kineticEnergy,
99  G4double cutEnergy) override;
100 
102  const G4Material* mat,
103  G4double kineticEnergy) override;
104 
106  const G4Material* mat,
107  G4double kineticEnergy) override;
108 
109  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
110  const G4DynamicParticle* dp,
111  G4double& eloss,
112  G4double&,
113  G4double length) override;
114 
115  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
116  const G4MaterialCutsCouple*,
117  const G4DynamicParticle*,
118  G4double tmin,
119  G4double maxEnergy) override;
120 
121 protected:
122 
124  G4double kinEnergy) override;
125 
126  inline G4double GetChargeSquareRatio() const;
127 
128  inline void SetChargeSquareRatio(G4double val);
129 
130 private:
131 
132  void SetupParameters();
133 
134  inline void SetParticle(const G4ParticleDefinition* p);
135 
136  // hide assignment operator
137  G4LindhardSorensenIonModel & operator=
138  (const G4LindhardSorensenIonModel &right) = delete;
140 
142 
148 
159 };
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162 
163 inline void
165 {
166  if(particle != p) {
167  particle = p;
168  SetupParameters();
169  }
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173 
175 {
176  return chargeSquare;
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180 
182 {
183  chargeSquare = val;
184 }
185 
186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187 
188 #endif
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4LindhardSorensenIonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="LindhardSorensen")
const char * p
Definition: xmltok.h:285
static G4LindhardSorensenData * lsdata
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4ParticleChangeForLoss * fParticleChange
Float_t Z
void SetParticle(const G4ParticleDefinition *p)
const G4ParticleDefinition * particle
double G4double
Definition: G4Types.hh:76
double A(double temperature)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
Float_t mat
int G4int
Definition: G4Types.hh:78
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &eloss, G4double &, G4double length) override
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override