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

Last change on this file since 991 was 961, checked in by garnier, 17 years ago

update processes

File size: 13.4 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.53 2008/07/15 16:56:38 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
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//
64// Class Description:
65//
66// A utility static class, responsable for the energy loss tables
67// for each particle
68//
69// Energy loss processes have to register their tables with this
70// class. The responsibility of creating and deleting the tables
71// remains with the energy loss classes.
72
73// -------------------------------------------------------------------
74//
75
76#ifndef G4LossTableManager_h
77#define G4LossTableManager_h 1
78
79#include <map>
80#include <vector>
81#include "globals.hh"
82#include "G4VEnergyLossProcess.hh"
83#include "G4EnergyLossTables.hh"
84
85class G4PhysicsTable;
86class G4MaterialCutsCouple;
87class G4EnergyLossMessenger;
88class G4ParticleDefinition;
89class G4VMultipleScattering;
90class G4VEmProcess;
91class G4EmCorrections;
92class G4EmSaturation;
93class G4LossTableBuilder;
94
95class G4LossTableManager
96{
97
98public:
99
100 static G4LossTableManager* Instance();
101
102 ~G4LossTableManager();
103
104 void Clear();
105
106 // get the DEDX or the range for a given particle/energy/material
107 inline G4double GetDEDX(
108 const G4ParticleDefinition *aParticle,
109 G4double kineticEnergy,
110 const G4MaterialCutsCouple *couple);
111
112 inline G4double GetSubDEDX(
113 const G4ParticleDefinition *aParticle,
114 G4double kineticEnergy,
115 const G4MaterialCutsCouple *couple);
116
117 inline G4double GetRange(
118 const G4ParticleDefinition *aParticle,
119 G4double kineticEnergy,
120 const G4MaterialCutsCouple *couple);
121
122 inline G4double GetCSDARange(
123 const G4ParticleDefinition *aParticle,
124 G4double kineticEnergy,
125 const G4MaterialCutsCouple *couple);
126
127 inline G4double GetRangeFromRestricteDEDX(
128 const G4ParticleDefinition *aParticle,
129 G4double kineticEnergy,
130 const G4MaterialCutsCouple *couple);
131
132 inline G4double GetEnergy(
133 const G4ParticleDefinition *aParticle,
134 G4double range,
135 const G4MaterialCutsCouple *couple);
136
137 inline G4double GetDEDXDispersion(
138 const G4MaterialCutsCouple *couple,
139 const G4DynamicParticle* dp,
140 G4double& length);
141
142 // to be called only by energy loss processes
143 void Register(G4VEnergyLossProcess* p);
144
145 void DeRegister(G4VEnergyLossProcess* p);
146
147 void Register(G4VMultipleScattering* p);
148
149 void DeRegister(G4VMultipleScattering* p);
150
151 void Register(G4VEmProcess* p);
152
153 void DeRegister(G4VEmProcess* p);
154
155 void Register(G4VEmModel* p);
156
157 void DeRegister(G4VEmModel* p);
158
159 void Register(G4VEmFluctuationModel* p);
160
161 void DeRegister(G4VEmFluctuationModel* p);
162
163 void EnergyLossProcessIsInitialised(const G4ParticleDefinition* aParticle,
164 G4VEnergyLossProcess* p);
165
166 void RegisterIon(const G4ParticleDefinition* aParticle,
167 G4VEnergyLossProcess* p);
168
169 void RegisterExtraParticle(const G4ParticleDefinition* aParticle,
170 G4VEnergyLossProcess* p);
171
172 void BuildPhysicsTable(const G4ParticleDefinition* aParticle,
173 G4VEnergyLossProcess* p);
174
175 void SetLossFluctuations(G4bool val);
176
177 void SetSubCutoff(G4bool val);
178
179 void SetIntegral(G4bool val);
180
181 void SetRandomStep(G4bool val);
182
183 void SetMinSubRange(G4double val);
184
185 void SetMinEnergy(G4double val);
186
187 void SetMaxEnergy(G4double val);
188
189 void SetMaxEnergyForCSDARange(G4double val);
190
191 void SetMaxEnergyForMuons(G4double val);
192
193 void SetDEDXBinning(G4int val);
194
195 void SetDEDXBinningForCSDARange(G4int val);
196
197 void SetLambdaBinning(G4int val);
198
199 void SetStepFunction(G4double v1, G4double v2);
200
201 void SetBuildCSDARange(G4bool val);
202
203 void SetLPMFlag(G4bool val);
204
205 void SetSplineFlag(G4bool val);
206
207 void SetLinearLossLimit(G4double val);
208
209 void SetBremsstrahlungTh(G4double val);
210
211 void SetVerbose(G4int val);
212
213 G4EnergyLossMessenger* GetMessenger();
214
215 G4bool BuildCSDARange() const;
216
217 G4bool LPMFlag() const;
218
219 G4bool SplineFlag() const;
220
221 G4double BremsstrahlungTh() const;
222
223 const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector();
224
225 const std::vector<G4VEmProcess*>& GetEmProcessVector();
226
227 const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector();
228
229 inline G4VEnergyLossProcess* GetEnergyLossProcess(const G4ParticleDefinition*);
230
231 G4EmCorrections* EmCorrections();
232
233 G4EmSaturation* EmSaturation();
234
235private:
236
237 G4LossTableManager();
238
239 G4VEnergyLossProcess* BuildTables(const G4ParticleDefinition* aParticle);
240
241 void CopyTables(const G4ParticleDefinition* aParticle,
242 G4VEnergyLossProcess*);
243
244 void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle);
245
246 void SetParameters(G4VEnergyLossProcess*);
247
248 void CopyDEDXTables();
249
250private:
251
252 static G4LossTableManager* theInstance;
253
254 typedef const G4ParticleDefinition* PD;
255 std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map;
256
257 std::vector<G4VEnergyLossProcess*> loss_vector;
258 std::vector<PD> part_vector;
259 std::vector<PD> base_part_vector;
260 std::vector<G4bool> tables_are_built;
261 std::vector<G4bool> isActive;
262 std::vector<G4PhysicsTable*> dedx_vector;
263 std::vector<G4PhysicsTable*> range_vector;
264 std::vector<G4PhysicsTable*> inv_range_vector;
265 std::vector<G4VMultipleScattering*> msc_vector;
266 std::vector<G4VEmProcess*> emp_vector;
267 std::vector<G4VEmModel*> mod_vector;
268 std::vector<G4VEmFluctuationModel*> fmod_vector;
269
270 // cash
271 G4VEnergyLossProcess* currentLoss;
272 PD currentParticle;
273 PD theElectron;
274
275 G4int n_loss;
276 G4int run;
277
278 G4bool all_tables_are_built;
279 // G4bool first_entry;
280 G4bool lossFluctuationFlag;
281 G4bool subCutoffFlag;
282 G4bool rndmStepFlag;
283 G4bool integral;
284 G4bool integralActive;
285 G4bool all_tables_are_stored;
286 G4bool buildCSDARange;
287 G4bool minEnergyActive;
288 G4bool maxEnergyActive;
289 G4bool maxEnergyForMuonsActive;
290 G4bool stepFunctionActive;
291 G4bool flagLPM;
292 G4bool splineFlag;
293
294 G4double minSubRange;
295 G4double maxRangeVariation;
296 G4double maxFinalStep;
297 G4double minKinEnergy;
298 G4double maxKinEnergy;
299 G4double maxKinEnergyForMuons;
300 G4double bremsTh;
301
302 G4LossTableBuilder* tableBuilder;
303 G4EnergyLossMessenger* theMessenger;
304 G4EmCorrections* emCorrections;
305 G4EmSaturation* emSaturation;
306
307 const G4ParticleDefinition* firstParticle;
308 G4int verbose;
309
310};
311
312//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
314
315inline G4double G4LossTableManager::GetDEDX(
316 const G4ParticleDefinition *aParticle,
317 G4double kineticEnergy,
318 const G4MaterialCutsCouple *couple)
319{
320 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
321 G4double x;
322 if(currentLoss) x = currentLoss->GetDEDX(kineticEnergy, couple);
323 else x = G4EnergyLossTables::GetDEDX(
324 currentParticle,kineticEnergy,couple,false);
325 return x;
326}
327
328//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
329
330inline G4double G4LossTableManager::GetSubDEDX(
331 const G4ParticleDefinition *aParticle,
332 G4double kineticEnergy,
333 const G4MaterialCutsCouple *couple)
334{
335 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
336 G4double x = 0.0;
337 if(currentLoss) x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple);
338 return x;
339}
340
341//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
342
343inline G4double G4LossTableManager::GetCSDARange(
344 const G4ParticleDefinition *aParticle,
345 G4double kineticEnergy,
346 const G4MaterialCutsCouple *couple)
347{
348 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
349 G4double x = DBL_MAX;
350 if(currentLoss) x = currentLoss->GetCSDARange(kineticEnergy, couple);
351 return x;
352}
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
355
356inline G4double G4LossTableManager::GetRangeFromRestricteDEDX(
357 const G4ParticleDefinition *aParticle,
358 G4double kineticEnergy,
359 const G4MaterialCutsCouple *couple)
360{
361 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
362 G4double x;
363 if(currentLoss) x = currentLoss->GetRangeForLoss(kineticEnergy, couple);
364 else
365 x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,couple,false);
366 return x;
367}
368
369//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
370
371inline G4double G4LossTableManager::GetRange(
372 const G4ParticleDefinition *aParticle,
373 G4double kineticEnergy,
374 const G4MaterialCutsCouple *couple)
375{
376 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
377 G4double x;
378 if(currentLoss) x = currentLoss->GetRange(kineticEnergy, couple);
379 else
380 x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,couple,false);
381 return x;
382}
383
384//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
385
386inline G4double G4LossTableManager::GetEnergy(
387 const G4ParticleDefinition *aParticle,
388 G4double range,
389 const G4MaterialCutsCouple *couple)
390{
391 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
392 G4double x;
393 if(currentLoss) x = currentLoss->GetKineticEnergy(range, couple);
394 else x = G4EnergyLossTables::GetPreciseEnergyFromRange(
395 currentParticle,range,couple,false);
396 return x;
397}
398
399//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
400
401inline G4double G4LossTableManager::GetDEDXDispersion(
402 const G4MaterialCutsCouple *couple,
403 const G4DynamicParticle* dp,
404 G4double& length)
405{
406 const G4ParticleDefinition* aParticle = dp->GetDefinition();
407 if(aParticle != currentParticle) {
408 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
409 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
410 currentParticle = aParticle;
411 currentLoss = (*pos).second;
412 } else {
413 ParticleHaveNoLoss(aParticle);
414 }
415 }
416 return currentLoss->GetDEDXDispersion(couple, dp, length);
417}
418
419//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
420
421inline G4VEnergyLossProcess* G4LossTableManager::GetEnergyLossProcess(
422 const G4ParticleDefinition *aParticle)
423{
424 if(aParticle != currentParticle) {
425 currentParticle = aParticle;
426 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
427 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
428 currentLoss = (*pos).second;
429 } else {
430 currentLoss = 0;
431 // ParticleHaveNoLoss(aParticle);
432 }
433 }
434 return currentLoss;
435}
436
437//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
438
439#endif
440
Note: See TracBrowser for help on using the repository browser.