Geant4  10.03
G4LossTableManager.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: G4LossTableManager.hh 99151 2016-09-07 08:03:17Z gcosmo $
27 //
28 //
29 // -------------------------------------------------------------------
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4LossTableManager
35 //
36 // Author: Vladimir Ivanchenko on base of G4LossTables class
37 // and Maria Grazia Pia ideas
38 //
39 // Creation date: 03.01.2002
40 //
41 // Modifications:
42 //
43 // 20-01-03 Migrade to cut per region (V.Ivanchenko)
44 // 17-02-03 Fix problem of store/restore tables for ions (V.Ivanchenko)
45 // 10-03-03 Add Ion registration (V.Ivanchenko)
46 // 25-03-03 Add deregistration (V.Ivanchenko)
47 // 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
48 // 02-04-03 Change messenger (V.Ivanchenko)
49 // 23-07-03 Add exchange with G4EnergyLossTables (V.Ivanchenko)
50 // 05-10-03 Add G4VEmProcesses registration (V.Ivanchenko)
51 // 17-10-03 Add SetParameters method (V.Ivanchenko)
52 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
53 // 14-01-04 Activate precise range calculation (V.Ivanchenko)
54 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
55 // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
56 // 20-01-06 Introduce GetSubDEDX method (VI)
57 // 26-01-06 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
58 // 10-05-06 Add methods SetMscStepLimitation, FacRange and MscFlag (VI)
59 // 22-05-06 Add methods Set/Get bremsTh (VI)
60 // 12-02-07 Add SetSkin, SetLinearLossLimit (V.Ivanchenko)
61 // 18-06-07 Move definition of msc parameters to G4EmProcessOptions (V.Ivanchenko)
62 // 12-04-10 Added PreparePhsyicsTables and BuildPhysicsTables entries (V.Ivanchenko)
63 // 04-06-13 Adaptation for MT mode, new method LocalPhysicsTables (V.Ivanchenko)
64 //
65 // Class Description:
66 //
67 // A utility static class, responsable for the energy loss tables
68 // for each particle
69 //
70 // Energy loss processes have to register their tables with this
71 // class. The responsibility of creating and deleting the tables
72 // remains with the energy loss classes.
73 
74 // -------------------------------------------------------------------
75 //
76 
77 #ifndef G4LossTableManager_h
78 #define G4LossTableManager_h 1
79 
80 #include <map>
81 #include <vector>
82 #include "globals.hh"
84 #include "G4VEnergyLossProcess.hh"
85 #include "G4EmCorrections.hh"
86 #include "G4LossTableBuilder.hh"
87 #include "G4VAtomDeexcitation.hh"
88 #include "G4VSubCutProducer.hh"
89 
90 class G4PhysicsTable;
93 class G4Region;
94 class G4EmSaturation;
95 class G4EmConfigurator;
96 class G4ElectronIonPair;
98 class G4VEmProcess;
99 
101 {
102 
104 
105 public:
106 
107  static G4LossTableManager* Instance();
108 
110 
111  //-------------------------------------------------
112  // called from destructor
113  //-------------------------------------------------
114 
115  void Clear();
116 
117  //-------------------------------------------------
118  // initialisation before a new run
119  //-------------------------------------------------
120 
121  void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
122  G4VEnergyLossProcess* p, G4bool theMaster);
123 
124  void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
125  G4VEmProcess* p, G4bool theMaster);
126 
127  void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
128  G4VMultipleScattering* p, G4bool theMaster);
129 
130  void BuildPhysicsTable(const G4ParticleDefinition* aParticle);
131 
132  void BuildPhysicsTable(const G4ParticleDefinition* aParticle,
134 
135  void LocalPhysicsTables(const G4ParticleDefinition* aParticle,
137 
138  //-------------------------------------------------
139  // Run time access to DEDX, range, energy for a given particle,
140  // energy, and G4MaterialCutsCouple
141  //-------------------------------------------------
142 
143  inline G4double GetDEDX(
144  const G4ParticleDefinition *aParticle,
145  G4double kineticEnergy,
146  const G4MaterialCutsCouple *couple);
147 
148  inline G4double GetSubDEDX(
149  const G4ParticleDefinition *aParticle,
150  G4double kineticEnergy,
151  const G4MaterialCutsCouple *couple);
152 
153  inline G4double GetRange(
154  const G4ParticleDefinition *aParticle,
155  G4double kineticEnergy,
156  const G4MaterialCutsCouple *couple);
157 
158  inline G4double GetCSDARange(
159  const G4ParticleDefinition *aParticle,
160  G4double kineticEnergy,
161  const G4MaterialCutsCouple *couple);
162 
164  const G4ParticleDefinition *aParticle,
165  G4double kineticEnergy,
166  const G4MaterialCutsCouple *couple);
167 
168  inline G4double GetEnergy(
169  const G4ParticleDefinition *aParticle,
170  G4double range,
171  const G4MaterialCutsCouple *couple);
172 
174  const G4MaterialCutsCouple *couple,
175  const G4DynamicParticle* dp,
176  G4double& length);
177 
178  //-------------------------------------------------
179  // Methods to be called only at initialisation
180  //-------------------------------------------------
181 
183 
185 
187 
189 
190  void Register(G4VEmProcess* p);
191 
192  void DeRegister(G4VEmProcess* p);
193 
194  void Register(G4VEmModel* p);
195 
196  void DeRegister(G4VEmModel* p);
197 
199 
201 
202  void RegisterExtraParticle(const G4ParticleDefinition* aParticle,
204 
205  void SetVerbose(G4int val);
206 
208 
210 
211  //-------------------------------------------------
212  // Access methods
213  //-------------------------------------------------
214 
215  inline G4bool IsMaster() const;
216 
217  inline
219 
220  const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector();
221 
222  const std::vector<G4VEmProcess*>& GetEmProcessVector();
223 
224  const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector();
225 
227 
229 
231 
232  inline G4EmCorrections* EmCorrections();
233 
235 
237 
239 
240 private:
241 
242  //-------------------------------------------------
243  // Private methods and members
244  //-------------------------------------------------
245 
247 
248  void ResetParameters();
249 
251 
252  void CopyTables(const G4ParticleDefinition* aParticle,
254 
255  void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle);
256 
257  void CopyDEDXTables();
258 
260 
263 
265 
266  typedef const G4ParticleDefinition* PD;
267 
268  std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map;
269 
270  std::vector<G4VEnergyLossProcess*> loss_vector;
271  std::vector<PD> part_vector;
272  std::vector<PD> base_part_vector;
273  std::vector<G4bool> tables_are_built;
274  std::vector<G4bool> isActive;
275  std::vector<G4PhysicsTable*> dedx_vector;
276  std::vector<G4PhysicsTable*> range_vector;
277  std::vector<G4PhysicsTable*> inv_range_vector;
278  std::vector<G4VMultipleScattering*> msc_vector;
279  std::vector<G4VEmProcess*> emp_vector;
280  std::vector<G4VEmModel*> mod_vector;
281  std::vector<G4VEmFluctuationModel*> fmod_vector;
282 
283  // cache
289 
292 
302 
304 
306 
307 };
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
310 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
311 
312 inline G4VEnergyLossProcess*
314 {
315  //G4cout << "G4LossTableManager::GetEnergyLossProcess: "
316  //<< aParticle << " " << currentParticle << " " << currentLoss << G4endl;
317  if(aParticle != currentParticle) {
318  currentParticle = aParticle;
319  std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
320  if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
321  currentLoss = (*pos).second;
322  } else {
323  currentLoss = 0;
324  if ((pos = loss_map.find(theGenericIon)) != loss_map.end()) {
325  currentLoss = (*pos).second;
326  }
327  }
328  }
329  return currentLoss;
330 }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333 
334 inline
336  G4double kineticEnergy,
337  const G4MaterialCutsCouple *couple)
338 {
339  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
340  G4double x = 0.0;
341  if(currentLoss) { x = currentLoss->GetDEDX(kineticEnergy, couple); }
342  return x;
343 }
344 
345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
346 
347 inline
349  G4double kineticEnergy,
350  const G4MaterialCutsCouple *couple)
351 {
352  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
353  G4double x = 0.0;
354  if(currentLoss) { x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple); }
355  return x;
356 }
357 
358 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
359 
360 inline
362  G4double kineticEnergy,
363  const G4MaterialCutsCouple *couple)
364 {
365  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
366  G4double x = DBL_MAX;
367  if(currentLoss) { x = currentLoss->GetCSDARange(kineticEnergy, couple); }
368  return x;
369 }
370 
371 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
372 
373 inline
375  const G4ParticleDefinition *aParticle,
376  G4double kineticEnergy,
377  const G4MaterialCutsCouple *couple)
378 {
379  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
380  G4double x = DBL_MAX;
381  if(currentLoss) { x = currentLoss->GetRangeForLoss(kineticEnergy, couple); }
382  return x;
383 }
384 
385 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
386 
387 inline
389  G4double kineticEnergy,
390  const G4MaterialCutsCouple *couple)
391 {
392  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
393  G4double x = DBL_MAX;
394  if(currentLoss) {
395  x = currentLoss->GetRange(kineticEnergy, couple);
396  }
397  return x;
398 }
399 
400 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
401 
402 inline
404  G4double range,
405  const G4MaterialCutsCouple *couple)
406 {
407  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
408  G4double x = 0;
409  if(currentLoss) { x = currentLoss->GetKineticEnergy(range, couple); }
410  return x;
411 }
412 
413 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
414 
415 inline
417  const G4MaterialCutsCouple *couple,
418  const G4DynamicParticle* dp,
419  G4double& length)
420 {
421  const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
422  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
423  G4double x = 0.0;
424  if(currentLoss) { currentLoss->GetDEDXDispersion(couple, dp, length); }
425  return x;
426 }
427 
428 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
429 
431 {
432  return isMaster;
433 }
434 
435 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
436 
438 {
439  return emCorrections;
440 }
441 
442 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
443 
445 {
446  return atomDeexcitation;
447 }
448 
449 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
450 
452 {
453  return subcutProducer;
454 }
455 
456 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
457 
459 {
460  return tableBuilder;
461 }
462 
463 #endif
464 
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4VSubCutProducer * SubCutProducer()
G4EmConfigurator * EmConfigurator()
static G4LossTableManager * Instance()
G4VAtomDeexcitation * atomDeexcitation
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4VSubCutProducer * subcutProducer
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4VEmModel * > mod_vector
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &length)
G4ElectronIonPair * emElectronIonPair
void DeRegister(G4VEnergyLossProcess *p)
G4double GetSubDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
static G4ThreadLocal G4LossTableManager * instance
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void CopyTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
std::vector< G4PhysicsTable * > range_vector
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
#define G4ThreadLocal
Definition: tls.hh:89
const std::vector< G4VEmProcess * > & GetEmProcessVector()
int G4int
Definition: G4Types.hh:78
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4LossTableBuilder * GetTableBuilder()
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
G4EmConfigurator * emConfigurator
G4VEnergyLossProcess * BuildTables(const G4ParticleDefinition *aParticle)
G4LossTableManager & operator=(const G4LossTableManager &right)=delete
void ParticleHaveNoLoss(const G4ParticleDefinition *aParticle)
void SetSubCutProducer(G4VSubCutProducer *)
std::vector< G4VEmFluctuationModel * > fmod_vector
bool G4bool
Definition: G4Types.hh:79
G4EmCorrections * EmCorrections()
G4VEnergyLossProcess * currentLoss
G4ElectronIonPair * ElectronIonPair()
G4EmSaturation * EmSaturation()
G4bool IsMaster() const
void Register(G4VEnergyLossProcess *p)
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetEnergy(const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
G4double GetRangeFromRestricteDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
std::vector< G4PhysicsTable * > inv_range_vector
const G4ParticleDefinition * PD
G4EmParameters * theParameters
G4EmCorrections * emCorrections
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4LossTableBuilder * tableBuilder
std::vector< G4PhysicsTable * > dedx_vector
std::vector< G4VEnergyLossProcess * > loss_vector
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
void PrintEWarning(G4String, G4double)
void SetVerbose(G4int val)
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
std::vector< PD > part_vector
std::vector< G4bool > isActive
G4VAtomDeexcitation * AtomDeexcitation()
std::vector< G4VEmProcess * > emp_vector
double G4double
Definition: G4Types.hh:76
std::vector< PD > base_part_vector
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
std::vector< G4VMultipleScattering * > msc_vector
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static const G4double pos
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4bool > tables_are_built