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.55 2009/10/29 19:25:28 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
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 | |
---|
85 | class G4PhysicsTable; |
---|
86 | class G4MaterialCutsCouple; |
---|
87 | class G4EnergyLossMessenger; |
---|
88 | class G4ParticleDefinition; |
---|
89 | class G4VMultipleScattering; |
---|
90 | class G4VEmProcess; |
---|
91 | class G4EmCorrections; |
---|
92 | class G4EmSaturation; |
---|
93 | class G4EmConfigurator; |
---|
94 | class G4LossTableBuilder; |
---|
95 | |
---|
96 | class G4LossTableManager |
---|
97 | { |
---|
98 | |
---|
99 | public: |
---|
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 SetFactorForAngleLimit(G4double val); |
---|
213 | |
---|
214 | void SetVerbose(G4int val); |
---|
215 | |
---|
216 | G4EnergyLossMessenger* GetMessenger(); |
---|
217 | |
---|
218 | G4bool BuildCSDARange() const; |
---|
219 | |
---|
220 | G4bool LPMFlag() const; |
---|
221 | |
---|
222 | G4bool SplineFlag() const; |
---|
223 | |
---|
224 | G4double BremsstrahlungTh() const; |
---|
225 | |
---|
226 | G4double FactorForAngleLimit() const; |
---|
227 | |
---|
228 | const std::vector<G4VEnergyLossProcess*>& GetEnergyLossProcessVector(); |
---|
229 | |
---|
230 | const std::vector<G4VEmProcess*>& GetEmProcessVector(); |
---|
231 | |
---|
232 | const std::vector<G4VMultipleScattering*>& GetMultipleScatteringVector(); |
---|
233 | |
---|
234 | inline G4VEnergyLossProcess* GetEnergyLossProcess(const G4ParticleDefinition*); |
---|
235 | |
---|
236 | G4EmCorrections* EmCorrections(); |
---|
237 | |
---|
238 | G4EmSaturation* EmSaturation(); |
---|
239 | |
---|
240 | G4EmConfigurator* EmConfigurator(); |
---|
241 | |
---|
242 | private: |
---|
243 | |
---|
244 | G4LossTableManager(); |
---|
245 | |
---|
246 | G4VEnergyLossProcess* BuildTables(const G4ParticleDefinition* aParticle); |
---|
247 | |
---|
248 | void CopyTables(const G4ParticleDefinition* aParticle, |
---|
249 | G4VEnergyLossProcess*); |
---|
250 | |
---|
251 | void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle); |
---|
252 | |
---|
253 | void SetParameters(G4VEnergyLossProcess*); |
---|
254 | |
---|
255 | void CopyDEDXTables(); |
---|
256 | |
---|
257 | private: |
---|
258 | |
---|
259 | static G4LossTableManager* theInstance; |
---|
260 | |
---|
261 | typedef const G4ParticleDefinition* PD; |
---|
262 | std::map<PD,G4VEnergyLossProcess*,std::less<PD> > loss_map; |
---|
263 | |
---|
264 | std::vector<G4VEnergyLossProcess*> loss_vector; |
---|
265 | std::vector<PD> part_vector; |
---|
266 | std::vector<PD> base_part_vector; |
---|
267 | std::vector<G4bool> tables_are_built; |
---|
268 | std::vector<G4bool> isActive; |
---|
269 | std::vector<G4PhysicsTable*> dedx_vector; |
---|
270 | std::vector<G4PhysicsTable*> range_vector; |
---|
271 | std::vector<G4PhysicsTable*> inv_range_vector; |
---|
272 | std::vector<G4VMultipleScattering*> msc_vector; |
---|
273 | std::vector<G4VEmProcess*> emp_vector; |
---|
274 | std::vector<G4VEmModel*> mod_vector; |
---|
275 | std::vector<G4VEmFluctuationModel*> fmod_vector; |
---|
276 | |
---|
277 | // cash |
---|
278 | G4VEnergyLossProcess* currentLoss; |
---|
279 | PD currentParticle; |
---|
280 | PD theElectron; |
---|
281 | |
---|
282 | G4int n_loss; |
---|
283 | G4int run; |
---|
284 | |
---|
285 | G4bool all_tables_are_built; |
---|
286 | // G4bool first_entry; |
---|
287 | G4bool lossFluctuationFlag; |
---|
288 | G4bool subCutoffFlag; |
---|
289 | G4bool rndmStepFlag; |
---|
290 | G4bool integral; |
---|
291 | G4bool integralActive; |
---|
292 | G4bool all_tables_are_stored; |
---|
293 | G4bool buildCSDARange; |
---|
294 | G4bool minEnergyActive; |
---|
295 | G4bool maxEnergyActive; |
---|
296 | G4bool maxEnergyForMuonsActive; |
---|
297 | G4bool stepFunctionActive; |
---|
298 | G4bool flagLPM; |
---|
299 | G4bool splineFlag; |
---|
300 | |
---|
301 | G4double minSubRange; |
---|
302 | G4double maxRangeVariation; |
---|
303 | G4double maxFinalStep; |
---|
304 | G4double minKinEnergy; |
---|
305 | G4double maxKinEnergy; |
---|
306 | G4double maxKinEnergyForMuons; |
---|
307 | G4double bremsTh; |
---|
308 | G4double factorForAngleLimit; |
---|
309 | |
---|
310 | G4LossTableBuilder* tableBuilder; |
---|
311 | G4EnergyLossMessenger* theMessenger; |
---|
312 | G4EmCorrections* emCorrections; |
---|
313 | G4EmSaturation* emSaturation; |
---|
314 | G4EmConfigurator* emConfigurator; |
---|
315 | |
---|
316 | const G4ParticleDefinition* firstParticle; |
---|
317 | G4int verbose; |
---|
318 | |
---|
319 | }; |
---|
320 | |
---|
321 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
322 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... |
---|
323 | |
---|
324 | inline G4double G4LossTableManager::GetDEDX( |
---|
325 | const G4ParticleDefinition *aParticle, |
---|
326 | G4double kineticEnergy, |
---|
327 | const G4MaterialCutsCouple *couple) |
---|
328 | { |
---|
329 | if(aParticle != currentParticle) GetEnergyLossProcess(aParticle); |
---|
330 | G4double x; |
---|
331 | if(currentLoss) x = currentLoss->GetDEDX(kineticEnergy, couple); |
---|
332 | else x = G4EnergyLossTables::GetDEDX( |
---|
333 | currentParticle,kineticEnergy,couple,false); |
---|
334 | return x; |
---|
335 | } |
---|
336 | |
---|
337 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... |
---|
338 | |
---|
339 | inline G4double G4LossTableManager::GetSubDEDX( |
---|
340 | const G4ParticleDefinition *aParticle, |
---|
341 | G4double kineticEnergy, |
---|
342 | const G4MaterialCutsCouple *couple) |
---|
343 | { |
---|
344 | if(aParticle != currentParticle) GetEnergyLossProcess(aParticle); |
---|
345 | G4double x = 0.0; |
---|
346 | if(currentLoss) x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple); |
---|
347 | return x; |
---|
348 | } |
---|
349 | |
---|
350 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... |
---|
351 | |
---|
352 | inline G4double G4LossTableManager::GetCSDARange( |
---|
353 | const G4ParticleDefinition *aParticle, |
---|
354 | G4double kineticEnergy, |
---|
355 | const G4MaterialCutsCouple *couple) |
---|
356 | { |
---|
357 | if(aParticle != currentParticle) GetEnergyLossProcess(aParticle); |
---|
358 | G4double x = DBL_MAX; |
---|
359 | if(currentLoss) x = currentLoss->GetCSDARange(kineticEnergy, couple); |
---|
360 | return x; |
---|
361 | } |
---|
362 | |
---|
363 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... |
---|
364 | |
---|
365 | inline G4double G4LossTableManager::GetRangeFromRestricteDEDX( |
---|
366 | const G4ParticleDefinition *aParticle, |
---|
367 | G4double kineticEnergy, |
---|
368 | const G4MaterialCutsCouple *couple) |
---|
369 | { |
---|
370 | if(aParticle != currentParticle) GetEnergyLossProcess(aParticle); |
---|
371 | G4double x; |
---|
372 | if(currentLoss) x = currentLoss->GetRangeForLoss(kineticEnergy, couple); |
---|
373 | else |
---|
374 | x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,couple,false); |
---|
375 | return x; |
---|
376 | } |
---|
377 | |
---|
378 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... |
---|
379 | |
---|
380 | inline G4double G4LossTableManager::GetRange( |
---|
381 | const G4ParticleDefinition *aParticle, |
---|
382 | G4double kineticEnergy, |
---|
383 | const G4MaterialCutsCouple *couple) |
---|
384 | { |
---|
385 | if(aParticle != currentParticle) GetEnergyLossProcess(aParticle); |
---|
386 | G4double x; |
---|
387 | if(currentLoss) x = currentLoss->GetRange(kineticEnergy, couple); |
---|
388 | else |
---|
389 | x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,couple,false); |
---|
390 | return x; |
---|
391 | } |
---|
392 | |
---|
393 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
394 | |
---|
395 | inline G4double G4LossTableManager::GetEnergy( |
---|
396 | const G4ParticleDefinition *aParticle, |
---|
397 | G4double range, |
---|
398 | const G4MaterialCutsCouple *couple) |
---|
399 | { |
---|
400 | if(aParticle != currentParticle) GetEnergyLossProcess(aParticle); |
---|
401 | G4double x; |
---|
402 | if(currentLoss) x = currentLoss->GetKineticEnergy(range, couple); |
---|
403 | else x = G4EnergyLossTables::GetPreciseEnergyFromRange( |
---|
404 | currentParticle,range,couple,false); |
---|
405 | return x; |
---|
406 | } |
---|
407 | |
---|
408 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
409 | |
---|
410 | inline G4double G4LossTableManager::GetDEDXDispersion( |
---|
411 | const G4MaterialCutsCouple *couple, |
---|
412 | const G4DynamicParticle* dp, |
---|
413 | G4double& length) |
---|
414 | { |
---|
415 | const G4ParticleDefinition* aParticle = dp->GetDefinition(); |
---|
416 | if(aParticle != currentParticle) { |
---|
417 | std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos; |
---|
418 | if ((pos = loss_map.find(aParticle)) != loss_map.end()) { |
---|
419 | currentParticle = aParticle; |
---|
420 | currentLoss = (*pos).second; |
---|
421 | } else { |
---|
422 | ParticleHaveNoLoss(aParticle); |
---|
423 | } |
---|
424 | } |
---|
425 | return currentLoss->GetDEDXDispersion(couple, dp, length); |
---|
426 | } |
---|
427 | |
---|
428 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
429 | |
---|
430 | inline G4VEnergyLossProcess* G4LossTableManager::GetEnergyLossProcess( |
---|
431 | const G4ParticleDefinition *aParticle) |
---|
432 | { |
---|
433 | if(aParticle != currentParticle) { |
---|
434 | currentParticle = aParticle; |
---|
435 | std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos; |
---|
436 | if ((pos = loss_map.find(aParticle)) != loss_map.end()) { |
---|
437 | currentLoss = (*pos).second; |
---|
438 | } else { |
---|
439 | currentLoss = 0; |
---|
440 | // ParticleHaveNoLoss(aParticle); |
---|
441 | } |
---|
442 | } |
---|
443 | return currentLoss; |
---|
444 | } |
---|
445 | |
---|
446 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
447 | |
---|
448 | #endif |
---|
449 | |
---|