Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// $Id$
28//
29//
30// --------------------------------------------------------------
31// GEANT 4 class implementation file/ History:
32// 5 Oct. 2002, H.Kuirashige : Structure created based on object model
33// --------------------------------------------------------------
34
36#include "G4ParticleTable.hh"
37#include "G4Material.hh"
38#include "G4MaterialTable.hh"
39#include "G4PhysicsLogVector.hh"
40
41#include "G4ios.hh"
42#include "G4SystemOfUnits.hh"
43
44// energy range
47
48// max energy cut
50
52 theParticle(0), theLossTable(0), NumberOfElements(0), TotBin(300),
53 verboseLevel(1)
54{
55 fMaxEnergyCut = 0.;
56}
57
58G4VRangeToEnergyConverter::G4VRangeToEnergyConverter(const G4VRangeToEnergyConverter& right) : theParticle(right.theParticle), theLossTable(0), TotBin(right.TotBin)
59{
60 fMaxEnergyCut = 0.;
61 *this = right;
62}
63
65{
66 if (this == &right) return *this;
67 if (theLossTable) {
69 delete theLossTable;
71 }
72
75 verboseLevel = right.verboseLevel;
76
77 // create the loss table
80 // fill the loss table
81 for (size_t j=0; j<size_t(NumberOfElements); j++){
82 G4LossVector* aVector= new
84 for (size_t i=0; i<size_t(TotBin); i++) {
85 G4double Value = (*((*right.theLossTable)[j]))[i];
86 aVector->PutValue(i,Value);
87 }
88 theLossTable->insert(aVector);
89 }
90
91 // clean up range vector store
92 for (size_t idx=0; idx<fRangeVectorStore.size(); idx++){
93 delete fRangeVectorStore.at(idx);
94 }
95 fRangeVectorStore.clear();
96
97 // copy range vector store
98 for (size_t j=0; j<((right.fRangeVectorStore).size()); j++){
99 G4RangeVector* vector = (right.fRangeVectorStore).at(j);
100 G4RangeVector* rangeVector = 0;
101 if (vector !=0 ) {
102 rangeVector = new G4RangeVector(LowestEnergy, MaxEnergyCut, TotBin);
104 for (size_t i=0; i<size_t(TotBin); i++) {
105 G4double Value = (*vector)[i];
106 rangeVector->PutValue(i,Value);
107 }
108 }
109 fRangeVectorStore.push_back(rangeVector);
110 }
111 return *this;
112}
113
114
116{
117 Reset();
118}
119
121{
122 return this == &right;
123}
124
126{
127 return this != &right;
128}
129
130
131// **********************************************************************
132// ************************* Convert ***********************************
133// **********************************************************************
135 const G4Material* material)
136{
137#ifdef G4VERBOSE
138 if (GetVerboseLevel()>3) {
139 G4cout << "G4VRangeToEnergyConverter::Convert() ";
140 G4cout << "Convert for " << material->GetName()
141 << " with Range Cut " << rangeCut/mm << "[mm]" << G4endl;
142 }
143#endif
144
145 G4double theKineticEnergyCuts = 0.;
146
149 // clear loss table and renge vectors
150 Reset();
151 }
152
153 // Build the energy loss table
155
156 // Build range vector for every material, convert cut into energy-cut,
157 // fill theKineticEnergyCuts and delete the range vector
158 G4double tune = 0.025*mm*g/cm3 ,lowen = 30.*keV ;
159
160 // check density
161 G4double density = material->GetDensity() ;
162 if(density <= 0.) {
163 #ifdef G4VERBOSE
164 if (GetVerboseLevel()>0) {
165 G4cout << "G4VRangeToEnergyConverter::Convert() ";
166 G4cout << material->GetName() << "has zero density "
167 << "( " << density << ")" << G4endl;
168 }
169#endif
170 return 0.;
171 }
172
173 // initialize RangeVectorStore
175 G4int ext_size = table->size() - fRangeVectorStore.size();
176 for (int i=0; i<ext_size; i++) fRangeVectorStore.push_back(0);
177
178 // Build Range Vector
179 G4int idx = material->GetIndex();
180 G4RangeVector* rangeVector = fRangeVectorStore.at(idx);
181 if (rangeVector == 0) {
182 rangeVector = new G4RangeVector(LowestEnergy, MaxEnergyCut, TotBin);
183 BuildRangeVector(material, rangeVector);
184 fRangeVectorStore.at(idx) = rangeVector;
185 }
186
187 // Convert Range Cut ro Kinetic Energy Cut
188 theKineticEnergyCuts = ConvertCutToKineticEnergy(rangeVector, rangeCut, idx);
189
190 if( ((theParticle->GetParticleName()=="e-")||(theParticle->GetParticleName()=="e+"))
191 && (theKineticEnergyCuts < lowen) ) {
192 // corr. should be switched on smoothly
193 theKineticEnergyCuts /= (1.+(1.-theKineticEnergyCuts/lowen)*
194 tune/(rangeCut*density));
195 }
196
197 if(theKineticEnergyCuts < LowestEnergy) {
198 theKineticEnergyCuts = LowestEnergy ;
199 } else if(theKineticEnergyCuts > MaxEnergyCut) {
200 theKineticEnergyCuts = MaxEnergyCut;
201 }
202
203 return theKineticEnergyCuts;
204}
205
206// **********************************************************************
207// ************************ SetEnergyRange *****************************
208// **********************************************************************
210 G4double highedge)
211{
212 // check LowestEnergy/ HighestEnergy
213 if ( (lowedge<0.0)||(highedge<=lowedge) ){
214#ifdef G4VERBOSE
215 G4cerr << "Error in G4VRangeToEnergyConverter::SetEnergyRange";
216 G4cerr << " : illegal energy range" << "(" << lowedge/GeV;
217 G4cerr << "," << highedge/GeV << ") [GeV]" << G4endl;
218#endif
219 G4Exception( "G4VRangeToEnergyConverter::SetEnergyRange()",
220 "ProcCuts101",
221 JustWarning, "Illegal energy range ");
222 } else {
223 LowestEnergy = lowedge;
224 HighestEnergy = highedge;
225 }
226}
227
228
230{
231 return LowestEnergy;
232}
233
234
236{
237 return HighestEnergy;
238}
239
240// **********************************************************************
241// ******************* Get/SetMaxEnergyCut *****************************
242// **********************************************************************
244{
245 return MaxEnergyCut;
246}
247
249{
250 MaxEnergyCut = value;
251}
252
253// **********************************************************************
254// ************************ Reset **************************************
255// **********************************************************************
257{
258 // delete loss table
259 if (theLossTable) {
261 delete theLossTable;
262 }
263 theLossTable=0;
265
266 //clear RangeVectorStore
267 for (size_t idx=0; idx<fRangeVectorStore.size(); idx++){
268 delete fRangeVectorStore.at(idx);
269 }
270 fRangeVectorStore.clear();
271}
272
273
274// **********************************************************************
275// ************************ BuildLossTable ******************************
276// **********************************************************************
277// create Energy Loss Table for charged particles
278// (cross section tabel for neutral )
280{
281 if (size_t(NumberOfElements) == G4Element::GetNumberOfElements()) return;
282
283 // clear Loss table and Range vectors
284 Reset();
285
286 // Build dE/dx tables for elements
290#ifdef G4VERBOSE
291 if (GetVerboseLevel()>3) {
292 G4cout << "G4VRangeToEnergyConverter::BuildLossTable() ";
293 G4cout << "Create theLossTable[" << theLossTable << "]";
294 G4cout << " NumberOfElements=" << NumberOfElements <<G4endl;
295 }
296#endif
297
298
299 // fill the loss table
300 for (size_t j=0; j<size_t(NumberOfElements); j++){
301 G4double Value;
302 G4LossVector* aVector= 0;
304 for (size_t i=0; i<size_t(TotBin); i++) {
305 Value = ComputeLoss( (*G4Element::GetElementTable())[j]->GetZ(),
306 aVector->GetLowEdgeEnergy(i)
307 );
308 aVector->PutValue(i,Value);
309 }
310 theLossTable->insert(aVector);
311 }
312}
313
314// **********************************************************************
315// ************************ BuildRangeVector ****************************
316// **********************************************************************
318 G4PhysicsLogVector* rangeVector)
319{
320 // create range vector for a material
321 const G4ElementVector* elementVector = aMaterial->GetElementVector();
322 const G4double* atomicNumDensityVector = aMaterial->GetAtomicNumDensityVector();
323 G4int NumEl = aMaterial->GetNumberOfElements();
324
325 // calculate parameters of the low energy part first
326 size_t i;
327 std::vector<G4double> lossV;
328 for ( size_t ib=0; ib<size_t(TotBin); ib++) {
329 G4double loss=0.;
330 for (i=0; i<size_t(NumEl); i++) {
331 G4int IndEl = (*elementVector)[i]->GetIndex();
332 loss += atomicNumDensityVector[i]*
333 (*((*theLossTable)[IndEl]))[ib];
334 }
335 lossV.push_back(loss);
336 }
337
338 // Integrate with Simpson formula with logarithmic binning
339 G4double ltt = std::log(MaxEnergyCut/LowestEnergy);
340 G4double dltau = ltt/TotBin;
341
342 G4double s0 = 0.;
343 G4double Value;
344 for ( i=0; i<size_t(TotBin); i++) {
345 G4double t = rangeVector->GetLowEdgeEnergy(i);
346 G4double q = t/lossV[i];
347 if (i==0) s0 += 0.5*q;
348 else s0 += q;
349
350 if (i==0) {
351 Value = (s0 + 0.5*q)*dltau ;
352 } else {
353 Value = (s0 - 0.5*q)*dltau ;
354 }
355 rangeVector->PutValue(i,Value);
356 }
357}
358
359// **********************************************************************
360// ****************** ConvertCutToKineticEnergy *************************
361// **********************************************************************
363 G4RangeVector* rangeVector,
364 G4double theCutInLength,
365#ifdef G4VERBOSE
366 size_t materialIndex
367#else
368 size_t
369#endif
370 ) const
371{
372 const G4double epsilon=0.01;
373
374 // find max. range and the corresponding energy (rmax,Tmax)
375 G4double rmax= -1.e10*mm;
376
378 G4double r1 =(*rangeVector)[0] ;
379
381
382 // check theCutInLength < r1
383 if ( theCutInLength <= r1 ) { return T1; }
384
385 // scan range vector to find nearest bin
386 // ( suppose that r(Ti) > r(Tj) if Ti >Tj )
387 for (size_t ibin=0; ibin<size_t(TotBin); ibin++) {
388 G4double T=rangeVector->GetLowEdgeEnergy(ibin);
389 G4double r=(*rangeVector)[ibin];
390 if ( r>rmax ) rmax=r;
391 if (r <theCutInLength ) {
392 T1 = T;
393 r1 = r;
394 } else if (r >theCutInLength ) {
395 T2 = T;
396 break;
397 }
398 }
399
400 // check cut in length is smaller than range max
401 if ( theCutInLength >= rmax ) {
402#ifdef G4VERBOSE
403 if (GetVerboseLevel()>2) {
404 G4cout << "G4VRangeToEnergyConverter::ConvertCutToKineticEnergy ";
405 G4cout << " for " << theParticle->GetParticleName() << G4endl;
406 G4cout << "The cut in range [" << theCutInLength/mm << " (mm)] ";
407 G4cout << " is too big " ;
408 G4cout << " for material idx=" << materialIndex <<G4endl;
409 }
410#endif
411 return MaxEnergyCut;
412 }
413
414 // convert range to energy
415 G4double T3 = std::sqrt(T1*T2);
416 G4double r3 = rangeVector->Value(T3);
417 while ( std::fabs(1.-r3/theCutInLength)>epsilon ) {
418 if ( theCutInLength <= r3 ) {
419 T2 = T3;
420 } else {
421 T1 = T3;
422 }
423 T3 = std::sqrt(T1*T2);
424 r3 = rangeVector->Value(T3);
425 }
426
427 return T3;
428}
429
std::vector< G4Element * > G4ElementVector
@ JustWarning
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
G4double GetDensity() const
Definition: G4Material.hh:179
static const G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:562
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
const G4String & GetName() const
Definition: G4Material.hh:177
size_t GetIndex() const
Definition: G4Material.hh:261
const G4String & GetParticleName() const
void clearAndDestroy()
void insert(G4PhysicsVector *)
G4double Value(G4double theEnergy)
virtual G4double GetLowEdgeEnergy(size_t binNumber) const
void PutValue(size_t index, G4double theValue)
std::vector< G4RangeVector * > fRangeVectorStore
virtual void BuildRangeVector(const G4Material *aMaterial, G4RangeVector *rangeVector)
G4int operator==(const G4VRangeToEnergyConverter &right) const
virtual G4double Convert(G4double rangeCut, const G4Material *material)
G4int operator!=(const G4VRangeToEnergyConverter &right) const
G4VRangeToEnergyConverter & operator=(const G4VRangeToEnergyConverter &right)
static void SetMaxEnergyCut(G4double value)
virtual G4double ComputeLoss(G4double AtomicNumber, G4double KineticEnergy) const =0
static void SetEnergyRange(G4double lowedge, G4double highedge)
const G4ParticleDefinition * theParticle
G4double ConvertCutToKineticEnergy(G4RangeVector *theRangeVector, G4double theCutInLength, size_t materialIndex) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41