Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VelocityTable Class Reference

#include <G4VelocityTable.hh>

Public Member Functions

 G4VelocityTable ()
 
 ~G4VelocityTable ()
 
G4double Value (G4double theEnergy)
 

Static Public Member Functions

static G4VelocityTableGetVelocityTable ()
 
static void SetVelocityTableProperties (G4double t_max, G4double t_min, G4int nbin)
 
static G4double GetMaxTOfVelocityTable ()
 
static G4double GetMinTOfVelocityTable ()
 
static G4int GetNbinOfVelocityTable ()
 

Detailed Description

Definition at line 51 of file G4VelocityTable.hh.

Constructor & Destructor Documentation

◆ G4VelocityTable()

G4VelocityTable::G4VelocityTable ( )

Definition at line 52 of file G4VelocityTable.cc.

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{
59 PrepareVelocityTable();
60}
#define DBL_MAX
Definition: templates.hh:83

◆ ~G4VelocityTable()

G4VelocityTable::~G4VelocityTable ( )

Definition at line 63 of file G4VelocityTable.cc.

65{
66 dataVector.clear();
67 binVector.clear();
68}

Member Function Documentation

◆ GetMaxTOfVelocityTable()

G4double G4VelocityTable::GetMaxTOfVelocityTable ( )
static

Definition at line 190 of file G4VelocityTable.cc.

192{ return theInstance->maxT; }

Referenced by G4Track::GetMaxTOfVelocityTable().

◆ GetMinTOfVelocityTable()

G4double G4VelocityTable::GetMinTOfVelocityTable ( )
static

Definition at line 195 of file G4VelocityTable.cc.

197{ return theInstance->minT; }

Referenced by G4Track::GetMinTOfVelocityTable().

◆ GetNbinOfVelocityTable()

G4int G4VelocityTable::GetNbinOfVelocityTable ( )
static

Definition at line 200 of file G4VelocityTable.cc.

202{ return theInstance->NbinT; }

Referenced by G4Track::GetNbinOfVelocityTable().

◆ GetVelocityTable()

G4VelocityTable * G4VelocityTable::GetVelocityTable ( )
static

Definition at line 160 of file G4VelocityTable.cc.

162{
163 return theInstance;
164}

Referenced by G4Track::G4Track(), and G4Track::SetVelocityTableProperties().

◆ SetVelocityTableProperties()

void G4VelocityTable::SetVelocityTableProperties ( G4double  t_max,
G4double  t_min,
G4int  nbin 
)
static

Definition at line 167 of file G4VelocityTable.cc.

169{
171 G4ApplicationState currentState = stateManager->GetCurrentState();
172
173 // check if state is outside event loop
174 if(!(currentState==G4State_Idle||currentState==G4State_PreInit)){
175 G4Exception("G4VelocityTable::SetVelocityTableProperties",
176 "Track101", JustWarning,
177 "Can modify only in PreInit or Idle state : Method ignored.");
178 return;
179 }
180
181 if (nbin > 100 ) theInstance->NbinT = nbin;
182 if ((t_min < t_max)&&(t_min>0.)) {
183 theInstance->minT = t_min;
184 theInstance->maxT = t_max;
185 }
186 theInstance->PrepareVelocityTable();
187}
G4ApplicationState
@ G4State_Idle
@ G4State_PreInit
@ JustWarning
G4ApplicationState GetCurrentState() const
static G4StateManager * GetStateManager()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4Track::SetVelocityTableProperties().

◆ Value()

G4double G4VelocityTable::Value ( G4double  theEnergy)

Definition at line 123 of file G4VelocityTable.cc.

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;
135 lastValue = Interpolation();
136
137 } else if( theEnergy <= edgeMin ) {
138 lastBin = 0;
139 lastEnergy = edgeMin;
140 lastValue = dataVector[0];
141
142 } else if( theEnergy >= edgeMax ){
143 lastBin = numberOfNodes-1;
144 lastEnergy = edgeMax;
145 lastValue = dataVector[lastBin];
146
147 } else {
148 lastBin = FindBinLocation(theEnergy);
149 lastEnergy = theEnergy;
150 lastValue = Interpolation();
151
152 }
153 return lastValue;
154}

Referenced by G4Track::CalculateVelocity().


The documentation for this class was generated from the following files: