Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
61class G4PhysicsTable;
64class G4Region;
65class G4EmSaturation;
70class G4VEmProcess;
71class G4EmCorrections;
75
77{
78
80
81public:
82
84
86
87 //-------------------------------------------------
88 // initialisation before a new run
89 //-------------------------------------------------
90
91 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
92 G4VEnergyLossProcess* p, G4bool theMaster);
93
94 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
95 G4VEmProcess* p, G4bool theMaster);
96
97 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
98 G4VMultipleScattering* p, G4bool theMaster);
99
100 void BuildPhysicsTable(const G4ParticleDefinition* aParticle);
101
102 void BuildPhysicsTable(const G4ParticleDefinition* aParticle,
104
105 void LocalPhysicsTables(const G4ParticleDefinition* aParticle,
107
108 void DumpHtml();
109
110 //-------------------------------------------------
111 // Run time access to DEDX, range, energy for a given particle,
112 // energy, and G4MaterialCutsCouple
113 //-------------------------------------------------
114
115 inline G4double GetDEDX(
116 const G4ParticleDefinition *aParticle,
117 G4double kineticEnergy,
118 const G4MaterialCutsCouple *couple);
119
120 inline G4double GetSubDEDX(
121 const G4ParticleDefinition *aParticle,
122 G4double kineticEnergy,
123 const G4MaterialCutsCouple *couple);
124
125 inline G4double GetRange(
126 const G4ParticleDefinition *aParticle,
127 G4double kineticEnergy,
128 const G4MaterialCutsCouple *couple);
129
130 inline G4double GetCSDARange(
131 const G4ParticleDefinition *aParticle,
132 G4double kineticEnergy,
133 const G4MaterialCutsCouple *couple);
134
136 const G4ParticleDefinition *aParticle,
137 G4double kineticEnergy,
138 const G4MaterialCutsCouple *couple);
139
140 inline G4double GetEnergy(
141 const G4ParticleDefinition *aParticle,
142 G4double range,
143 const G4MaterialCutsCouple *couple);
144
146 const G4MaterialCutsCouple *couple,
147 const G4DynamicParticle* dp,
148 G4double& length);
149
150 //-------------------------------------------------
151 // Methods to be called only at initialisation
152 // and at the end of the job
153 //-------------------------------------------------
154
156
158
160
162
163 void Register(G4VEmProcess* p);
164
165 void DeRegister(G4VEmProcess* p);
166
167 void Register(G4VProcess* p);
168
169 void DeRegister(G4VProcess* p);
170
171 void Register(G4VEmModel* p);
172
173 void DeRegister(G4VEmModel* p);
174
176
178
179 void RegisterExtraParticle(const G4ParticleDefinition* aParticle,
181
182 void SetVerbose(G4int val);
183
184 void ResetParameters();
185
187
189
191
192 //-------------------------------------------------
193 // Access methods
194 //-------------------------------------------------
195
196 inline G4bool IsMaster() const;
197
199
200 const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector();
201
202 const std::vector<G4VEmProcess*>& GetEmProcessVector();
203
204 const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector();
205
207
209
211
213
215
217
219
221
223
225
227
229
231
233
234private:
235
236 //-------------------------------------------------
237 // Private methods and members
238 //-------------------------------------------------
239
241
242 void Clear();
243
244 G4VEnergyLossProcess* BuildTables(const G4ParticleDefinition* aParticle);
245
246 void CopyTables(const G4ParticleDefinition* aParticle,
248
249 void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle);
250
251 void CopyDEDXTables();
252
253 void PrintEWarning(G4String, G4double);
254
256 G4LossTableManager & operator=(const G4LossTableManager &right) = delete;
257
258 static G4ThreadLocal G4LossTableManager* instance;
259
260 typedef const G4ParticleDefinition* PD;
261
262 std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map;
263
264 std::vector<G4VEnergyLossProcess*> loss_vector;
265 std::vector<PD> part_vector;
266 std::vector<PD> base_part_vector;
267 std::vector<G4bool> tables_are_built;
268 std::vector<G4bool> isActive;
269 std::vector<G4PhysicsTable*> dedx_vector;
270 std::vector<G4PhysicsTable*> range_vector;
271 std::vector<G4PhysicsTable*> inv_range_vector;
272 std::vector<G4VMultipleScattering*> msc_vector;
273 std::vector<G4VEmProcess*> emp_vector;
274 std::vector<G4VEmModel*> mod_vector;
275 std::vector<G4VEmFluctuationModel*> fmod_vector;
276 std::vector<G4VProcess*> p_vector;
277
278 // cache
279 G4VEnergyLossProcess* currentLoss;
280 PD currentParticle;
281 PD theElectron;
282 PD theGenericIon;
283 PD firstParticle;
284
285 G4int n_loss;
286 G4int run;
287
288 G4bool all_tables_are_built;
289 G4bool startInitialisation;
290 G4bool isMaster;
291 G4LossTableBuilder* tableBuilder;
292 G4EmCorrections* emCorrections;
293 G4EmConfigurator* emConfigurator;
294 G4ElectronIonPair* emElectronIonPair;
295 G4NIELCalculator* nielCalculator;
296 G4VAtomDeexcitation* atomDeexcitation;
297 G4VSubCutProducer* subcutProducer;
298
299 G4EmParameters* theParameters;
300 G4VEmProcess* gGeneral;
301 G4VEmProcess* eGeneral;
302 G4VEmProcess* pGeneral;
303
304 G4int verbose;
305};
306
307//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
308//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
309
310inline
312 G4double kineticEnergy,
313 const G4MaterialCutsCouple *couple)
314{
315 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
316 return currentLoss ? currentLoss->GetDEDX(kineticEnergy, couple) : 0.0;
317}
318
319//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
320
321inline
323 G4double kineticEnergy,
324 const G4MaterialCutsCouple *couple)
325{
326 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
327 return currentLoss ? currentLoss->GetDEDXForSubsec(kineticEnergy, couple) : 0.0;
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
331
332inline
334 G4double kineticEnergy,
335 const G4MaterialCutsCouple *couple)
336{
337 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
338 return currentLoss ? currentLoss->GetCSDARange(kineticEnergy, couple) : DBL_MAX;
339}
340
341//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
342
343inline
345 const G4ParticleDefinition *aParticle,
346 G4double kineticEnergy,
347 const G4MaterialCutsCouple *couple)
348{
349 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
350 return currentLoss ? currentLoss->GetRangeForLoss(kineticEnergy, couple) : DBL_MAX;
351}
352
353//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
354
355inline
357 G4double kineticEnergy,
358 const G4MaterialCutsCouple *couple)
359{
360 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
361 return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
362}
363
364//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
365
366inline
368 G4double range,
369 const G4MaterialCutsCouple *couple)
370{
371 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
372 return currentLoss ? currentLoss->GetKineticEnergy(range, couple) : 0.0;
373}
374
375//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
376
377inline
379 const G4MaterialCutsCouple *couple,
380 const G4DynamicParticle* dp,
381 G4double& length)
382{
383 const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
384 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
385 return currentLoss ? currentLoss->GetDEDXDispersion(couple, dp, length) : 0.0;
386}
387
388//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
389
391{
392 return isMaster;
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
396
398{
399 return emCorrections;
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
403
405{
406 return atomDeexcitation;
407}
408
409//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
410
412{
413 return subcutProducer;
414}
415
416//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
417
419{
420 return tableBuilder;
421}
422
423//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
424
426{
427 gGeneral = ptr;
428}
429
430//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
431
433{
434 return gGeneral;
435}
436
437//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
438
440{
441 eGeneral = ptr;
442}
443
444//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
445
447{
448 return eGeneral;
449}
450
451//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
452
454{
455 pGeneral = ptr;
456}
457
458//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
459
461{
462 return pGeneral;
463}
464//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
465
466#endif
467
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4ParticleDefinition * GetParticleDefinition() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const std::vector< G4VEmProcess * > & GetEmProcessVector()
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void SetGammaGeneralProcess(G4VEmProcess *)
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4LossTableBuilder * GetTableBuilder()
G4double GetSubDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
void SetPositronGeneralProcess(G4VEmProcess *)
G4VSubCutProducer * SubCutProducer()
G4double GetEnergy(const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4VEmProcess * GetGammaGeneralProcess()
void SetVerbose(G4int val)
void DeRegister(G4VEnergyLossProcess *p)
G4NIELCalculator * NIELCalculator()
void SetNIELCalculator(G4NIELCalculator *)
G4VEmProcess * GetPositronGeneralProcess()
G4EmConfigurator * EmConfigurator()
void Register(G4VEnergyLossProcess *p)
G4ElectronIonPair * ElectronIonPair()
void SetElectronGeneralProcess(G4VEmProcess *)
G4EmSaturation * EmSaturation()
G4VEmProcess * GetElectronGeneralProcess()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
G4VAtomDeexcitation * AtomDeexcitation()
G4double GetRangeFromRestricteDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4EmCorrections * EmCorrections()
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)
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetRangeForLoss(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
G4double GetDEDXForSubsec(G4double kineticEnergy, 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