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

Last change on this file since 1127 was 1055, checked in by garnier, 17 years ago

maj sur la beta de geant 4.9.3

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