source: trunk/source/processes/electromagnetic/utils/include/G4LossTableManager.hh@ 1337

Last change on this file since 1337 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 14.9 KB
Line 
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,v 1.58 2010/04/27 16:59:52 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
28//
29//
30// -------------------------------------------------------------------
31//
32// GEANT4 Class header file
33//
34//
35// File name: G4LossTableManager
36//
37// Author: Vladimir Ivanchenko on base of G4LossTables class
38// and Maria Grazia Pia ideas
39//
40// Creation date: 03.01.2002
41//
42// Modifications:
43//
44// 20-01-03 Migrade to cut per region (V.Ivanchenko)
45// 17-02-03 Fix problem of store/restore tables for ions (V.Ivanchenko)
46// 10-03-03 Add Ion registration (V.Ivanchenko)
47// 25-03-03 Add deregistration (V.Ivanchenko)
48// 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
49// 02-04-03 Change messenger (V.Ivanchenko)
50// 23-07-03 Add exchange with G4EnergyLossTables (V.Ivanchenko)
51// 05-10-03 Add G4VEmProcesses registration (V.Ivanchenko)
52// 17-10-03 Add SetParameters method (V.Ivanchenko)
53// 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
54// 14-01-04 Activate precise range calculation (V.Ivanchenko)
55// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
56// 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
57// 20-01-06 Introduce GetSubDEDX method (VI)
58// 26-01-06 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
59// 10-05-06 Add methods SetMscStepLimitation, FacRange and MscFlag (VI)
60// 22-05-06 Add methods Set/Get bremsTh (VI)
61// 12-02-07 Add SetSkin, SetLinearLossLimit (V.Ivanchenko)
62// 18-06-07 Move definition of msc parameters to G4EmProcessOptions (V.Ivanchenko)
63// 12-04-10 Added PreparePhsyicsTables and BuildPhysicsTables entries (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"
83#include "G4VEnergyLossProcess.hh"
84#include "G4EnergyLossTables.hh"
85
86class G4PhysicsTable;
87class G4MaterialCutsCouple;
88class G4EnergyLossMessenger;
89class G4ParticleDefinition;
90class G4VMultipleScattering;
91class G4VEmProcess;
92class G4EmCorrections;
93class G4EmSaturation;
94class G4EmConfigurator;
95class G4LossTableBuilder;
96class G4Region;
97
98class G4LossTableManager
99{
100
101public:
102
103 static G4LossTableManager* Instance();
104
105 ~G4LossTableManager();
106
107 //-------------------------------------------------
108 // called from destructor
109 //-------------------------------------------------
110
111 void Clear();
112
113 //-------------------------------------------------
114 // initialisation before a new run
115 //-------------------------------------------------
116
117 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
118 G4VEnergyLossProcess* p);
119 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
120 G4VEmProcess* p);
121 void PreparePhysicsTable(const G4ParticleDefinition* aParticle,
122 G4VMultipleScattering* p);
123 void BuildPhysicsTable(const G4ParticleDefinition* aParticle);
124 void BuildPhysicsTable(const G4ParticleDefinition* aParticle,
125 G4VEnergyLossProcess* p);
126
127 //-------------------------------------------------
128 // Run time access to DEDX, range, energy for a given particle,
129 // energy, and G4MaterialCutsCouple
130 //-------------------------------------------------
131
132 inline G4double GetDEDX(
133 const G4ParticleDefinition *aParticle,
134 G4double kineticEnergy,
135 const G4MaterialCutsCouple *couple);
136
137 inline G4double GetSubDEDX(
138 const G4ParticleDefinition *aParticle,
139 G4double kineticEnergy,
140 const G4MaterialCutsCouple *couple);
141
142 inline G4double GetRange(
143 const G4ParticleDefinition *aParticle,
144 G4double kineticEnergy,
145 const G4MaterialCutsCouple *couple);
146
147 inline G4double GetCSDARange(
148 const G4ParticleDefinition *aParticle,
149 G4double kineticEnergy,
150 const G4MaterialCutsCouple *couple);
151
152 inline G4double GetRangeFromRestricteDEDX(
153 const G4ParticleDefinition *aParticle,
154 G4double kineticEnergy,
155 const G4MaterialCutsCouple *couple);
156
157 inline G4double GetEnergy(
158 const G4ParticleDefinition *aParticle,
159 G4double range,
160 const G4MaterialCutsCouple *couple);
161
162 inline G4double GetDEDXDispersion(
163 const G4MaterialCutsCouple *couple,
164 const G4DynamicParticle* dp,
165 G4double& length);
166
167 //-------------------------------------------------
168 // Methods to be called only at initialisation
169 //-------------------------------------------------
170
171 void Register(G4VEnergyLossProcess* p);
172
173 void DeRegister(G4VEnergyLossProcess* p);
174
175 void Register(G4VMultipleScattering* p);
176
177 void DeRegister(G4VMultipleScattering* p);
178
179 void Register(G4VEmProcess* p);
180
181 void DeRegister(G4VEmProcess* p);
182
183 void Register(G4VEmModel* p);
184
185 void DeRegister(G4VEmModel* p);
186
187 void Register(G4VEmFluctuationModel* p);
188
189 void DeRegister(G4VEmFluctuationModel* p);
190
191 void RegisterIon(const G4ParticleDefinition* aParticle,
192 G4VEnergyLossProcess* p);
193
194 void RegisterExtraParticle(const G4ParticleDefinition* aParticle,
195 G4VEnergyLossProcess* p);
196
197 void SetLossFluctuations(G4bool val);
198
199 void SetSubCutoff(G4bool val, const G4Region* r=0);
200
201 void SetIntegral(G4bool val);
202
203 void SetRandomStep(G4bool val);
204
205 void SetMinSubRange(G4double val);
206
207 void SetMinEnergy(G4double val);
208
209 void SetMaxEnergy(G4double val);
210
211 void SetMaxEnergyForCSDARange(G4double val);
212
213 void SetMaxEnergyForMuons(G4double val);
214
215 void SetDEDXBinning(G4int val);
216
217 void SetDEDXBinningForCSDARange(G4int val);
218
219 void SetLambdaBinning(G4int val);
220
221 G4int GetNumberOfBinsPerDecade() const;
222
223 void SetStepFunction(G4double v1, G4double v2);
224
225 void SetBuildCSDARange(G4bool val);
226
227 void SetLPMFlag(G4bool val);
228
229 void SetSplineFlag(G4bool val);
230
231 void SetLinearLossLimit(G4double val);
232
233 void SetBremsstrahlungTh(G4double val);
234
235 void SetFactorForAngleLimit(G4double val);
236
237 void SetVerbose(G4int val);
238
239 //-------------------------------------------------
240 // Access methods
241 //-------------------------------------------------
242
243 G4EnergyLossMessenger* GetMessenger();
244
245 G4bool BuildCSDARange() const;
246
247 G4bool LPMFlag() const;
248
249 G4bool SplineFlag() const;
250
251 G4double BremsstrahlungTh() const;
252
253 G4double FactorForAngleLimit() const;
254
255 const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector();
256
257 const std::vector<G4VEmProcess*>& GetEmProcessVector();
258
259 const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector();
260
261 inline G4VEnergyLossProcess* GetEnergyLossProcess(const G4ParticleDefinition*);
262
263 G4EmCorrections* EmCorrections();
264
265 G4EmSaturation* EmSaturation();
266
267 G4EmConfigurator* EmConfigurator();
268
269private:
270
271 //-------------------------------------------------
272 // Private methods and members
273 //-------------------------------------------------
274
275 G4LossTableManager();
276
277 G4VEnergyLossProcess* BuildTables(const G4ParticleDefinition* aParticle);
278
279 void CopyTables(const G4ParticleDefinition* aParticle,
280 G4VEnergyLossProcess*);
281
282 void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle);
283
284 void SetParameters(const G4ParticleDefinition* aParticle,
285 G4VEnergyLossProcess*);
286
287 void CopyDEDXTables();
288
289 static G4LossTableManager* theInstance;
290
291 typedef const G4ParticleDefinition* PD;
292
293 std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map;
294
295 std::vector<G4VEnergyLossProcess*> loss_vector;
296 std::vector<PD> part_vector;
297 std::vector<PD> base_part_vector;
298 std::vector<G4bool> tables_are_built;
299 std::vector<G4bool> isActive;
300 std::vector<G4PhysicsTable*> dedx_vector;
301 std::vector<G4PhysicsTable*> range_vector;
302 std::vector<G4PhysicsTable*> inv_range_vector;
303 std::vector<G4VMultipleScattering*> msc_vector;
304 std::vector<G4VEmProcess*> emp_vector;
305 std::vector<G4VEmModel*> mod_vector;
306 std::vector<G4VEmFluctuationModel*> fmod_vector;
307
308 // cash
309 G4VEnergyLossProcess* currentLoss;
310 PD currentParticle;
311 PD theElectron;
312 PD firstParticle;
313
314 G4int n_loss;
315 G4int run;
316
317 G4bool all_tables_are_built;
318 G4bool startInitialisation;
319
320 G4bool lossFluctuationFlag;
321 G4bool subCutoffFlag;
322 G4bool rndmStepFlag;
323 G4bool integral;
324 G4bool integralActive;
325 G4bool buildCSDARange;
326 G4bool minEnergyActive;
327 G4bool maxEnergyActive;
328 G4bool maxEnergyForMuonsActive;
329 G4bool stepFunctionActive;
330 G4bool flagLPM;
331 G4bool splineFlag;
332
333 G4double minSubRange;
334 G4double maxRangeVariation;
335 G4double maxFinalStep;
336 G4double minKinEnergy;
337 G4double maxKinEnergy;
338 G4double maxKinEnergyForMuons;
339 G4double bremsTh;
340 G4double factorForAngleLimit;
341
342 G4LossTableBuilder* tableBuilder;
343 G4EnergyLossMessenger* theMessenger;
344 G4EmCorrections* emCorrections;
345 G4EmSaturation* emSaturation;
346 G4EmConfigurator* emConfigurator;
347
348 G4int nbinsLambda;
349 G4int nbinsPerDecade;
350 G4int verbose;
351
352};
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
355//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
356
357inline G4VEnergyLossProcess* G4LossTableManager::GetEnergyLossProcess(
358 const G4ParticleDefinition *aParticle)
359{
360 if(aParticle != currentParticle) {
361 currentParticle = aParticle;
362 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
363 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
364 currentLoss = (*pos).second;
365 } else {
366 currentLoss = 0;
367 // ParticleHaveNoLoss(aParticle);
368 }
369 }
370 return currentLoss;
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374
375inline G4double G4LossTableManager::GetDEDX(
376 const G4ParticleDefinition *aParticle,
377 G4double kineticEnergy,
378 const G4MaterialCutsCouple *couple)
379{
380 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
381 G4double x;
382 if(currentLoss) { x = currentLoss->GetDEDX(kineticEnergy, couple); }
383 else { x = G4EnergyLossTables::GetDEDX(currentParticle,
384 kineticEnergy,couple,false); }
385 return x;
386}
387
388//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
389
390inline G4double G4LossTableManager::GetSubDEDX(
391 const G4ParticleDefinition *aParticle,
392 G4double kineticEnergy,
393 const G4MaterialCutsCouple *couple)
394{
395 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
396 G4double x = 0.0;
397 if(currentLoss) { x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple); }
398 return x;
399}
400
401//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
402
403inline G4double G4LossTableManager::GetCSDARange(
404 const G4ParticleDefinition *aParticle,
405 G4double kineticEnergy,
406 const G4MaterialCutsCouple *couple)
407{
408 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
409 G4double x = DBL_MAX;
410 if(currentLoss) { x = currentLoss->GetCSDARange(kineticEnergy, couple); }
411 return x;
412}
413
414//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
415
416inline G4double G4LossTableManager::GetRangeFromRestricteDEDX(
417 const G4ParticleDefinition *aParticle,
418 G4double kineticEnergy,
419 const G4MaterialCutsCouple *couple)
420{
421 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
422 G4double x;
423 if(currentLoss) { x = currentLoss->GetRangeForLoss(kineticEnergy, couple); }
424 else { x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,
425 couple,false); }
426 return x;
427}
428
429//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
430
431inline G4double G4LossTableManager::GetRange(
432 const G4ParticleDefinition *aParticle,
433 G4double kineticEnergy,
434 const G4MaterialCutsCouple *couple)
435{
436 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
437 G4double x;
438 if(currentLoss) { x = currentLoss->GetRange(kineticEnergy, couple); }
439 else { x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,
440 couple,false); }
441 return x;
442}
443
444//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
445
446inline G4double G4LossTableManager::GetEnergy(
447 const G4ParticleDefinition *aParticle,
448 G4double range,
449 const G4MaterialCutsCouple *couple)
450{
451 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
452 G4double x;
453 if(currentLoss) { x = currentLoss->GetKineticEnergy(range, couple); }
454 else { x = G4EnergyLossTables::GetPreciseEnergyFromRange(currentParticle,range,
455 couple,false); }
456 return x;
457}
458
459//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
460
461inline G4double G4LossTableManager::GetDEDXDispersion(
462 const G4MaterialCutsCouple *couple,
463 const G4DynamicParticle* dp,
464 G4double& length)
465{
466 const G4ParticleDefinition* aParticle = dp->GetDefinition();
467 if(aParticle != currentParticle) {
468 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
469 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
470 currentParticle = aParticle;
471 currentLoss = (*pos).second;
472 } else {
473 ParticleHaveNoLoss(aParticle);
474 }
475 }
476 return currentLoss->GetDEDXDispersion(couple, dp, length);
477}
478
479//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
480
481#endif
482
Note: See TracBrowser for help on using the repository browser.