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