Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id$
27 //
28 //---------------------------------------------------------------------------
29 //
30 // ClassName: G4EnergyLossForExtrapolator
31 //
32 // Description: This class provide calculation of energy loss, fluctuation,
33 // and msc angle
34 //
35 // Author: 09.12.04 V.Ivanchenko
36 //
37 // Modification:
38 // 08-04-05 Rename Propogator -> Extrapolator
39 // 16-03-06 Add muon tables
40 // 21-03-06 Add verbosity defined in the constructor and Initialisation
41 // start only when first public method is called (V.Ivanchenko)
42 // 03-05-06 Remove unused pointer G4Material* from number of methods (VI)
43 // 28-07-07 Add maxEnergyTransfer for computation of energy loss (VI)
44 //
45 //----------------------------------------------------------------------------
46 //
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 
50 #ifndef G4EnergyLossForExtrapolator_h
51 #define G4EnergyLossForExtrapolator_h 1
52 
53 #include <vector>
55 
56 #include "globals.hh"
57 #include "G4PhysicsTable.hh"
58 
60 class G4Material;
62 class G4ProductionCuts;
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65 
67 {
68 public:
70 
72 
74 
76 
78 
79  G4double EnergyAfterStep(G4double kinEnergy, G4double step,
80  const G4Material*, const G4ParticleDefinition*);
81 
82  G4double EnergyBeforeStep(G4double kinEnergy, G4double step,
83  const G4Material*, const G4ParticleDefinition*);
84 
85  G4double TrueStepLength(G4double kinEnergy, G4double step,
86  const G4Material*, const G4ParticleDefinition* part);
87 
88  inline G4double EnergyAfterStep(G4double kinEnergy, G4double step,
89  const G4Material*, const G4String& particleName);
90 
91  inline G4double EnergyBeforeStep(G4double kinEnergy, G4double step,
92  const G4Material*, const G4String& particleName);
93 
94  inline G4double AverageScatteringAngle(G4double kinEnergy, G4double step,
95  const G4Material*,
96  const G4ParticleDefinition* part);
97 
98  inline G4double AverageScatteringAngle(G4double kinEnergy, G4double step,
99  const G4Material*,
100  const G4String& particleName);
101 
102  inline G4double ComputeTrueStep(const G4Material*, const G4ParticleDefinition* part,
103  G4double kinEnergy, G4double stepLength);
104 
105  inline G4double EnergyDispersion(G4double kinEnergy, G4double step,
106  const G4Material*, const G4ParticleDefinition*);
107 
108  inline G4double EnergyDispersion(G4double kinEnergy, G4double step,
109  const G4Material*, const G4String& particleName);
110 
111  inline void SetVerbose(G4int val);
112 
113  inline void SetMinKinEnergy(G4double);
114 
115  inline void SetMaxKinEnergy(G4double);
116 
117  inline void SetMaxEnergyTransfer(G4double);
118 
119 private:
120 
121  void Initialisation();
122 
123  G4bool SetupKinematics(const G4ParticleDefinition*, const G4Material*,
124  G4double kinEnergy);
125 
126  G4PhysicsTable* PrepareTable();
127 
128  const G4ParticleDefinition* FindParticle(const G4String& name);
129 
130  void ComputeElectronDEDX(const G4ParticleDefinition* part, G4PhysicsTable* table);
131 
132  void ComputeMuonDEDX(const G4ParticleDefinition* part, G4PhysicsTable* table);
133 
134  void ComputeProtonDEDX(const G4ParticleDefinition* part, G4PhysicsTable* table);
135 
136  void ComputeTrasportXS(const G4ParticleDefinition* part, G4PhysicsTable* table);
137 
138  inline G4double ComputeValue(G4double x, const G4PhysicsTable* table);
139 
140  // hide assignment operator
143 
144  const G4ParticleDefinition* currentParticle;
145  const G4ParticleDefinition* electron;
146  const G4ParticleDefinition* positron;
147  const G4ParticleDefinition* muonPlus;
148  const G4ParticleDefinition* muonMinus;
149  const G4ParticleDefinition* proton;
150 
151  G4ProductionCuts* cuts;
152  std::vector<const G4MaterialCutsCouple*> couples;
153 
154  G4String currentParticleName;
155 
156  G4PhysicsTable* dedxElectron;
157  G4PhysicsTable* dedxPositron;
158  G4PhysicsTable* dedxMuon;
159  G4PhysicsTable* dedxProton;
160  G4PhysicsTable* rangeElectron;
161  G4PhysicsTable* rangePositron;
162  G4PhysicsTable* rangeMuon;
163  G4PhysicsTable* rangeProton;
164  G4PhysicsTable* invRangeElectron;
165  G4PhysicsTable* invRangePositron;
166  G4PhysicsTable* invRangeMuon;
167  G4PhysicsTable* invRangeProton;
168  G4PhysicsTable* mscElectron;
169 
170  const G4Material* currentMaterial;
171  G4int index;
172  G4double electronDensity;
173  G4double radLength;
174  G4double mass;
175  G4double charge2;
176  G4double kineticEnergy;
177  G4double gam;
178  G4double bg2;
179  G4double beta2;
180  G4double tmax;
181 
182  G4double linLossLimit;
183  G4double emin;
184  G4double emax;
185  G4double maxEnergyTransfer;
186 
187  G4int nbins;
188  G4int nmat;
189  G4int verbose;
190  G4bool isInitialised;
191 };
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194 
196  G4double step,
197  const G4Material* mat,
198  const G4String& name)
199 {
200  return EnergyAfterStep(kinEnergy,step,mat,FindParticle(name));
201 }
202 
203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
204 
206  G4double step,
207  const G4Material* mat,
208  const G4String& name)
209 {
210  return EnergyBeforeStep(kinEnergy,step,mat,FindParticle(name));
211 }
212 
213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214 
216  G4double step,
217  const G4Material* mat,
218  const G4String& name)
219 {
220  return AverageScatteringAngle(kinEnergy,step,mat,FindParticle(name));
221 }
222 
223 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224 
225 inline G4double
227  G4double step,
228  const G4Material* mat,
229  const G4String& name)
230 {
231  return EnergyDispersion(kinEnergy,step,mat,FindParticle(name));
232 }
233 
234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
235 
236 inline G4double
238  G4double stepLength,
239  const G4Material* mat,
240  const G4ParticleDefinition* part)
241 {
242  G4double theta = 0.0;
243  if(SetupKinematics(part, mat, kinEnergy)) {
244  G4double t = stepLength/radLength;
245  G4double y = std::max(0.001, t);
246  theta = 19.23*CLHEP::MeV*std::sqrt(charge2*t)*(1.0 + 0.038*std::log(y))/(beta2*gam*mass);
247  }
248  return theta;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252 
253 inline G4double
255  const G4ParticleDefinition* part,
256  G4double kinEnergy,
257  G4double stepLength)
258 {
259  G4double theta = AverageScatteringAngle(kinEnergy,stepLength,mat,part);
260  return stepLength*std::sqrt(1.0 + 0.625*theta*theta);
261 }
262 
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
264 
265 inline
267  G4double stepLength,
268  const G4Material* mat,
269  const G4ParticleDefinition* part)
270 {
271  G4double sig2 = 0.0;
272  if(SetupKinematics(part, mat, kinEnergy)) {
273  G4double step = ComputeTrueStep(mat,part,kinEnergy,stepLength);
274  sig2 = (1.0/beta2 - 0.5)*CLHEP::twopi_mc2_rcl2*tmax*step*electronDensity*charge2;
275  }
276  return sig2;
277 }
278 
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
280 
281 inline G4double G4EnergyLossForExtrapolator::ComputeValue(G4double x,
282  const G4PhysicsTable* table)
283 {
284  G4double res = 0.0;
285  G4bool b;
286  if(table) res = ((*table)[index])->GetValue(x, b);
287  return res;
288 }
289 
290 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
291 
293 {
294  verbose = val;
295 }
296 
297 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
298 
300 {
301  emin = val;
302 }
303 
304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
305 
307 {
308  emax = val;
309 }
310 
311 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
312 
314 {
315  maxEnergyTransfer = val;
316 }
317 
318 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
319 
320 #endif
321