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: G4BetheBlochModel.hh,v 1.23 2010/05/27 14:26:17 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-04-beta-01 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class header file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4BetheBlochModel |
---|
35 | // |
---|
36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code |
---|
37 | // |
---|
38 | // Creation date: 03.01.2002 |
---|
39 | // |
---|
40 | // Modifications: |
---|
41 | // |
---|
42 | // 23-12-02 V.Ivanchenko change interface in order to moveto cut per region |
---|
43 | // 24-01-03 Make models region aware (V.Ivanchenko) |
---|
44 | // 13-02-03 Add name (V.Ivanchenko) |
---|
45 | // 12-11-03 Fix for GenericIons (V.Ivanchenko) |
---|
46 | // 24-03-05 Add G4EmCorrections (V.Ivanchenko) |
---|
47 | // 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko) |
---|
48 | // 11-04-04 Move MaxSecondaryEnergy to models (V.Ivanchenko) |
---|
49 | // 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma) |
---|
50 | // 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio, |
---|
51 | // CorrectionsAlongStep needed for ions(V.Ivanchenko) |
---|
52 | |
---|
53 | // |
---|
54 | // Class Description: |
---|
55 | // |
---|
56 | // Implementation of Bethe-Bloch model of energy loss and |
---|
57 | // delta-electron production by heavy charged particles |
---|
58 | |
---|
59 | // ------------------------------------------------------------------- |
---|
60 | // |
---|
61 | |
---|
62 | #ifndef G4BetheBlochModel_h |
---|
63 | #define G4BetheBlochModel_h 1 |
---|
64 | |
---|
65 | #include "G4VEmModel.hh" |
---|
66 | #include "G4NistManager.hh" |
---|
67 | |
---|
68 | class G4EmCorrections; |
---|
69 | class G4ParticleChangeForLoss; |
---|
70 | |
---|
71 | class G4BetheBlochModel : public G4VEmModel |
---|
72 | { |
---|
73 | |
---|
74 | public: |
---|
75 | |
---|
76 | G4BetheBlochModel(const G4ParticleDefinition* p = 0, |
---|
77 | const G4String& nam = "BetheBloch"); |
---|
78 | |
---|
79 | virtual ~G4BetheBlochModel(); |
---|
80 | |
---|
81 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
82 | |
---|
83 | virtual G4double MinEnergyCut(const G4ParticleDefinition*, |
---|
84 | const G4MaterialCutsCouple*); |
---|
85 | |
---|
86 | virtual G4double ComputeCrossSectionPerElectron( |
---|
87 | const G4ParticleDefinition*, |
---|
88 | G4double kineticEnergy, |
---|
89 | G4double cutEnergy, |
---|
90 | G4double maxEnergy); |
---|
91 | |
---|
92 | virtual G4double ComputeCrossSectionPerAtom( |
---|
93 | const G4ParticleDefinition*, |
---|
94 | G4double kineticEnergy, |
---|
95 | G4double Z, G4double A, |
---|
96 | G4double cutEnergy, |
---|
97 | G4double maxEnergy); |
---|
98 | |
---|
99 | virtual G4double CrossSectionPerVolume(const G4Material*, |
---|
100 | const G4ParticleDefinition*, |
---|
101 | G4double kineticEnergy, |
---|
102 | G4double cutEnergy, |
---|
103 | G4double maxEnergy); |
---|
104 | |
---|
105 | virtual G4double ComputeDEDXPerVolume(const G4Material*, |
---|
106 | const G4ParticleDefinition*, |
---|
107 | G4double kineticEnergy, |
---|
108 | G4double cutEnergy); |
---|
109 | |
---|
110 | virtual G4double GetChargeSquareRatio(const G4ParticleDefinition* p, |
---|
111 | const G4Material* mat, |
---|
112 | G4double kineticEnergy); |
---|
113 | |
---|
114 | virtual G4double GetParticleCharge(const G4ParticleDefinition* p, |
---|
115 | const G4Material* mat, |
---|
116 | G4double kineticEnergy); |
---|
117 | |
---|
118 | virtual void CorrectionsAlongStep(const G4MaterialCutsCouple* couple, |
---|
119 | const G4DynamicParticle* dp, |
---|
120 | G4double& eloss, |
---|
121 | G4double&, |
---|
122 | G4double length); |
---|
123 | |
---|
124 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
125 | const G4MaterialCutsCouple*, |
---|
126 | const G4DynamicParticle*, |
---|
127 | G4double tmin, |
---|
128 | G4double maxEnergy); |
---|
129 | |
---|
130 | protected: |
---|
131 | |
---|
132 | virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*, |
---|
133 | G4double kinEnergy); |
---|
134 | |
---|
135 | inline G4double GetChargeSquareRatio() const; |
---|
136 | |
---|
137 | inline void SetChargeSquareRatio(G4double val); |
---|
138 | |
---|
139 | private: |
---|
140 | |
---|
141 | void SetupParameters(); |
---|
142 | |
---|
143 | inline void SetParticle(const G4ParticleDefinition* p); |
---|
144 | |
---|
145 | inline void SetGenericIon(const G4ParticleDefinition* p); |
---|
146 | |
---|
147 | // hide assignment operator |
---|
148 | G4BetheBlochModel & operator=(const G4BetheBlochModel &right); |
---|
149 | G4BetheBlochModel(const G4BetheBlochModel&); |
---|
150 | |
---|
151 | const G4ParticleDefinition* particle; |
---|
152 | const G4Material* currentMaterial; |
---|
153 | G4ParticleDefinition* theElectron; |
---|
154 | G4EmCorrections* corr; |
---|
155 | G4ParticleChangeForLoss* fParticleChange; |
---|
156 | G4NistManager* nist; |
---|
157 | |
---|
158 | G4double mass; |
---|
159 | G4double tlimit; |
---|
160 | G4double spin; |
---|
161 | G4double magMoment2; |
---|
162 | G4double chargeSquare; |
---|
163 | G4double ratio; |
---|
164 | G4double formfact; |
---|
165 | G4double twoln10; |
---|
166 | G4double bg2lim; |
---|
167 | G4double taulim; |
---|
168 | G4double corrFactor; |
---|
169 | G4bool isIon; |
---|
170 | G4bool isInitialised; |
---|
171 | }; |
---|
172 | |
---|
173 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
174 | |
---|
175 | inline void G4BetheBlochModel::SetParticle(const G4ParticleDefinition* p) |
---|
176 | { |
---|
177 | if(particle != p) { |
---|
178 | particle = p; |
---|
179 | if (p->GetPDGCharge()/eplus > 1.5 && p->GetBaryonNumber() > 2) |
---|
180 | { isIon = true; } |
---|
181 | SetupParameters(); |
---|
182 | } |
---|
183 | } |
---|
184 | |
---|
185 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
186 | |
---|
187 | inline void G4BetheBlochModel::SetGenericIon(const G4ParticleDefinition* p) |
---|
188 | { |
---|
189 | if(p && particle != p) { |
---|
190 | if(p->GetParticleName() == "GenericIon") { isIon = true; } |
---|
191 | } |
---|
192 | } |
---|
193 | |
---|
194 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
195 | |
---|
196 | inline G4double G4BetheBlochModel::GetChargeSquareRatio() const |
---|
197 | { |
---|
198 | return chargeSquare; |
---|
199 | } |
---|
200 | |
---|
201 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
202 | |
---|
203 | inline void G4BetheBlochModel::SetChargeSquareRatio(G4double val) |
---|
204 | { |
---|
205 | chargeSquare = val; |
---|
206 | } |
---|
207 | |
---|
208 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
209 | |
---|
210 | #endif |
---|