Geant4  10.01.p02
G4VelocityTable.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: G4VelocityTable.cc 69005 2013-04-15 09:26:47Z gcosmo $
28 //
29 //
30 //---------------------------------------------------------------
31 //
32 // G4VelocityTable.cc
33 //
34 // class description:
35 // This class keeps a table of velocity as a function of
36 // the ratio kinetic erngy and mass
37 //
38 //---------------------------------------------------------------
39 // created 17.Aug. 2011 H.Kurashige
40 //
41 
42 #include "G4VelocityTable.hh"
43 #include "G4PhysicalConstants.hh"
44 #include "G4StateManager.hh"
45 #include "G4ApplicationState.hh"
46 
47 #include "G4ios.hh"
48 
50 
54  : edgeMin(0.), edgeMax(0.), numberOfNodes(0),
55  dBin(0.), baseBin(0.),
56  lastEnergy(-DBL_MAX), lastValue(0.), lastBin(0),
57  maxT( 1000.0 ), minT( 0.0001 ), NbinT( 500 )
58 {
60 }
61 
65 {
66  dataVector.clear();
67  binVector.clear();
68 }
69 
70 
74 {
75  dataVector.clear();
76  binVector.clear();
77 
78  dBin = std::log10(maxT/minT)/NbinT;
79  baseBin = std::log10(minT)/dBin;
80 
81  numberOfNodes = NbinT + 1;
82  dataVector.reserve(numberOfNodes);
83  binVector.reserve(numberOfNodes);
84 
85  const G4double g4log10 = std::log(10.);
86 
87  binVector.push_back(minT);
88  dataVector.push_back(0.0);
89 
90  for (size_t i=1; i<numberOfNodes-1; i++){
91  binVector.push_back(std::exp(g4log10*(baseBin+i)*dBin));
92  dataVector.push_back(0.0);
93  }
94  binVector.push_back(maxT);
95  dataVector.push_back(0.0);
96 
97  edgeMin = binVector[0];
98  edgeMax = binVector[numberOfNodes-1];
99 
100  for (G4int i=0; i<=NbinT; i++){
101  G4double T = binVector[i];
102  dataVector[i]= c_light*std::sqrt(T*(T+2.))/(T+1.0);
103  }
104 
105  return;
106 }
107 
109 {
110  // For G4PhysicsLogVector, FindBinLocation is implemented using
111  // a simple arithmetic calculation.
112  //
113  // Because this is a virtual function, it is accessed through a
114  // pointer to the G4PhyiscsVector object for most usages. In this
115  // case, 'inline' will not be invoked. However, there is a possibility
116  // that the user access to the G4PhysicsLogVector object directly and
117  // not through pointers or references. In this case, the 'inline' will
118  // be invoked. (See R.B.Murray, "C++ Strategies and Tactics", Chap.6.6)
119 
120  return size_t( std::log10(theEnergy)/dBin - baseBin );
121 }
122 
124 {
125  // Use cache for speed up - check if the value 'theEnergy' is same as the
126  // last call. If it is same, then use the last bin location. Also the
127  // value 'theEnergy' lies between the last energy and low edge of of the
128  // bin of last call, then the last bin location is used.
129 
130  if( theEnergy == lastEnergy ) {
131 
132  } else if( theEnergy < lastEnergy
133  && theEnergy >= binVector[lastBin]) {
134  lastEnergy = theEnergy;
136 
137  } else if( theEnergy <= edgeMin ) {
138  lastBin = 0;
140  lastValue = dataVector[0];
141 
142  } else if( theEnergy >= edgeMax ){
143  lastBin = numberOfNodes-1;
146 
147  } else {
148  lastBin = FindBinLocation(theEnergy);
149  lastEnergy = theEnergy;
151 
152  }
153  return lastValue;
154 }
155 
157 // Static methods
158 
162 {
163  if (!theInstance) { theInstance = new G4VelocityTable(); }
164  return theInstance;
165 }
166 
170 {
171  if (theInstance == 0) { theInstance = new G4VelocityTable(); }
172 
174  G4ApplicationState currentState = stateManager->GetCurrentState();
175 
176  // check if state is outside event loop
177  if(!(currentState==G4State_Idle||currentState==G4State_PreInit)){
178  G4Exception("G4VelocityTable::SetVelocityTableProperties",
179  "Track101", JustWarning,
180  "Can modify only in PreInit or Idle state : Method ignored.");
181  return;
182  }
183 
184  if (nbin > 100 ) theInstance->NbinT = nbin;
185  if ((t_min < t_max)&&(t_min>0.)) {
186  theInstance->minT = t_min;
187  theInstance->maxT = t_max;
188  }
190 }
191 
195 {
196  return GetVelocityTable()->maxT;
197 }
198 
201 {
203  return GetVelocityTable()->minT;
204 }
205 
208 {
210  return GetVelocityTable()->NbinT;
211 }
static G4VelocityTable * GetVelocityTable()
static G4double GetMinTOfVelocityTable()
void PrepareVelocityTable()
G4VTDataVector binVector
G4VTDataVector dataVector
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
static G4StateManager * GetStateManager()
size_t FindBinLocation(G4double theEnergy) const
static G4int GetNbinOfVelocityTable()
G4ApplicationState GetCurrentState() const
static void SetVelocityTableProperties(G4double t_max, G4double t_min, G4int nbin)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ThreadLocal G4VelocityTable * theInstance
double G4double
Definition: G4Types.hh:76
static G4double GetMaxTOfVelocityTable()
#define DBL_MAX
Definition: templates.hh:83
G4double Value(G4double theEnergy)
G4double Interpolation() const
G4ApplicationState