Geant4-11
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// -------------------------------------------------------------------
27//
28// GEANT4 Class header file
29//
30//
31// File name: G4LossTableManager
32//
33// Author: Vladimir Ivanchenko on base of G4LossTables class
34// and Maria Grazia Pia ideas
35//
36// Creation date: 03.01.2002
37//
38// Modifications by V.Ivanchenko
39//
40// Class Description:
41//
42// A utility static class, responsable for the energy loss tables
43// for each particle
44//
45// Energy loss processes have to register their tables with this
46// class. The responsibility of creating and deleting the tables
47// remains with the energy loss classes.
48
49// -------------------------------------------------------------------
50//
51
52#ifndef G4LossTableManager_h
53#define G4LossTableManager_h 1
54
55#include <map>
56#include <vector>
57#include "globals.hh"
60#include "G4EmParameters.hh"
61
62class G4PhysicsTable;
65class G4Region;
66class G4EmSaturation;
71class G4VEmProcess;
72class G4EmCorrections;
76
78{
79
81
82public:
83
85
87
88 //-------------------------------------------------
89 // initialisation before a new run
90 //-------------------------------------------------
91
92 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
93 G4VEnergyLossProcess* p, G4bool theMaster);
94
95 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
96 G4VEmProcess* p, G4bool theMaster);
97
98 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
99 G4VMultipleScattering* p, G4bool theMaster);
100
101 void BuildPhysicsTable(const G4ParticleDefinition* aParticle);
102
103 void BuildPhysicsTable(const G4ParticleDefinition* aParticle,
105
106 void LocalPhysicsTables(const G4ParticleDefinition* aParticle,
108
109 void DumpHtml();
110
111 //-------------------------------------------------
112 // Run time access to DEDX, range, energy for a given particle,
113 // energy, and G4MaterialCutsCouple
114 //-------------------------------------------------
115
116 inline G4double GetDEDX(
117 const G4ParticleDefinition *aParticle,
118 G4double kineticEnergy,
119 const G4MaterialCutsCouple *couple);
120
121 inline G4double GetRange(
122 const G4ParticleDefinition *aParticle,
123 G4double kineticEnergy,
124 const G4MaterialCutsCouple *couple);
125
126 inline G4double GetCSDARange(
127 const G4ParticleDefinition *aParticle,
128 G4double kineticEnergy,
129 const G4MaterialCutsCouple *couple);
130
132 const G4ParticleDefinition *aParticle,
133 G4double kineticEnergy,
134 const G4MaterialCutsCouple *couple);
135
136 inline G4double GetEnergy(
137 const G4ParticleDefinition *aParticle,
138 G4double range,
139 const G4MaterialCutsCouple *couple);
140
142 const G4MaterialCutsCouple *couple,
143 const G4DynamicParticle* dp,
144 G4double& length);
145
146 //-------------------------------------------------
147 // Methods to be called only at initialisation
148 // and at the end of the job
149 //-------------------------------------------------
150
152
154
156
158
159 void Register(G4VEmProcess* p);
160
161 void DeRegister(G4VEmProcess* p);
162
163 void Register(G4VProcess* p);
164
165 void DeRegister(G4VProcess* p);
166
167 void Register(G4VEmModel* p);
168
169 void DeRegister(G4VEmModel* p);
170
172
174
175 void RegisterExtraParticle(const G4ParticleDefinition* aParticle,
177
178 void SetVerbose(G4int val);
179
180 void ResetParameters();
181
183
185
187
188 //-------------------------------------------------
189 // Access methods
190 //-------------------------------------------------
191
192 inline G4bool IsMaster() const;
193
195
196 const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector();
197
198 const std::vector<G4VEmProcess*>& GetEmProcessVector();
199
200 const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector();
201
203
205
207
209
211
213
215
217
219
221
223
225
227
229
232
233private:
234
235 //-------------------------------------------------
236 // Private methods and members
237 //-------------------------------------------------
238
240
241 void Clear();
242
244
245 void CopyTables(const G4ParticleDefinition* aParticle,
247
248 void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle);
249
251
253
255
256 typedef const G4ParticleDefinition* PD;
257
258 std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map;
259
260 std::vector<G4VEnergyLossProcess*> loss_vector;
261 std::vector<PD> part_vector;
262 std::vector<PD> base_part_vector;
263 std::vector<G4bool> tables_are_built;
264 std::vector<G4bool> isActive;
265 std::vector<G4PhysicsTable*> dedx_vector;
266 std::vector<G4PhysicsTable*> range_vector;
267 std::vector<G4PhysicsTable*> inv_range_vector;
268 std::vector<G4VMultipleScattering*> msc_vector;
269 std::vector<G4VEmProcess*> emp_vector;
270 std::vector<G4VEmModel*> mod_vector;
271 std::vector<G4VEmFluctuationModel*> fmod_vector;
272 std::vector<G4VProcess*> p_vector;
273
274 // cache
280
288
293
297
301};
302
303//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
305
306inline
308 G4double kineticEnergy,
309 const G4MaterialCutsCouple *couple)
310{
311 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
312 return currentLoss ? currentLoss->GetDEDX(kineticEnergy, couple) : 0.0;
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
316
317inline
319 G4double kineticEnergy,
320 const G4MaterialCutsCouple *couple)
321{
322 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
323 return currentLoss ? currentLoss->GetCSDARange(kineticEnergy, couple) : DBL_MAX;
324}
325
326//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
327
328inline
330 const G4ParticleDefinition *aParticle,
331 G4double kineticEnergy,
332 const G4MaterialCutsCouple *couple)
333{
334 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
335 return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
336}
337
338//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
339
340inline
342 G4double kineticEnergy,
343 const G4MaterialCutsCouple *couple)
344{
345 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
346 return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
350
351inline
353 G4double range,
354 const G4MaterialCutsCouple *couple)
355{
356 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
357 return currentLoss ? currentLoss->GetKineticEnergy(range, couple) : 0.0;
358}
359
360//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
361
362inline
364 const G4MaterialCutsCouple *couple,
365 const G4DynamicParticle* dp,
366 G4double& length)
367{
368 const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
369 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
370 return currentLoss ? currentLoss->GetDEDXDispersion(couple, dp, length) : 0.0;
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
374
376{
377 return isMaster;
378}
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
381
383{
384 return emCorrections;
385}
386
387//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
388
390{
391 return atomDeexcitation;
392}
393
394//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
395
397{
398 return subcutProducer;
399}
400
401//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
402
404{
405 return tableBuilder;
406}
407
408//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
409
411{
412 gGeneral = ptr;
413}
414
415//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
416
418{
419 return gGeneral;
420}
421
422//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
423
425{
426 eGeneral = ptr;
427}
428
429//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
430
432{
433 return eGeneral;
434}
435
436//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
437
439{
440 pGeneral = ptr;
441}
442
443//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
444
446{
447 return pGeneral;
448}
449//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
450
451#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4ParticleDefinition * GetParticleDefinition() const
std::vector< PD > base_part_vector
std::vector< G4PhysicsTable * > range_vector
void SetAtomDeexcitation(G4VAtomDeexcitation *)
const G4ParticleDefinition * PD
std::vector< G4bool > isActive
std::vector< G4PhysicsTable * > dedx_vector
static G4LossTableManager * Instance()
const std::vector< G4VEmProcess * > & GetEmProcessVector()
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void CopyTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
G4ElectronIonPair * emElectronIonPair
std::vector< G4VEmProcess * > emp_vector
void ParticleHaveNoLoss(const G4ParticleDefinition *aParticle)
G4VSubCutProducer * subcutProducer
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
void SetGammaGeneralProcess(G4VEmProcess *)
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4VEnergyLossProcess * currentLoss
G4EmParameters * theParameters
G4LossTableBuilder * GetTableBuilder()
std::vector< G4VMultipleScattering * > msc_vector
G4EmConfigurator * emConfigurator
std::vector< G4VEmModel * > mod_vector
static G4ThreadLocal G4LossTableManager * instance
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
void SetPositronGeneralProcess(G4VEmProcess *)
std::vector< G4VEnergyLossProcess * > loss_vector
G4VSubCutProducer * SubCutProducer()
std::vector< G4VProcess * > p_vector
G4double GetEnergy(const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4VEmProcess * GetGammaGeneralProcess()
void SetVerbose(G4int val)
G4LossTableBuilder * tableBuilder
void DeRegister(G4VEnergyLossProcess *p)
G4NIELCalculator * NIELCalculator()
void SetNIELCalculator(G4NIELCalculator *)
G4EmCorrections * emCorrections
G4VEmProcess * GetPositronGeneralProcess()
G4EmConfigurator * EmConfigurator()
void Register(G4VEnergyLossProcess *p)
G4ElectronIonPair * ElectronIonPair()
void SetElectronGeneralProcess(G4VEmProcess *)
std::vector< G4VEmFluctuationModel * > fmod_vector
void PrintEWarning(G4String, G4double)
G4EmSaturation * EmSaturation()
G4VEmProcess * GetElectronGeneralProcess()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
G4VAtomDeexcitation * AtomDeexcitation()
std::vector< G4bool > tables_are_built
G4double GetRangeFromRestricteDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4LossTableManager & operator=(const G4LossTableManager &right)=delete
std::vector< PD > part_vector
G4NIELCalculator * nielCalculator
G4EmCorrections * EmCorrections()
G4VEnergyLossProcess * BuildTables(const G4ParticleDefinition *aParticle)
std::vector< G4PhysicsTable * > inv_range_vector
void SetSubCutProducer(G4VSubCutProducer *)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &length)
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4LossTableManager(G4LossTableManager &)=delete
G4VAtomDeexcitation * atomDeexcitation
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
G4double GetCSDARange(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77