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

Last change on this file since 889 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 13.2 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.48 2007/11/07 18:38:49 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-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
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 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 EnergyLossProcessIsInitialised(const G4ParticleDefinition* aParticle,
156 G4VEnergyLossProcess* p);
157
158 void RegisterIon(const G4ParticleDefinition* aParticle,
159 G4VEnergyLossProcess* p);
160
161 void RegisterExtraParticle(const G4ParticleDefinition* aParticle,
162 G4VEnergyLossProcess* p);
163
164 void BuildPhysicsTable(const G4ParticleDefinition* aParticle,
165 G4VEnergyLossProcess* p);
166
167 void SetLossFluctuations(G4bool val);
168
169 void SetSubCutoff(G4bool val);
170
171 void SetIntegral(G4bool val);
172
173 void SetRandomStep(G4bool val);
174
175 void SetMinSubRange(G4double val);
176
177 void SetMinEnergy(G4double val);
178
179 void SetMaxEnergy(G4double val);
180
181 void SetMaxEnergyForCSDARange(G4double val);
182
183 void SetMaxEnergyForMuons(G4double val);
184
185 void SetDEDXBinning(G4int val);
186
187 void SetDEDXBinningForCSDARange(G4int val);
188
189 void SetLambdaBinning(G4int val);
190
191 void SetStepFunction(G4double v1, G4double v2);
192
193 void SetBuildCSDARange(G4bool val);
194
195 void SetLPMFlag(G4bool val);
196
197 void SetLinearLossLimit(G4double val);
198
199 void SetBremsstrahlungTh(G4double val);
200
201 void SetVerbose(G4int val);
202
203 G4EnergyLossMessenger* GetMessenger();
204
205 G4bool BuildCSDARange() const;
206
207 G4bool LPMFlag() const;
208
209 G4double BremsstrahlungTh() const;
210
211 const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector();
212
213 const std::vector<G4VEmProcess*>& GetEmProcessVector();
214
215 const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector();
216
217 inline G4VEnergyLossProcess* GetEnergyLossProcess(const G4ParticleDefinition*);
218
219 inline G4EmCorrections* EmCorrections();
220
221private:
222
223 G4LossTableManager();
224
225 G4VEnergyLossProcess* BuildTables(const G4ParticleDefinition* aParticle);
226
227 void CopyTables(const G4ParticleDefinition* aParticle,
228 G4VEnergyLossProcess*);
229
230 void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle);
231
232 void SetParameters(G4VEnergyLossProcess*);
233
234 void CopyDEDXTables();
235
236private:
237
238 static G4LossTableManager* theInstance;
239
240 typedef const G4ParticleDefinition* PD;
241 std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map;
242
243 std::vector<G4VEnergyLossProcess*> loss_vector;
244 std::vector<PD> part_vector;
245 std::vector<PD> base_part_vector;
246 std::vector<G4bool> tables_are_built;
247 std::vector<G4bool> isActive;
248 std::vector<G4PhysicsTable*> dedx_vector;
249 std::vector<G4PhysicsTable*> range_vector;
250 std::vector<G4PhysicsTable*> inv_range_vector;
251 std::vector<G4VMultipleScattering*> msc_vector;
252 std::vector<G4VEmProcess*> emp_vector;
253
254 // cash
255 G4VEnergyLossProcess* currentLoss;
256 PD currentParticle;
257 PD theElectron;
258
259 G4int n_loss;
260 G4int run;
261
262 G4bool all_tables_are_built;
263 // G4bool first_entry;
264 G4bool lossFluctuationFlag;
265 G4bool subCutoffFlag;
266 G4bool rndmStepFlag;
267 G4bool integral;
268 G4bool integralActive;
269 G4bool all_tables_are_stored;
270 G4bool buildCSDARange;
271 G4bool minEnergyActive;
272 G4bool maxEnergyActive;
273 G4bool maxEnergyForMuonsActive;
274 G4bool stepFunctionActive;
275 G4bool flagLPM;
276
277 G4double minSubRange;
278 G4double maxRangeVariation;
279 G4double maxFinalStep;
280 G4double minKinEnergy;
281 G4double maxKinEnergy;
282 G4double maxKinEnergyForMuons;
283 G4double bremsTh;
284
285 G4LossTableBuilder* tableBuilder;
286 G4EnergyLossMessenger* theMessenger;
287 G4EmCorrections* emCorrections;
288 const G4ParticleDefinition* firstParticle;
289 G4int verbose;
290
291};
292
293//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
294//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
295
296inline G4double G4LossTableManager::GetDEDX(
297 const G4ParticleDefinition *aParticle,
298 G4double kineticEnergy,
299 const G4MaterialCutsCouple *couple)
300{
301 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
302 G4double x;
303 if(currentLoss) x = currentLoss->GetDEDX(kineticEnergy, couple);
304 else x = G4EnergyLossTables::GetDEDX(
305 currentParticle,kineticEnergy,couple,false);
306 return x;
307}
308
309//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
310
311inline G4double G4LossTableManager::GetSubDEDX(
312 const G4ParticleDefinition *aParticle,
313 G4double kineticEnergy,
314 const G4MaterialCutsCouple *couple)
315{
316 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
317 G4double x = 0.0;
318 if(currentLoss) x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple);
319 return x;
320}
321
322//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
323
324inline G4double G4LossTableManager::GetCSDARange(
325 const G4ParticleDefinition *aParticle,
326 G4double kineticEnergy,
327 const G4MaterialCutsCouple *couple)
328{
329 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
330 G4double x = DBL_MAX;
331 if(currentLoss) x = currentLoss->GetCSDARange(kineticEnergy, couple);
332 return x;
333}
334
335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
336
337inline G4double G4LossTableManager::GetRangeFromRestricteDEDX(
338 const G4ParticleDefinition *aParticle,
339 G4double kineticEnergy,
340 const G4MaterialCutsCouple *couple)
341{
342 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
343 G4double x;
344 if(currentLoss) x = currentLoss->GetRangeForLoss(kineticEnergy, couple);
345 else
346 x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,couple,false);
347 return x;
348}
349
350//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
351
352inline G4double G4LossTableManager::GetRange(
353 const G4ParticleDefinition *aParticle,
354 G4double kineticEnergy,
355 const G4MaterialCutsCouple *couple)
356{
357 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
358 G4double x;
359 if(currentLoss) x = currentLoss->GetRange(kineticEnergy, couple);
360 else
361 x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,couple,false);
362 return x;
363}
364
365//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
366
367inline G4double G4LossTableManager::GetEnergy(
368 const G4ParticleDefinition *aParticle,
369 G4double range,
370 const G4MaterialCutsCouple *couple)
371{
372 if(aParticle != currentParticle) GetEnergyLossProcess(aParticle);
373 G4double x;
374 if(currentLoss) x = currentLoss->GetKineticEnergy(range, couple);
375 else x = G4EnergyLossTables::GetPreciseEnergyFromRange(
376 currentParticle,range,couple,false);
377 return x;
378}
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
381
382inline G4double G4LossTableManager::GetDEDXDispersion(
383 const G4MaterialCutsCouple *couple,
384 const G4DynamicParticle* dp,
385 G4double& length)
386{
387 const G4ParticleDefinition* aParticle = dp->GetDefinition();
388 if(aParticle != currentParticle) {
389 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
390 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
391 currentParticle = aParticle;
392 currentLoss = (*pos).second;
393 } else {
394 ParticleHaveNoLoss(aParticle);
395 }
396 }
397 return currentLoss->GetDEDXDispersion(couple, dp, length);
398}
399
400//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
401
402inline G4VEnergyLossProcess* G4LossTableManager::GetEnergyLossProcess(
403 const G4ParticleDefinition *aParticle)
404{
405 if(aParticle != currentParticle) {
406 currentParticle = aParticle;
407 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
408 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
409 currentLoss = (*pos).second;
410 } else {
411 currentLoss = 0;
412 // ParticleHaveNoLoss(aParticle);
413 }
414 }
415 return currentLoss;
416}
417
418//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
419
420inline G4EmCorrections* G4LossTableManager::EmCorrections()
421{
422 return emCorrections;
423}
424
425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
426
427#endif
428
Note: See TracBrowser for help on using the repository browser.