Geant4-11
G4VRangeToEnergyConverter.cc
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// G4VRangeToEnergyConverter class implementation
27//
28// Author: H.Kurashige, 05 October 2002 - First implementation
29// --------------------------------------------------------------------
30
32#include "G4ParticleTable.hh"
33#include "G4Element.hh"
34#include "G4SystemOfUnits.hh"
35#include "G4Log.hh"
36#include "G4Exp.hh"
37
38#ifdef G4MULTITHREADED
39G4Mutex G4VRangeToEnergyConverter::theMutex = G4MUTEX_INITIALIZER;
40#endif
41
44
45std::vector<G4double>* G4VRangeToEnergyConverter::Energy = nullptr;
46
49
50// --------------------------------------------------------------------
52{
53 if(nullptr == Energy)
54 {
55#ifdef G4MULTITHREADED
56 G4MUTEXLOCK(&theMutex);
57 if(nullptr == Energy)
58 {
59#endif
60 isFirstInstance = true;
61 Energy = new std::vector<G4double>(Nbin + 1);
62#ifdef G4MULTITHREADED
63 }
64 G4MUTEXUNLOCK(&theMutex);
65#endif
66 }
67 // this method defines lock itself
69 {
71 }
72}
73
74// --------------------------------------------------------------------
76{
78 {
79 delete Energy;
80 Energy = nullptr;
81 Emin = 0.;
82 Emax = 0.;
83 }
84}
85
86// --------------------------------------------------------------------
88 const G4Material* material)
89{
90#ifdef G4VERBOSE
91 if (GetVerboseLevel()>3)
92 {
93 G4cout << "G4VRangeToEnergyConverter::Convert() - ";
94 G4cout << "Convert for " << material->GetName()
95 << " with Range Cut " << rangeCut/mm << "[mm]" << G4endl;
96 }
97#endif
98
99 G4double cut = 0.0;
100 if(fPDG == 22)
101 {
102 cut = ConvertForGamma(rangeCut, material);
103 }
104 else
105 {
106 cut = ConvertForElectron(rangeCut, material);
107
108 const G4double tune = 0.025*CLHEP::mm*CLHEP::g/CLHEP::cm3;
109 const G4double lowen = 30.*CLHEP::keV;
110 if(cut < lowen)
111 {
112 // corr. should be switched on smoothly
113 cut /= (1.+(1.-cut/lowen)*tune/(rangeCut*material->GetDensity()));
114 }
115 }
116
117 cut = std::max(Emin, std::min(cut, Emax));
118 return cut;
119}
120
121// --------------------------------------------------------------------
123 const G4double highedge)
124{
125 G4double ehigh = std::min(Emax, highedge);
126 if(ehigh > lowedge)
127 {
128 FillEnergyVector(lowedge, ehigh);
129 }
130}
131
132// --------------------------------------------------------------------
134{
135 return Emin;
136}
137
138// --------------------------------------------------------------------
140{
141 return Emax;
142}
143
144// --------------------------------------------------------------------
145
147{
148 return Emax;
149}
150
151// --------------------------------------------------------------------
153{
154 if(value > Emin)
155 {
156 FillEnergyVector(Emin, value);
157 }
158}
159
160// --------------------------------------------------------------------
162 const G4double emax)
163{
164 if(emin != Emin || emax != Emax)
165 {
166#ifdef G4MULTITHREADED
167 G4MUTEXLOCK(&theMutex);
168 if(emin != Emin || emax != Emax)
169 {
170#endif
171 Emin = emin;
172 Emax = emax;
173 Nbin = NbinPerDecade*static_cast<G4int>(std::log10(emax/emin));
174 Energy->resize(Nbin + 1);
175 (*Energy)[0] = emin;
176 (*Energy)[Nbin] = emax;
177 G4double fact = G4Log(emax/emin)/Nbin;
178 for(G4int i=1; i<Nbin; ++i) { (*Energy)[i] = emin*G4Exp(i * fact); }
179#ifdef G4MULTITHREADED
180 }
181 G4MUTEXUNLOCK(&theMutex);
182#endif
183 }
184}
185
186// --------------------------------------------------------------------
189 const G4Material* material)
190{
191 const G4ElementVector* elm = material->GetElementVector();
192 const G4double* dens = material->GetAtomicNumDensityVector();
193
194 // fill absorption length vector
195 G4int nelm = material->GetNumberOfElements();
196 G4double range1 = 0.0;
197 G4double range2 = 0.0;
198 G4double e1 = 0.0;
199 G4double e2 = 0.0;
200 for (G4int i=0; i<Nbin; ++i)
201 {
202 e2 = (*Energy)[i];
203 G4double sig = 0.;
204
205 for (G4int j=0; j<nelm; ++j)
206 {
207 sig += dens[j]*ComputeValue((*elm)[j]->GetZasInt(), e2);
208 }
209 range2 = (sig > 0.0) ? 5./sig : DBL_MAX;
210 if(i == 0 || range2 < rangeCut)
211 {
212 e1 = e2;
213 range1 = range2;
214 }
215 else
216 {
217 break;
218 }
219 }
220 return LiniearInterpolation(e1, e2, range1, range2, rangeCut);
221}
222
223// --------------------------------------------------------------------
226 const G4Material* material)
227{
228 const G4ElementVector* elm = material->GetElementVector();
229 const G4double* dens = material->GetAtomicNumDensityVector();
230
231 // fill absorption length vector
232 G4int nelm = material->GetNumberOfElements();
233 G4double dedx1 = 0.0;
234 G4double dedx2 = 0.0;
235 G4double range1 = 0.0;
236 G4double range2 = 0.0;
237 G4double e1 = 0.0;
238 G4double e2 = 0.0;
239 G4double range = 0.;
240 for (G4int i=0; i<Nbin; ++i)
241 {
242 e2 = (*Energy)[i];
243 dedx2 = 0.0;
244 for (G4int j=0; j<nelm; ++j)
245 {
246 dedx2 += dens[j]*ComputeValue((*elm)[j]->GetZasInt(), e2);
247 }
248 range += (dedx1 + dedx2 > 0.0) ? 2*(e2 - e1)/(dedx1 + dedx2) : 0.0;
249 range2 = range;
250 if(range2 < rangeCut)
251 {
252 e1 = e2;
253 dedx1 = dedx2;
254 range1 = range2;
255 }
256 else
257 {
258 break;
259 }
260 }
261 return LiniearInterpolation(e1, e2, range1, range2, rangeCut);
262}
263
264// --------------------------------------------------------------------
static const G4double e1[44]
static const G4double e2[44]
static const G4double emax
std::vector< const G4Element * > G4ElementVector
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double mm
Definition: G4SIunits.hh:95
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static void SetMaxEnergyCut(const G4double value)
virtual G4double ComputeValue(const G4int Z, const G4double kinEnergy)=0
G4double LiniearInterpolation(const G4double e1, const G4double e2, const G4double r1, const G4double r2, const G4double r)
G4double ConvertForGamma(const G4double rangeCut, const G4Material *material)
static std::vector< G4double > * Energy
G4double ConvertForElectron(const G4double rangeCut, const G4Material *material)
static void FillEnergyVector(const G4double emin, const G4double emax)
virtual G4double Convert(const G4double rangeCut, const G4Material *material)
static void SetEnergyRange(const G4double lowedge, const G4double highedge)
static constexpr double mm
Definition: SystemOfUnits.h:96
static constexpr double cm3
static constexpr double GeV
static constexpr double keV
static constexpr double g
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62