Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EnergyLossTables.hh
Go to the documentation of this file.
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//
27//
28//
29// -----------------------------------------------------------------------------
30
31#ifndef included_G4EnergyLossTables
32#define included_G4EnergyLossTables
33
34#include <map>
35#include "globals.hh"
36
37#include "G4PhysicsTable.hh"
39#include "G4Material.hh"
40#include "G4ios.hh"
41
42//------------------------------------------------------------------------------
43// A utility class, containing the energy loss tables
44// for each particle
45//
46// Energy loss processes have to register their tables with this
47// class. The responsibility of creating and deleting the tables
48// remains with the energy loss classes.
49// -----------------------------------------------------------------------------
50//
51// P. Urban, 06/04/1998
52// L. Urban, 27/05/1988 , modifications + new functions added
53// L.Urban , 13/10/98 , revision
54// L.Urban, 26/10/98 , revision, Interpolate removed
55// L.Urban , 08/02/99, cache mechanism
56// L.Urban , 12/04/99 , bug fixed
57// don't use the helper class.
58// It can't be hidden for Rogue Wave uses it.
59// 10.11.99: moved from RWT hash dictionary to STL map, G.Barrand, M.Maire
60// 26.10.01: all static functions movev from .icc to .cc file (mma)
61// 15.01.03 Add interfaces required for "cut per region" (V.Ivanchenko)
62// 12.03.03 Add warnings to obsolete interfaces (V.Ivanchenko)
63// 12.04.03 move exception to new method (V.Ivanchenko)
64//
65// -----------------------------------------------------------------------------
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
70
71friend class G4EnergyLossTables;
72 // the only instances are within the class G4EnergyLossTables
73
74public:
76
77private:
79 const G4PhysicsTable* aRangeTable,
80 const G4PhysicsTable* anInverseRangeTable,
81 const G4PhysicsTable* aLabTimeTable,
82 const G4PhysicsTable* aProperTimeTable,
83 G4double aLowestKineticEnergy,
84 G4double aHighestKineticEnergy,
85 G4double aMassRatio,
86 G4int aNumberOfBins);
87 // data to be stored in the dictionary
88 const G4PhysicsTable* theDEDXTable;
89 const G4PhysicsTable* theRangeTable;
90 const G4PhysicsTable* theInverseRangeTable;
91 const G4PhysicsTable* theLabTimeTable;
92 const G4PhysicsTable* theProperTimeTable;
93 G4double theLowestKineticEnergy;
94 G4double theHighestKineticEnergy;
95 G4double theMassRatio;
96 G4int theNumberOfBins;
97};
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
102
104
105public:
106
107 // get the table for a given particle
108 // (0 if the table was not found)
109 static const G4PhysicsTable* GetDEDXTable(
110 const G4ParticleDefinition* p);
111 static const G4PhysicsTable* GetRangeTable(
112 const G4ParticleDefinition* p);
114 const G4ParticleDefinition* p);
115 static const G4PhysicsTable* GetLabTimeTable(
116 const G4ParticleDefinition* p);
118 const G4ParticleDefinition* p);
119
120 // get the DEDX or the range for a given particle/energy/material
121 static G4double GetDEDX(
122 const G4ParticleDefinition *aParticle,
123 G4double KineticEnergy,
124 const G4Material *aMaterial);
125 static G4double GetRange(
126 const G4ParticleDefinition *aParticle,
127 G4double KineticEnergy,
128 const G4Material *aMaterial);
129 static G4double GetLabTime(
130 const G4ParticleDefinition *aParticle,
131 G4double KineticEnergy,
132 const G4Material *aMaterial);
134 const G4ParticleDefinition *aParticle,
135 G4double KineticEnergyStart,
136 G4double KineticEnergyEnd,
137 const G4Material *aMaterial);
138 static G4double GetProperTime(
139 const G4ParticleDefinition *aParticle,
140 G4double KineticEnergy,
141 const G4Material *aMaterial);
143 const G4ParticleDefinition *aParticle,
144 G4double KineticEnergyStart,
145 G4double KineticEnergyEnd,
146 const G4Material *aMaterial);
147
149 const G4ParticleDefinition *aParticle,
150 G4double KineticEnergy,
151 const G4Material *aMaterial);
153 const G4ParticleDefinition *aParticle,
154 G4double KineticEnergy,
155 const G4Material *aMaterial);
157 const G4ParticleDefinition *aParticle,
158 G4double range,
159 const G4Material *aMaterial);
160
161 // get the DEDX or the range for a given particle/energy/materialCutsCouple
162 static G4double GetDEDX(
163 const G4ParticleDefinition *aParticle,
164 G4double KineticEnergy,
165 const G4MaterialCutsCouple *couple,
166 G4bool check = true);
167 static G4double GetRange(
168 const G4ParticleDefinition *aParticle,
169 G4double KineticEnergy,
170 const G4MaterialCutsCouple *couple,
171 G4bool check = true);
172
174 const G4ParticleDefinition *aParticle,
175 G4double KineticEnergy,
176 const G4MaterialCutsCouple *couple);
178 const G4ParticleDefinition *aParticle,
179 G4double KineticEnergy,
180 const G4MaterialCutsCouple *couple);
182 const G4ParticleDefinition *aParticle,
183 G4double range,
184 const G4MaterialCutsCouple *couple,
185 G4bool check = true);
186
187 // to be called only by energy loss processes
188 static void Register(
189 const G4ParticleDefinition* p,
190 const G4PhysicsTable* tDEDX,
191 const G4PhysicsTable* tRange,
192 const G4PhysicsTable* tInverseRange,
193 const G4PhysicsTable* tLabTime,
194 const G4PhysicsTable* tProperTime,
195 G4double lowestKineticEnergy,
196 G4double highestKineticEnergy,
197 G4double massRatio,
198 G4int NumberOfBins);
199
200public:
201 typedef const G4ParticleDefinition* K;
202
203private:
204
205 static void CPRWarning();
206 static void ParticleHaveNoLoss(const G4ParticleDefinition* aParticle,
207 const G4String&);
208
209 /*
210 typedef std::map<K,G4EnergyLossTablesHelper,std::less<K> > helper_map;
211 static G4ThreadLocal helper_map *dict;
212
213 static G4EnergyLossTablesHelper GetTables(const G4ParticleDefinition* p);
214
215 static G4ThreadLocal G4EnergyLossTablesHelper *t;
216 static G4ThreadLocal G4EnergyLossTablesHelper *null_loss;
217 static G4ThreadLocal G4ParticleDefinition* lastParticle ;
218 static G4ThreadLocal G4double QQPositron ;
219 static G4ThreadLocal G4double Chargesquare ;
220 static G4ThreadLocal G4int oldIndex ;
221 static G4ThreadLocal G4double rmin,rmax,Thigh ;
222 static G4ThreadLocal G4int let_counter;
223 static G4ThreadLocal G4int let_max_num_warnings;
224 static G4ThreadLocal G4bool first_loss;
225 */
226 typedef std::map<K,G4EnergyLossTablesHelper,std::less<K> > helper_map;
227 static helper_map *dict;
228
229 static G4EnergyLossTablesHelper GetTables(const G4ParticleDefinition* p);
230
231 static G4EnergyLossTablesHelper *t;
232 static G4EnergyLossTablesHelper *null_loss;
233 static G4ParticleDefinition* lastParticle ;
234 static G4double QQPositron ;
235 static G4double Chargesquare ;
236 static G4int oldIndex ;
237 static G4double rmin,rmax,Thigh ;
238 static G4int let_counter;
239 static G4int let_max_num_warnings;
240 static G4bool first_loss;
241
242};
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245
246#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static G4double GetRange(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
static const G4PhysicsTable * GetLabTimeTable(const G4ParticleDefinition *p)
static void Register(const G4ParticleDefinition *p, const G4PhysicsTable *tDEDX, const G4PhysicsTable *tRange, const G4PhysicsTable *tInverseRange, const G4PhysicsTable *tLabTime, const G4PhysicsTable *tProperTime, G4double lowestKineticEnergy, G4double highestKineticEnergy, G4double massRatio, G4int NumberOfBins)
static G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
static G4double GetDeltaProperTime(const G4ParticleDefinition *aParticle, G4double KineticEnergyStart, G4double KineticEnergyEnd, const G4Material *aMaterial)
static G4double GetLabTime(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
static G4double GetPreciseDEDX(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
static G4double GetPreciseEnergyFromRange(const G4ParticleDefinition *aParticle, G4double range, const G4Material *aMaterial)
static const G4PhysicsTable * GetDEDXTable(const G4ParticleDefinition *p)
static const G4PhysicsTable * GetRangeTable(const G4ParticleDefinition *p)
static const G4PhysicsTable * GetInverseRangeTable(const G4ParticleDefinition *p)
static G4double GetDeltaLabTime(const G4ParticleDefinition *aParticle, G4double KineticEnergyStart, G4double KineticEnergyEnd, const G4Material *aMaterial)
static G4double GetPreciseRangeFromEnergy(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
static G4double GetProperTime(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
static const G4PhysicsTable * GetProperTimeTable(const G4ParticleDefinition *p)
const G4ParticleDefinition * K