source: trunk/source/processes/electromagnetic/utils/src/G4EnergyLossTables.cc@ 991

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

update processes

File size: 30.7 KB
RevLine 
[819]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//
[961]27// $Id: G4EnergyLossTables.cc,v 1.34 2008/07/08 10:57:22 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
[819]29//
30// -------------------------------------------------------------------
31// first version created by P.Urban , 06/04/1998
32// modifications + "precise" functions added by L.Urban , 27/05/98
33// modifications , TOF functions , 26/10/98, L.Urban
34// cache mechanism in order to gain time, 11/02/99, L.Urban
35// bug fixed , 12/04/99 , L.Urban
36// 10.11.99: moved from RWT hash dictionary to STL map, G.Barrand, M.Maire
37// 27.09.01 L.Urban , bug fixed (negative energy deposit)
38// 26.10.01 all static functions moved from .icc files (mma)
39// 15.01.03 Add interfaces required for "cut per region" (V.Ivanchenko)
40// 12.03.03 Add warnings to obsolete interfaces (V.Ivanchenko)
41// 10.04.03 Add call to G4LossTableManager is particle is not registered (V.Ivanchenko)
42//
43// -------------------------------------------------------------------
44
45#include "G4EnergyLossTables.hh"
46#include "G4MaterialCutsCouple.hh"
47#include "G4RegionStore.hh"
48#include "G4LossTableManager.hh"
49
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
53G4EnergyLossTablesHelper G4EnergyLossTables::t ;
54G4EnergyLossTablesHelper G4EnergyLossTables::null_loss ;
55const G4ParticleDefinition* G4EnergyLossTables::lastParticle = 0;
56G4double G4EnergyLossTables::QQPositron = eplus*eplus ;
57G4double G4EnergyLossTables::Chargesquare ;
58G4int G4EnergyLossTables::oldIndex = -1 ;
59G4double G4EnergyLossTables::rmin = 0. ;
60G4double G4EnergyLossTables::rmax = 0. ;
61G4double G4EnergyLossTables::Thigh = 0. ;
62G4int G4EnergyLossTables::let_counter = 0;
63G4int G4EnergyLossTables::let_max_num_warnings = 100;
64G4bool G4EnergyLossTables::first_loss = true;
65
66G4EnergyLossTables::helper_map G4EnergyLossTables::dict;
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
70G4EnergyLossTablesHelper::G4EnergyLossTablesHelper(
71 const G4PhysicsTable* aDEDXTable,
72 const G4PhysicsTable* aRangeTable,
73 const G4PhysicsTable* anInverseRangeTable,
74 const G4PhysicsTable* aLabTimeTable,
75 const G4PhysicsTable* aProperTimeTable,
76 G4double aLowestKineticEnergy,
77 G4double aHighestKineticEnergy,
78 G4double aMassRatio,
79 G4int aNumberOfBins)
80 :
81 theDEDXTable(aDEDXTable), theRangeTable(aRangeTable),
82 theInverseRangeTable(anInverseRangeTable),
83 theLabTimeTable(aLabTimeTable),
84 theProperTimeTable(aProperTimeTable),
85 theLowestKineticEnergy(aLowestKineticEnergy),
86 theHighestKineticEnergy(aHighestKineticEnergy),
87 theMassRatio(aMassRatio),
88 theNumberOfBins(aNumberOfBins)
89{ }
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
93G4EnergyLossTablesHelper::G4EnergyLossTablesHelper()
94{ }
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
98void G4EnergyLossTables::Register(
99 const G4ParticleDefinition* p,
100 const G4PhysicsTable* tDEDX,
101 const G4PhysicsTable* tRange,
102 const G4PhysicsTable* tInverseRange,
103 const G4PhysicsTable* tLabTime,
104 const G4PhysicsTable* tProperTime,
105 G4double lowestKineticEnergy,
106 G4double highestKineticEnergy,
107 G4double massRatio,
108 G4int NumberOfBins)
109{
110 dict[p]= G4EnergyLossTablesHelper(tDEDX, tRange,tInverseRange,
111 tLabTime,tProperTime,lowestKineticEnergy,
112 highestKineticEnergy, massRatio,NumberOfBins);
113
114 t = GetTables(p) ; // important for cache !!!!!
115 lastParticle = p ;
116 Chargesquare = (p->GetPDGCharge())*(p->GetPDGCharge())/
117 QQPositron ;
118 if (first_loss ) {
119 null_loss = G4EnergyLossTablesHelper(0, 0, 0, 0, 0, 0.0, 0.0, 0.0, 0);
120 first_loss = false;
121 }
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
126const G4PhysicsTable* G4EnergyLossTables::GetDEDXTable(
127 const G4ParticleDefinition* p)
128{
129 helper_map::iterator it;
130 if((it=dict.find(p))==dict.end()) return 0;
131 return (*it).second.theDEDXTable;
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135
136const G4PhysicsTable* G4EnergyLossTables::GetRangeTable(
137 const G4ParticleDefinition* p)
138{
139 helper_map::iterator it;
140 if((it=dict.find(p))==dict.end()) return 0;
141 return (*it).second.theRangeTable;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
146const G4PhysicsTable* G4EnergyLossTables::GetInverseRangeTable(
147 const G4ParticleDefinition* p)
148{
149 helper_map::iterator it;
150 if((it=dict.find(p))==dict.end()) return 0;
151 return (*it).second.theInverseRangeTable;
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155
156const G4PhysicsTable* G4EnergyLossTables::GetLabTimeTable(
157 const G4ParticleDefinition* p)
158{
159 helper_map::iterator it;
160 if((it=dict.find(p))==dict.end()) return 0;
161 return (*it).second.theLabTimeTable;
162}
163
164//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165
166const G4PhysicsTable* G4EnergyLossTables::GetProperTimeTable(
167 const G4ParticleDefinition* p)
168{
169 helper_map::iterator it;
170 if((it=dict.find(p))==dict.end()) return 0;
171 return (*it).second.theProperTimeTable;
172}
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175
176G4EnergyLossTablesHelper G4EnergyLossTables::GetTables(
177 const G4ParticleDefinition* p)
178{
179 helper_map::iterator it;
180 if ((it=dict.find(p))==dict.end()) {
181// G4cout << "Table is not found out for " << p->GetParticleName() << G4endl;
182// G4Exception("G4EnergyLossTables::GetTables: table not found!");
183// exit(1);
184 return null_loss;
185 }
186 return (*it).second;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
191G4double G4EnergyLossTables::GetDEDX(
192 const G4ParticleDefinition *aParticle,
193 G4double KineticEnergy,
194 const G4Material *aMaterial)
195{
196 CPRWarning();
197 if(aParticle != lastParticle)
198 {
199 t= GetTables(aParticle);
200 lastParticle = aParticle ;
201 Chargesquare = (aParticle->GetPDGCharge())*
202 (aParticle->GetPDGCharge())/
203 QQPositron ;
204 oldIndex = -1 ;
205 }
206 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
207 if (!dEdxTable) ParticleHaveNoLoss(aParticle,"dEdx");
208
209 G4int materialIndex = aMaterial->GetIndex();
210 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
211 G4double dEdx;
212 G4bool isOut;
213
214 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
215
216 dEdx =(*dEdxTable)(materialIndex)->GetValue(
217 t.theLowestKineticEnergy,isOut)
218 *std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy);
219
220 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
221
222 dEdx = (*dEdxTable)(materialIndex)->GetValue(
223 t.theHighestKineticEnergy,isOut);
224
225 } else {
226
227 dEdx = (*dEdxTable)(materialIndex)->GetValue(
228 scaledKineticEnergy,isOut);
229
230 }
231
232 return dEdx*Chargesquare;
233}
234
235//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
236
237G4double G4EnergyLossTables::GetLabTime(
238 const G4ParticleDefinition *aParticle,
239 G4double KineticEnergy,
240 const G4Material *aMaterial)
241{
242 CPRWarning();
243 if(aParticle != lastParticle)
244 {
245 t= GetTables(aParticle);
246 lastParticle = aParticle ;
247 oldIndex = -1 ;
248 }
249 const G4PhysicsTable* labtimeTable= t.theLabTimeTable;
250 if (!labtimeTable) ParticleHaveNoLoss(aParticle,"LabTime");
251
252 const G4double parlowen=0.4 , ppar=0.5-parlowen ;
253 G4int materialIndex = aMaterial->GetIndex();
254 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
255 G4double time;
256 G4bool isOut;
257
258 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
259
260 time = std::exp(ppar*std::log(scaledKineticEnergy/t.theLowestKineticEnergy))*
261 (*labtimeTable)(materialIndex)->GetValue(
262 t.theLowestKineticEnergy,isOut);
263
264
265 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
266
267 time = (*labtimeTable)(materialIndex)->GetValue(
268 t.theHighestKineticEnergy,isOut);
269
270 } else {
271
272 time = (*labtimeTable)(materialIndex)->GetValue(
273 scaledKineticEnergy,isOut);
274
275 }
276
277 return time/t.theMassRatio ;
278}
279
280//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
281
282G4double G4EnergyLossTables::GetDeltaLabTime(
283 const G4ParticleDefinition *aParticle,
284 G4double KineticEnergyStart,
285 G4double KineticEnergyEnd,
286 const G4Material *aMaterial)
287{
288 CPRWarning();
289 if(aParticle != lastParticle)
290 {
291 t= GetTables(aParticle);
292 lastParticle = aParticle ;
293 oldIndex = -1 ;
294 }
295 const G4PhysicsTable* labtimeTable= t.theLabTimeTable;
296 if (!labtimeTable) ParticleHaveNoLoss(aParticle,"LabTime");
297
298 const G4double parlowen=0.4 , ppar=0.5-parlowen ;
299 const G4double dToverT = 0.05 , facT = 1. -dToverT ;
300 G4double timestart,timeend,deltatime,dTT;
301 G4bool isOut;
302
303 G4int materialIndex = aMaterial->GetIndex();
304 G4double scaledKineticEnergy = KineticEnergyStart*t.theMassRatio;
305
306 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
307
308 timestart = std::exp(ppar*std::log(scaledKineticEnergy/t.theLowestKineticEnergy))*
309 (*labtimeTable)(materialIndex)->GetValue(
310 t.theLowestKineticEnergy,isOut);
311
312
313 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
314
315 timestart = (*labtimeTable)(materialIndex)->GetValue(
316 t.theHighestKineticEnergy,isOut);
317
318 } else {
319
320 timestart = (*labtimeTable)(materialIndex)->GetValue(
321 scaledKineticEnergy,isOut);
322
323 }
324
325 dTT = (KineticEnergyStart - KineticEnergyEnd)/KineticEnergyStart ;
326
327 if( dTT < dToverT )
328 scaledKineticEnergy = facT*KineticEnergyStart*t.theMassRatio;
329 else
330 scaledKineticEnergy = KineticEnergyEnd*t.theMassRatio;
331
332 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
333
334 timeend = std::exp(ppar*std::log(scaledKineticEnergy/t.theLowestKineticEnergy))*
335 (*labtimeTable)(materialIndex)->GetValue(
336 t.theLowestKineticEnergy,isOut);
337
338
339 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
340
341 timeend = (*labtimeTable)(materialIndex)->GetValue(
342 t.theHighestKineticEnergy,isOut);
343
344 } else {
345
346 timeend = (*labtimeTable)(materialIndex)->GetValue(
347 scaledKineticEnergy,isOut);
348
349 }
350
351 deltatime = timestart - timeend ;
352
353 if( dTT < dToverT )
354 deltatime *= dTT/dToverT;
355
356 return deltatime/t.theMassRatio ;
357}
358
359//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
360
361G4double G4EnergyLossTables::GetProperTime(
362 const G4ParticleDefinition *aParticle,
363 G4double KineticEnergy,
364 const G4Material *aMaterial)
365{
366 CPRWarning();
367 if(aParticle != lastParticle)
368 {
369 t= GetTables(aParticle);
370 lastParticle = aParticle ;
371 oldIndex = -1 ;
372 }
373 const G4PhysicsTable* propertimeTable= t.theProperTimeTable;
374 if (!propertimeTable) ParticleHaveNoLoss(aParticle,"ProperTime");
375
376 const G4double parlowen=0.4 , ppar=0.5-parlowen ;
377 G4int materialIndex = aMaterial->GetIndex();
378 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
379 G4double time;
380 G4bool isOut;
381
382 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
383
384 time = std::exp(ppar*std::log(scaledKineticEnergy/t.theLowestKineticEnergy))*
385 (*propertimeTable)(materialIndex)->GetValue(
386 t.theLowestKineticEnergy,isOut);
387
388
389 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
390
391 time = (*propertimeTable)(materialIndex)->GetValue(
392 t.theHighestKineticEnergy,isOut);
393
394 } else {
395
396 time = (*propertimeTable)(materialIndex)->GetValue(
397 scaledKineticEnergy,isOut);
398
399 }
400
401 return time/t.theMassRatio ;
402}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
405
406G4double G4EnergyLossTables::GetDeltaProperTime(
407 const G4ParticleDefinition *aParticle,
408 G4double KineticEnergyStart,
409 G4double KineticEnergyEnd,
410 const G4Material *aMaterial)
411{
412 CPRWarning();
413 if(aParticle != lastParticle)
414 {
415 t= GetTables(aParticle);
416 lastParticle = aParticle ;
417 oldIndex = -1 ;
418 }
419 const G4PhysicsTable* propertimeTable= t.theProperTimeTable;
420 if (!propertimeTable) ParticleHaveNoLoss(aParticle,"ProperTime");
421
422 const G4double parlowen=0.4 , ppar=0.5-parlowen ;
423 const G4double dToverT = 0.05 , facT = 1. -dToverT ;
424 G4double timestart,timeend,deltatime,dTT;
425 G4bool isOut;
426
427 G4int materialIndex = aMaterial->GetIndex();
428 G4double scaledKineticEnergy = KineticEnergyStart*t.theMassRatio;
429
430 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
431
432 timestart = std::exp(ppar*std::log(scaledKineticEnergy/t.theLowestKineticEnergy))*
433 (*propertimeTable)(materialIndex)->GetValue(
434 t.theLowestKineticEnergy,isOut);
435
436
437 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
438
439 timestart = (*propertimeTable)(materialIndex)->GetValue(
440 t.theHighestKineticEnergy,isOut);
441
442 } else {
443
444 timestart = (*propertimeTable)(materialIndex)->GetValue(
445 scaledKineticEnergy,isOut);
446
447 }
448
449 dTT = (KineticEnergyStart - KineticEnergyEnd)/KineticEnergyStart ;
450
451 if( dTT < dToverT )
452 scaledKineticEnergy = facT*KineticEnergyStart*t.theMassRatio;
453 else
454 scaledKineticEnergy = KineticEnergyEnd*t.theMassRatio;
455
456 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
457
458 timeend = std::exp(ppar*std::log(scaledKineticEnergy/t.theLowestKineticEnergy))*
459 (*propertimeTable)(materialIndex)->GetValue(
460 t.theLowestKineticEnergy,isOut);
461
462
463 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
464
465 timeend = (*propertimeTable)(materialIndex)->GetValue(
466 t.theHighestKineticEnergy,isOut);
467
468 } else {
469
470 timeend = (*propertimeTable)(materialIndex)->GetValue(
471 scaledKineticEnergy,isOut);
472
473 }
474
475 deltatime = timestart - timeend ;
476
477 if( dTT < dToverT )
478 deltatime *= dTT/dToverT ;
479
480 return deltatime/t.theMassRatio ;
481}
482
483//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
484
485G4double G4EnergyLossTables::GetRange(
486 const G4ParticleDefinition *aParticle,
487 G4double KineticEnergy,
488 const G4Material *aMaterial)
489{
490 CPRWarning();
491 if(aParticle != lastParticle)
492 {
493 t= GetTables(aParticle);
494 lastParticle = aParticle ;
495 Chargesquare = (aParticle->GetPDGCharge())*
496 (aParticle->GetPDGCharge())/
497 QQPositron ;
498 oldIndex = -1 ;
499 }
500 const G4PhysicsTable* rangeTable= t.theRangeTable;
501 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
502 if (!rangeTable) ParticleHaveNoLoss(aParticle,"Range");
503
504 G4int materialIndex = aMaterial->GetIndex();
505 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
506 G4double Range;
507 G4bool isOut;
508
509 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
510
511 Range = std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy)*
512 (*rangeTable)(materialIndex)->GetValue(
513 t.theLowestKineticEnergy,isOut);
514
515 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
516
517 Range = (*rangeTable)(materialIndex)->GetValue(
518 t.theHighestKineticEnergy,isOut)+
519 (scaledKineticEnergy-t.theHighestKineticEnergy)/
520 (*dEdxTable)(materialIndex)->GetValue(
521 t.theHighestKineticEnergy,isOut);
522
523 } else {
524
525 Range = (*rangeTable)(materialIndex)->GetValue(
526 scaledKineticEnergy,isOut);
527
528 }
529
530 return Range/(Chargesquare*t.theMassRatio);
531}
532
533//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
534
535G4double G4EnergyLossTables::GetPreciseEnergyFromRange(
536 const G4ParticleDefinition *aParticle,
537 G4double range,
538 const G4Material *aMaterial)
539// it returns the value of the kinetic energy for a given range
540{
541 CPRWarning();
542 if( aParticle != lastParticle)
543 {
544 t= GetTables(aParticle);
545 lastParticle = aParticle;
546 Chargesquare = (aParticle->GetPDGCharge())*
547 (aParticle->GetPDGCharge())/
548 QQPositron ;
549 oldIndex = -1 ;
550 }
551 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
552 const G4PhysicsTable* inverseRangeTable= t.theInverseRangeTable;
553 if (!inverseRangeTable) ParticleHaveNoLoss(aParticle,"InverseRange");
554
555 G4double scaledrange,scaledKineticEnergy ;
556 G4bool isOut ;
557
558 G4int materialIndex = aMaterial->GetIndex() ;
559
560 if(materialIndex != oldIndex)
561 {
562 oldIndex = materialIndex ;
563 rmin = (*inverseRangeTable)(materialIndex)->
564 GetLowEdgeEnergy(0) ;
565 rmax = (*inverseRangeTable)(materialIndex)->
566 GetLowEdgeEnergy(t.theNumberOfBins-2) ;
567 Thigh = (*inverseRangeTable)(materialIndex)->
568 GetValue(rmax,isOut) ;
569 }
570
571 scaledrange = range*Chargesquare*t.theMassRatio ;
572
573 if(scaledrange < rmin)
574 {
575 scaledKineticEnergy = t.theLowestKineticEnergy*
576 scaledrange*scaledrange/(rmin*rmin) ;
577 }
578 else
579 {
580 if(scaledrange < rmax)
581 {
582 scaledKineticEnergy = (*inverseRangeTable)(materialIndex)->
583 GetValue( scaledrange,isOut) ;
584 }
585 else
586 {
587 scaledKineticEnergy = Thigh +
588 (scaledrange-rmax)*
589 (*dEdxTable)(materialIndex)->
590 GetValue(Thigh,isOut) ;
591 }
592 }
593
594 return scaledKineticEnergy/t.theMassRatio ;
595}
596
597//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
598
599 G4double G4EnergyLossTables::GetPreciseDEDX(
600 const G4ParticleDefinition *aParticle,
601 G4double KineticEnergy,
602 const G4Material *aMaterial)
603{
604 CPRWarning();
605 if( aParticle != lastParticle)
606 {
607 t= GetTables(aParticle);
608 lastParticle = aParticle;
609 Chargesquare = (aParticle->GetPDGCharge())*
610 (aParticle->GetPDGCharge())/
611 QQPositron ;
612 oldIndex = -1 ;
613 }
614 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
615 if (!dEdxTable) ParticleHaveNoLoss(aParticle,"dEdx");
616
617 G4int materialIndex = aMaterial->GetIndex();
618 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
619 G4double dEdx;
620 G4bool isOut;
621
622 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
623
624 dEdx = std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy)
625 *(*dEdxTable)(materialIndex)->GetValue(
626 t.theLowestKineticEnergy,isOut);
627
628 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
629
630 dEdx = (*dEdxTable)(materialIndex)->GetValue(
631 t.theHighestKineticEnergy,isOut);
632
633 } else {
634
635 dEdx = (*dEdxTable)(materialIndex)->GetValue(
636 scaledKineticEnergy,isOut) ;
637
638 }
639
640 return dEdx*Chargesquare;
641}
642
643//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
644
645 G4double G4EnergyLossTables::GetPreciseRangeFromEnergy(
646 const G4ParticleDefinition *aParticle,
647 G4double KineticEnergy,
648 const G4Material *aMaterial)
649{
650 CPRWarning();
651 if( aParticle != lastParticle)
652 {
653 t= GetTables(aParticle);
654 lastParticle = aParticle;
655 Chargesquare = (aParticle->GetPDGCharge())*
656 (aParticle->GetPDGCharge())/
657 QQPositron ;
658 oldIndex = -1 ;
659 }
660 const G4PhysicsTable* rangeTable= t.theRangeTable;
661 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
662 if (!rangeTable) ParticleHaveNoLoss(aParticle,"Range");
663
664 G4int materialIndex = aMaterial->GetIndex();
665
666 G4double Thighr = t.theHighestKineticEnergy*t.theLowestKineticEnergy/
667 (*rangeTable)(materialIndex)->
668 GetLowEdgeEnergy(1) ;
669
670 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
671 G4double Range;
672 G4bool isOut;
673
674 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
675
676 Range = std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy)*
677 (*rangeTable)(materialIndex)->GetValue(
678 t.theLowestKineticEnergy,isOut);
679
680 } else if (scaledKineticEnergy>Thighr) {
681
682 Range = (*rangeTable)(materialIndex)->GetValue(
683 Thighr,isOut)+
684 (scaledKineticEnergy-Thighr)/
685 (*dEdxTable)(materialIndex)->GetValue(
686 Thighr,isOut);
687
688 } else {
689
690 Range = (*rangeTable)(materialIndex)->GetValue(
691 scaledKineticEnergy,isOut) ;
692
693 }
694
695 return Range/(Chargesquare*t.theMassRatio);
696}
697
698//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
699
700G4double G4EnergyLossTables::GetDEDX(
701 const G4ParticleDefinition *aParticle,
702 G4double KineticEnergy,
703 const G4MaterialCutsCouple *couple,
704 G4bool check)
705{
706 if(aParticle != lastParticle)
707 {
708 t= GetTables(aParticle);
709 lastParticle = aParticle ;
710 Chargesquare = (aParticle->GetPDGCharge())*
711 (aParticle->GetPDGCharge())/
712 QQPositron ;
713 oldIndex = -1 ;
714 }
715 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
716
717 if (!dEdxTable ) {
718 if (check) return G4LossTableManager::Instance()->GetDEDX(aParticle,KineticEnergy,couple);
719 else ParticleHaveNoLoss(aParticle, "dEdx");
720 }
721
722 G4int materialIndex = couple->GetIndex();
723 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
724 G4double dEdx;
725 G4bool isOut;
726
727 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
728
729 dEdx =(*dEdxTable)(materialIndex)->GetValue(
730 t.theLowestKineticEnergy,isOut)
731 *std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy);
732
733 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
734
735 dEdx = (*dEdxTable)(materialIndex)->GetValue(
736 t.theHighestKineticEnergy,isOut);
737
738 } else {
739
740 dEdx = (*dEdxTable)(materialIndex)->GetValue(
741 scaledKineticEnergy,isOut);
742
743 }
744
745 return dEdx*Chargesquare;
746}
747
748//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
749
750G4double G4EnergyLossTables::GetRange(
751 const G4ParticleDefinition *aParticle,
752 G4double KineticEnergy,
753 const G4MaterialCutsCouple *couple,
754 G4bool check)
755{
756 if(aParticle != lastParticle)
757 {
758 t= GetTables(aParticle);
759 lastParticle = aParticle ;
760 Chargesquare = (aParticle->GetPDGCharge())*
761 (aParticle->GetPDGCharge())/
762 QQPositron ;
763 oldIndex = -1 ;
764 }
765 const G4PhysicsTable* rangeTable= t.theRangeTable;
766 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
767 if (!rangeTable) {
768 if(check) return G4LossTableManager::Instance()->GetRange(aParticle,KineticEnergy,couple);
769 else return DBL_MAX;
770 //ParticleHaveNoLoss(aParticle,"Range");
771 }
772
773 G4int materialIndex = couple->GetIndex();
774 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
775 G4double Range;
776 G4bool isOut;
777
778 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
779
780 Range = std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy)*
781 (*rangeTable)(materialIndex)->GetValue(
782 t.theLowestKineticEnergy,isOut);
783
784 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
785
786 Range = (*rangeTable)(materialIndex)->GetValue(
787 t.theHighestKineticEnergy,isOut)+
788 (scaledKineticEnergy-t.theHighestKineticEnergy)/
789 (*dEdxTable)(materialIndex)->GetValue(
790 t.theHighestKineticEnergy,isOut);
791
792 } else {
793
794 Range = (*rangeTable)(materialIndex)->GetValue(
795 scaledKineticEnergy,isOut);
796
797 }
798
799 return Range/(Chargesquare*t.theMassRatio);
800}
801
802//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
803
804G4double G4EnergyLossTables::GetPreciseEnergyFromRange(
805 const G4ParticleDefinition *aParticle,
806 G4double range,
807 const G4MaterialCutsCouple *couple,
808 G4bool check)
809// it returns the value of the kinetic energy for a given range
810{
811 if( aParticle != lastParticle)
812 {
813 t= GetTables(aParticle);
814 lastParticle = aParticle;
815 Chargesquare = (aParticle->GetPDGCharge())*
816 (aParticle->GetPDGCharge())/
817 QQPositron ;
818 oldIndex = -1 ;
819 }
820 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
821 const G4PhysicsTable* inverseRangeTable= t.theInverseRangeTable;
822
823 if (!inverseRangeTable) {
824 if(check) return G4LossTableManager::Instance()->GetEnergy(aParticle,range,couple);
825 else return DBL_MAX;
826 // else ParticleHaveNoLoss(aParticle,"InverseRange");
827 }
828
829 G4double scaledrange,scaledKineticEnergy ;
830 G4bool isOut ;
831
832 G4int materialIndex = couple->GetIndex() ;
833
834 if(materialIndex != oldIndex)
835 {
836 oldIndex = materialIndex ;
837 rmin = (*inverseRangeTable)(materialIndex)->
838 GetLowEdgeEnergy(0) ;
839 rmax = (*inverseRangeTable)(materialIndex)->
840 GetLowEdgeEnergy(t.theNumberOfBins-2) ;
841 Thigh = (*inverseRangeTable)(materialIndex)->
842 GetValue(rmax,isOut) ;
843 }
844
845 scaledrange = range*Chargesquare*t.theMassRatio ;
846
847 if(scaledrange < rmin)
848 {
849 scaledKineticEnergy = t.theLowestKineticEnergy*
850 scaledrange*scaledrange/(rmin*rmin) ;
851 }
852 else
853 {
854 if(scaledrange < rmax)
855 {
856 scaledKineticEnergy = (*inverseRangeTable)(materialIndex)->
857 GetValue( scaledrange,isOut) ;
858 }
859 else
860 {
861 scaledKineticEnergy = Thigh +
862 (scaledrange-rmax)*
863 (*dEdxTable)(materialIndex)->
864 GetValue(Thigh,isOut) ;
865 }
866 }
867
868 return scaledKineticEnergy/t.theMassRatio ;
869}
870
871//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
872
873G4double G4EnergyLossTables::GetPreciseDEDX(
874 const G4ParticleDefinition *aParticle,
875 G4double KineticEnergy,
876 const G4MaterialCutsCouple *couple)
877{
878
879 if( aParticle != lastParticle)
880 {
881 t= GetTables(aParticle);
882 lastParticle = aParticle;
883 Chargesquare = (aParticle->GetPDGCharge())*
884 (aParticle->GetPDGCharge())/
885 QQPositron ;
886 oldIndex = -1 ;
887 }
888 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
889 if ( !dEdxTable )
890 return G4LossTableManager::Instance()->GetDEDX(aParticle,KineticEnergy,couple);
891
892 G4int materialIndex = couple->GetIndex();
893 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
894 G4double dEdx;
895 G4bool isOut;
896
897 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
898
899 dEdx = std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy)
900 *(*dEdxTable)(materialIndex)->GetValue(
901 t.theLowestKineticEnergy,isOut);
902
903 } else if (scaledKineticEnergy>t.theHighestKineticEnergy) {
904
905 dEdx = (*dEdxTable)(materialIndex)->GetValue(
906 t.theHighestKineticEnergy,isOut);
907
908 } else {
909
910 dEdx = (*dEdxTable)(materialIndex)->GetValue(
911 scaledKineticEnergy,isOut) ;
912
913 }
914
915 return dEdx*Chargesquare;
916}
917
918//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
919
920G4double G4EnergyLossTables::GetPreciseRangeFromEnergy(
921 const G4ParticleDefinition *aParticle,
922 G4double KineticEnergy,
923 const G4MaterialCutsCouple *couple)
924{
925 if( aParticle != lastParticle)
926 {
927 t= GetTables(aParticle);
928 lastParticle = aParticle;
929 Chargesquare = (aParticle->GetPDGCharge())*
930 (aParticle->GetPDGCharge())/
931 QQPositron ;
932 oldIndex = -1 ;
933 }
934 const G4PhysicsTable* rangeTable= t.theRangeTable;
935 const G4PhysicsTable* dEdxTable= t.theDEDXTable;
936 if ( !dEdxTable || !rangeTable)
937 return G4LossTableManager::Instance()->GetDEDX(aParticle,KineticEnergy,couple);
938
939 G4int materialIndex = couple->GetIndex();
940
941 G4double Thighr = t.theHighestKineticEnergy*t.theLowestKineticEnergy/
942 (*rangeTable)(materialIndex)->
943 GetLowEdgeEnergy(1) ;
944
945 G4double scaledKineticEnergy = KineticEnergy*t.theMassRatio;
946 G4double Range;
947 G4bool isOut;
948
949 if (scaledKineticEnergy<t.theLowestKineticEnergy) {
950
951 Range = std::sqrt(scaledKineticEnergy/t.theLowestKineticEnergy)*
952 (*rangeTable)(materialIndex)->GetValue(
953 t.theLowestKineticEnergy,isOut);
954
955 } else if (scaledKineticEnergy>Thighr) {
956
957 Range = (*rangeTable)(materialIndex)->GetValue(
958 Thighr,isOut)+
959 (scaledKineticEnergy-Thighr)/
960 (*dEdxTable)(materialIndex)->GetValue(
961 Thighr,isOut);
962
963 } else {
964
965 Range = (*rangeTable)(materialIndex)->GetValue(
966 scaledKineticEnergy,isOut) ;
967
968 }
969
970 return Range/(Chargesquare*t.theMassRatio);
971}
972
973//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
974
975void G4EnergyLossTables::CPRWarning()
976{
977 if (let_counter < let_max_num_warnings) {
978
979 G4cout << G4endl;
980 G4cout << "##### G4EnergyLossTable WARNING: The obsolete interface is used!" << G4endl;
981 G4cout << "##### RESULTS ARE NOT GARANTEED!" << G4endl;
982 G4cout << "##### Please, substitute G4Material by G4MaterialCutsCouple" << G4endl;
983 G4cout << "##### Obsolete interface will be removed soon" << G4endl;
984 G4cout << G4endl;
985 let_counter++;
986// if ((G4RegionStore::GetInstance())->size() > 1) {
987// G4Exception("G4EnergyLossTables:: More than 1 region - table can't be accessed with obsolete interface");
988// exit(1);
989// }
990
991 } else if (let_counter == let_max_num_warnings) {
992
993 G4cout << "##### G4EnergyLossTable WARNING closed" << G4endl;
994 let_counter++;
995 }
996}
997
998//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
999
[961]1000void G4EnergyLossTables::ParticleHaveNoLoss(const G4ParticleDefinition* aParticle,
1001 const G4String& q)
[819]1002{
[961]1003 G4String s = " " + q + " table not found for "
1004 + aParticle->GetParticleName() + " !";
1005 G4Exception("G4EnergyLossTables::ParticleHaveNoLoss", "EM01",
1006 FatalException, s);
[819]1007}
1008
1009//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.