Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4LossTableManager.hh
이 파일의 문서화 페이지로 가기
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 110572 2018-05-30 13:08:12Z 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 class G4GammaShark;
100 
102 {
103 
105 
106 public:
107 
108  static G4LossTableManager* Instance();
109 
111 
112  //-------------------------------------------------
113  // called from destructor
114  //-------------------------------------------------
115 
116  void Clear();
117 
118  //-------------------------------------------------
119  // initialisation before a new run
120  //-------------------------------------------------
121 
122  void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
123  G4VEnergyLossProcess* p, G4bool theMaster);
124 
125  void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
126  G4VEmProcess* p, G4bool theMaster);
127 
128  void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
129  G4VMultipleScattering* p, G4bool theMaster);
130 
131  void BuildPhysicsTable(const G4ParticleDefinition* aParticle);
132 
133  void BuildPhysicsTable(const G4ParticleDefinition* aParticle,
135 
136  void LocalPhysicsTables(const G4ParticleDefinition* aParticle,
138 
139  void DumpHtml();
140 
141  //-------------------------------------------------
142  // Run time access to DEDX, range, energy for a given particle,
143  // energy, and G4MaterialCutsCouple
144  //-------------------------------------------------
145 
146  inline G4double GetDEDX(
147  const G4ParticleDefinition *aParticle,
149  const G4MaterialCutsCouple *couple);
150 
151  inline G4double GetSubDEDX(
152  const G4ParticleDefinition *aParticle,
154  const G4MaterialCutsCouple *couple);
155 
156  inline G4double GetRange(
157  const G4ParticleDefinition *aParticle,
159  const G4MaterialCutsCouple *couple);
160 
161  inline G4double GetCSDARange(
162  const G4ParticleDefinition *aParticle,
164  const G4MaterialCutsCouple *couple);
165 
167  const G4ParticleDefinition *aParticle,
169  const G4MaterialCutsCouple *couple);
170 
171  inline G4double GetEnergy(
172  const G4ParticleDefinition *aParticle,
173  G4double range,
174  const G4MaterialCutsCouple *couple);
175 
177  const G4MaterialCutsCouple *couple,
178  const G4DynamicParticle* dp,
179  G4double& length);
180 
181  //-------------------------------------------------
182  // Methods to be called only at initialisation
183  //-------------------------------------------------
184 
186 
188 
190 
192 
193  void Register(G4VEmProcess* p);
194 
195  void DeRegister(G4VEmProcess* p);
196 
197  void Register(G4VEmModel* p);
198 
199  void DeRegister(G4VEmModel* p);
200 
202 
204 
205  void RegisterExtraParticle(const G4ParticleDefinition* aParticle,
207 
208  void SetVerbose(G4int val);
209 
211 
213 
214  //-------------------------------------------------
215  // Access methods
216  //-------------------------------------------------
217 
218  inline G4bool IsMaster() const;
219 
221 
222  const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector();
223 
224  const std::vector<G4VEmProcess*>& GetEmProcessVector();
225 
226  const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector();
227 
229 
231 
233 
234  inline G4EmCorrections* EmCorrections();
235 
237 
239 
241 
242  inline void SetGammaShark(G4GammaShark*);
243 
244  inline G4GammaShark* GetGammaShark();
245 
246 private:
247 
248  //-------------------------------------------------
249  // Private methods and members
250  //-------------------------------------------------
251 
253 
254  void ResetParameters();
255 
257 
258  void CopyTables(const G4ParticleDefinition* aParticle,
260 
261  void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle);
262 
263  void CopyDEDXTables();
264 
266 
269 
271 
272  typedef const G4ParticleDefinition* PD;
273 
274  std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map;
275 
276  std::vector<G4VEnergyLossProcess*> loss_vector;
277  std::vector<PD> part_vector;
278  std::vector<PD> base_part_vector;
279  std::vector<G4bool> tables_are_built;
280  std::vector<G4bool> isActive;
281  std::vector<G4PhysicsTable*> dedx_vector;
282  std::vector<G4PhysicsTable*> range_vector;
283  std::vector<G4PhysicsTable*> inv_range_vector;
284  std::vector<G4VMultipleScattering*> msc_vector;
285  std::vector<G4VEmProcess*> emp_vector;
286  std::vector<G4VEmModel*> mod_vector;
287  std::vector<G4VEmFluctuationModel*> fmod_vector;
288 
289  // cache
295 
298 
308 
310  G4GammaShark* gammaShark;
311 
313 
314 };
315 
316 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
318 
319 inline
322  const G4MaterialCutsCouple *couple)
323 {
324  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
325  return currentLoss ? currentLoss->GetDEDX(kineticEnergy, couple) : 0.0;
326 }
327 
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
329 
330 inline
333  const G4MaterialCutsCouple *couple)
334 {
335  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
336  return currentLoss ? currentLoss->GetDEDXForSubsec(kineticEnergy, couple) : 0.0;
337 }
338 
339 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
340 
341 inline
344  const G4MaterialCutsCouple *couple)
345 {
346  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
347  return currentLoss ? currentLoss->GetCSDARange(kineticEnergy, couple) : DBL_MAX;
348 }
349 
350 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
351 
352 inline
354  const G4ParticleDefinition *aParticle,
356  const G4MaterialCutsCouple *couple)
357 {
358  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
359  return currentLoss ? currentLoss->GetRangeForLoss(kineticEnergy, couple) : DBL_MAX;
360 }
361 
362 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
363 
364 inline
367  const G4MaterialCutsCouple *couple)
368 {
369  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
370  return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
371 }
372 
373 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374 
375 inline
377  G4double range,
378  const G4MaterialCutsCouple *couple)
379 {
380  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
381  return currentLoss ? currentLoss->GetKineticEnergy(range, couple) : 0.0;
382 }
383 
384 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
385 
386 inline
388  const G4MaterialCutsCouple *couple,
389  const G4DynamicParticle* dp,
390  G4double& length)
391 {
392  const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
393  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
394  return currentLoss ? currentLoss->GetDEDXDispersion(couple, dp, length) : 0.0;
395 }
396 
397 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
398 
400 {
401  return isMaster;
402 }
403 
404 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
405 
407 {
408  return emCorrections;
409 }
410 
411 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
412 
414 {
415  return atomDeexcitation;
416 }
417 
418 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
419 
421 {
422  return subcutProducer;
423 }
424 
425 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
426 
428 {
429  return tableBuilder;
430 }
431 
432 inline void G4LossTableManager::SetGammaShark(G4GammaShark* ptr)
433 {
434  gammaShark = ptr;
435 }
436 
437 inline G4GammaShark* G4LossTableManager::GetGammaShark()
438 {
439  return gammaShark;
440 }
441 
442 #endif
443 
std::vector< G4VMultipleScattering * > msc_vector
G4LossTableBuilder * tableBuilder
G4VSubCutProducer * subcutProducer
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4VEmModel * > mod_vector
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4VAtomDeexcitation * AtomDeexcitation()
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double GetSubDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
const char * p
Definition: xmltok.h:285
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
G4bool IsMaster() const
G4LossTableManager & operator=(const G4LossTableManager &right)=delete
std::vector< G4PhysicsTable * > inv_range_vector
static G4ThreadLocal G4LossTableManager * instance
G4EmCorrections * emCorrections
const G4ParticleDefinition * PD
void SetGammaShark(G4GammaShark *)
void SetVerbose(G4int val)
#define G4ThreadLocal
Definition: tls.hh:69
std::vector< G4VEnergyLossProcess * > loss_vector
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
std::vector< G4VEmFluctuationModel * > fmod_vector
G4ElectronIonPair * emElectronIonPair
void ParticleHaveNoLoss(const G4ParticleDefinition *aParticle)
std::vector< G4PhysicsTable * > dedx_vector
G4ElectronIonPair * ElectronIonPair()
void CopyTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
std::vector< G4PhysicsTable * > range_vector
void Register(G4VEnergyLossProcess *p)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4VAtomDeexcitation * atomDeexcitation
G4EmCorrections * EmCorrections()
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
const std::vector< G4VEmProcess * > & GetEmProcessVector()
void SetSubCutProducer(G4VSubCutProducer *)
G4VEnergyLossProcess * currentLoss
const G4ParticleDefinition const G4Material *G4double range
G4EmConfigurator * emConfigurator
std::vector< G4bool > tables_are_built
G4VSubCutProducer * SubCutProducer()
void DeRegister(G4VEnergyLossProcess *p)
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void PrintEWarning(G4String, G4double)
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4LossTableBuilder * GetTableBuilder()
G4VEnergyLossProcess * BuildTables(const G4ParticleDefinition *aParticle)
G4GammaShark * gammaShark
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4GammaShark * GetGammaShark()
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
int G4int
Definition: G4Types.hh:78
G4EmSaturation * EmSaturation()
G4double GetEnergy(const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &length)
G4EmParameters * theParameters
static G4LossTableManager * Instance()
void SetAtomDeexcitation(G4VAtomDeexcitation *)
G4double GetRangeFromRestricteDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4ParticleDefinition * GetParticleDefinition() const
std::vector< G4VEmProcess * > emp_vector
std::vector< G4bool > isActive
#define DBL_MAX
Definition: templates.hh:83
G4EmConfigurator * EmConfigurator()
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
std::vector< PD > base_part_vector
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
std::vector< PD > part_vector
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)