Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EnergyLossForExtrapolator.hh
Go to the documentation of this file.
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//
27//---------------------------------------------------------------------------
28//
29// ClassName: G4EnergyLossForExtrapolator
30//
31// Description: This class provide calculation of energy loss, fluctuation,
32// and msc angle
33//
34// Author: 09.12.04 V.Ivanchenko
35//
36// Modification:
37// 08-04-05 Rename Propogator -> Extrapolator
38// 16-03-06 Add muon tables
39// 21-03-06 Add verbosity defined in the constructor and Initialisation
40// start only when first public method is called (V.Ivanchenko)
41// 03-05-06 Remove unused pointer G4Material* from number of methods (VI)
42// 28-07-07 Add maxEnergyTransfer for computation of energy loss (VI)
43//
44//----------------------------------------------------------------------------
45//
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
49#ifndef G4EnergyLossForExtrapolator_h
50#define G4EnergyLossForExtrapolator_h 1
51
52#include <vector>
54
55#include "globals.hh"
56#include "G4PhysicsTable.hh"
58#include "G4Log.hh"
59#include "G4Threading.hh"
60
62class G4Material;
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68{
69public:
70
71 explicit G4EnergyLossForExtrapolator(G4int verb = 1);
72
74
76
78
80
82 const G4Material*, const G4ParticleDefinition*);
83
85 const G4Material*, const G4ParticleDefinition*);
86
88 const G4Material*, const G4ParticleDefinition* part);
89
90 inline G4double EnergyAfterStep(G4double kinEnergy, G4double step,
91 const G4Material*,
92 const G4String& particleName);
93
94 inline G4double EnergyBeforeStep(G4double kinEnergy, G4double step,
95 const G4Material*,
96 const G4String& particleName);
97
98 inline G4double AverageScatteringAngle(G4double kinEnergy, G4double step,
99 const G4Material*,
100 const G4ParticleDefinition* part);
101
102 inline G4double AverageScatteringAngle(G4double kinEnergy, G4double step,
103 const G4Material*,
104 const G4String& particleName);
105
106 inline G4double ComputeTrueStep(const G4Material*,
107 const G4ParticleDefinition* part,
108 G4double kinEnergy, G4double stepLength);
109
110 inline G4double EnergyDispersion(G4double kinEnergy, G4double step,
111 const G4Material*,
112 const G4ParticleDefinition*);
113
114 inline G4double EnergyDispersion(G4double kinEnergy, G4double step,
115 const G4Material*,
116 const G4String& particleName);
117
118 inline void SetVerbose(G4int val);
119
120 inline void SetMinKinEnergy(G4double);
121
122 inline void SetMaxKinEnergy(G4double);
123
124 inline void SetMaxEnergyTransfer(G4double);
125
126private:
127
128 void Initialisation();
129
130 void BuildTables();
131
132 G4bool SetupKinematics(const G4ParticleDefinition*, const G4Material*,
133 G4double kinEnergy);
134
135 const G4ParticleDefinition* FindParticle(const G4String& name);
136
137 inline G4double ComputeValue(G4double x, const G4PhysicsTable* table,
138 size_t idx);
139
140 inline const G4PhysicsTable* GetPhysicsTable(ExtTableType type) const;
141
142 // hide assignment operator
145
146#ifdef G4MULTITHREADED
147 static G4Mutex extrapolatorMutex;
148#endif
149 static G4TablesForExtrapolator* tables;
150
151 const G4ParticleDefinition* currentParticle;
152 const G4ParticleDefinition* electron;
153 const G4ParticleDefinition* positron;
154 const G4ParticleDefinition* muonPlus;
155 const G4ParticleDefinition* muonMinus;
156 const G4ParticleDefinition* proton;
157
158 G4String currentParticleName;
159
160 size_t idxDedxElectron;
161 size_t idxDedxPositron;
162 size_t idxDedxMuon;
163 size_t idxDedxProton;
164 size_t idxRangeElectron;
165 size_t idxRangePositron;
166 size_t idxRangeMuon;
167 size_t idxRangeProton;
168 size_t idxInvRangeElectron;
169 size_t idxInvRangePositron;
170 size_t idxInvRangeMuon;
171 size_t idxInvRangeProton;
172 size_t idxMscElectron;
173
174 const G4Material* currentMaterial;
175 G4int index;
176
177 G4double electronDensity;
178 G4double radLength;
179 G4double mass;
180 G4double charge2;
181 G4double kineticEnergy;
182 G4double gam;
183 G4double bg2;
184 G4double beta2;
185 G4double tmax;
186
187 G4double linLossLimit;
188 G4double emin;
189 G4double emax;
190 G4double maxEnergyTransfer;
191
192 G4int nbins;
193 G4int nmat;
194 G4int verbose;
195};
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198
199inline const G4PhysicsTable*
200G4EnergyLossForExtrapolator::GetPhysicsTable(ExtTableType type) const
201{
202 return tables->GetPhysicsTable(type);
203}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
206
207inline G4double
209 G4double step,
210 const G4Material* mat,
211 const G4String& name)
212{
213 return EnergyAfterStep(kinEnergy,step,mat,FindParticle(name));
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217
218inline G4double
220 G4double step,
221 const G4Material* mat,
222 const G4String& name)
223{
224 return EnergyBeforeStep(kinEnergy,step,mat,FindParticle(name));
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
228
229inline G4double
231 G4double step,
232 const G4Material* mat,
233 const G4String& name)
234{
235 return AverageScatteringAngle(kinEnergy,step,mat,FindParticle(name));
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
239
240inline G4double
242 G4double step,
243 const G4Material* mat,
244 const G4String& name)
245{
246 return EnergyDispersion(kinEnergy,step,mat,FindParticle(name));
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
250
252 G4double kinEnergy,
253 G4double stepLength,
254 const G4Material* mat,
255 const G4ParticleDefinition* part)
256{
257 G4double theta = 0.0;
258 if(SetupKinematics(part, mat, kinEnergy)) {
259 G4double t = stepLength/radLength;
260 G4double y = std::max(0.001, t);
261 theta = 19.23*CLHEP::MeV*std::sqrt(charge2*t)*(1.0 + 0.038*G4Log(y))
262 /(beta2*gam*mass);
263 }
264 return theta;
265}
266
267//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
268
269inline G4double
271 const G4ParticleDefinition* part,
272 G4double kinEnergy,
273 G4double stepLength)
274{
275 G4double theta = AverageScatteringAngle(kinEnergy,stepLength,mat,part);
276 return stepLength*std::sqrt(1.0 + 0.625*theta*theta);
277}
278
279//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
280
281inline G4double
283 G4double stepLength,
284 const G4Material* mat,
285 const G4ParticleDefinition* part)
286{
287 G4double sig2 = 0.0;
288 if(SetupKinematics(part, mat, kinEnergy)) {
289 G4double step = ComputeTrueStep(mat,part,kinEnergy,stepLength);
290 sig2 = (1.0/beta2 - 0.5)
291 *CLHEP::twopi_mc2_rcl2*tmax*step*electronDensity*charge2;
292 }
293 return sig2;
294}
295
296//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
297
298inline G4double
299G4EnergyLossForExtrapolator::ComputeValue(G4double x,
300 const G4PhysicsTable* table,
301 size_t idx)
302{
303 G4double res = 0.0;
304 if(table) { res = ((*table)[index])->Value(x, idx); }
305 return res;
306}
307
308//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
309
311{
312 verbose = val;
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
316
318{
319 emin = val;
320}
321
322//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
323
325{
326 emax = val;
327}
328
329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
330
332{
333 maxEnergyTransfer = val;
334}
335
336//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
337
338#endif
339
G4double G4Log(G4double x)
Definition: G4Log.hh:226
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double TrueStepLength(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
G4double EnergyDispersion(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4double ComputeRange(G4double kinEnergy, const G4ParticleDefinition *)
G4double EnergyBeforeStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4double AverageScatteringAngle(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *)
G4double ComputeEnergy(G4double range, const G4ParticleDefinition *)
G4double ComputeTrueStep(const G4Material *, const G4ParticleDefinition *part, G4double kinEnergy, G4double stepLength)
G4double EnergyAfterStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const