Geant4  10.02.p03
G4PhysicsVector Class Reference

#include <G4PhysicsVector.hh>

Inheritance diagram for G4PhysicsVector:
Collaboration diagram for G4PhysicsVector:

Public Member Functions

 G4PhysicsVector (G4bool spline=false)
 
 G4PhysicsVector (const G4PhysicsVector &)
 
G4PhysicsVectoroperator= (const G4PhysicsVector &)
 
virtual ~G4PhysicsVector ()
 
G4double Value (G4double theEnergy, size_t &lastidx) const
 
G4double Value (G4double theEnergy) const
 
G4double GetValue (G4double theEnergy, G4bool &isOutRange) const
 
G4int operator== (const G4PhysicsVector &right) const
 
G4int operator!= (const G4PhysicsVector &right) const
 
G4double operator[] (const size_t binNumber) const
 
G4double operator() (const size_t binNumber) const
 
void PutValue (size_t index, G4double theValue)
 
virtual void ScaleVector (G4double factorE, G4double factorV)
 
G4double Energy (size_t index) const
 
G4double GetMaxEnergy () const
 
G4double GetLowEdgeEnergy (size_t binNumber) const
 
size_t GetVectorLength () const
 
size_t FindBin (G4double energy, size_t idx) const
 
void FillSecondDerivatives ()
 
void ComputeSecDerivatives ()
 
void ComputeSecondDerivatives (G4double firstPointDerivative, G4double endPointDerivative)
 
G4double FindLinearEnergy (G4double rand) const
 
G4bool IsFilledVectorExist () const
 
G4PhysicsVectorType GetType () const
 
void SetSpline (G4bool)
 
virtual G4bool Store (std::ofstream &fOut, G4bool ascii=false)
 
virtual G4bool Retrieve (std::ifstream &fIn, G4bool ascii=false)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel (G4int)
 

Protected Member Functions

void DeleteData ()
 
void CopyData (const G4PhysicsVector &vec)
 

Protected Attributes

G4PhysicsVectorType type
 
G4double edgeMin
 
G4double edgeMax
 
size_t numberOfNodes
 
G4PVDataVector dataVector
 
G4PVDataVector binVector
 
G4PVDataVector secDerivative
 
G4double dBin
 
G4double baseBin
 
G4int verboseLevel
 

Private Member Functions

G4bool SplinePossible ()
 
G4double LinearInterpolation (size_t idx, G4double energy) const
 
G4double SplineInterpolation (size_t idx, G4double energy) const
 
G4double Interpolation (size_t idx, G4double energy) const
 
size_t FindBinLocation (G4double theEnergy) const
 

Private Attributes

G4bool useSpline
 

Friends

std::ostream & operator<< (std::ostream &, const G4PhysicsVector &)
 

Detailed Description

Definition at line 76 of file G4PhysicsVector.hh.

Constructor & Destructor Documentation

◆ G4PhysicsVector() [1/2]

G4PhysicsVector::G4PhysicsVector ( G4bool  spline = false)

Definition at line 61 of file G4PhysicsVector.cc.

63  edgeMin(0.), edgeMax(0.), numberOfNodes(0),
64  useSpline(false),
65  dBin(0.), baseBin(0.),
66  verboseLevel(0)
67 {
68 }
G4PhysicsVectorType type

◆ G4PhysicsVector() [2/2]

G4PhysicsVector::G4PhysicsVector ( const G4PhysicsVector right)

Definition at line 78 of file G4PhysicsVector.cc.

79 {
80  dBin = right.dBin;
81  baseBin = right.baseBin;
82  verboseLevel = right.verboseLevel;
83 
84  DeleteData();
85  CopyData(right);
86 }
void CopyData(const G4PhysicsVector &vec)
Here is the call graph for this function:

◆ ~G4PhysicsVector()

G4PhysicsVector::~G4PhysicsVector ( )
virtual

Definition at line 72 of file G4PhysicsVector.cc.

73 {
74 }

Member Function Documentation

◆ ComputeSecDerivatives()

void G4PhysicsVector::ComputeSecDerivatives ( )

Definition at line 431 of file G4PhysicsVector.cc.

433 {
434  if(3 > numberOfNodes) // cannot compute derivatives for less than 4 bins
435  {
436  useSpline = false;
437  return;
438  }
439 
440  if(!SplinePossible()) { return; }
441 
442  useSpline = true;
443 
444  size_t n = numberOfNodes-1;
445 
446  for(size_t i=1; i<n; ++i)
447  {
448  secDerivative[i] =
449  3.0*((dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i]) -
450  (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]))
451  /(binVector[i+1]-binVector[i-1]);
452  }
453  secDerivative[n] = secDerivative[n-1];
455 }
G4PVDataVector dataVector
G4PVDataVector binVector
Char_t n[5]
G4PVDataVector secDerivative
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeSecondDerivatives()

void G4PhysicsVector::ComputeSecondDerivatives ( G4double  firstPointDerivative,
G4double  endPointDerivative 
)

Definition at line 294 of file G4PhysicsVector.cc.

300 {
301  if(4 > numberOfNodes) // cannot compute derivatives for less than 4 bins
302  {
304  return;
305  }
306 
307  if(!SplinePossible()) { return; }
308 
309  useSpline = true;
310 
311  G4int n = numberOfNodes-1;
312 
313  G4double* u = new G4double [n];
314 
315  G4double p, sig, un;
316 
317  u[0] = (6.0/(binVector[1]-binVector[0]))
318  * ((dataVector[1]-dataVector[0])/(binVector[1]-binVector[0])
319  - firstPointDerivative);
320 
321  secDerivative[0] = - 0.5;
322 
323  // Decomposition loop for tridiagonal algorithm. secDerivative[i]
324  // and u[i] are used for temporary storage of the decomposed factors.
325 
326  for(G4int i=1; i<n; ++i)
327  {
328  sig = (binVector[i]-binVector[i-1]) / (binVector[i+1]-binVector[i-1]);
329  p = sig*(secDerivative[i-1]) + 2.0;
330  secDerivative[i] = (sig - 1.0)/p;
331  u[i] = (dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i])
332  - (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]);
333  u[i] = 6.0*u[i]/(binVector[i+1]-binVector[i-1]) - sig*u[i-1]/p;
334  }
335 
336  sig = (binVector[n-1]-binVector[n-2]) / (binVector[n]-binVector[n-2]);
337  p = sig*secDerivative[n-2] + 2.0;
338  un = (6.0/(binVector[n]-binVector[n-1]))
339  *(endPointDerivative -
340  (dataVector[n]-dataVector[n-1])/(binVector[n]-binVector[n-1])) - u[n-1]/p;
341  secDerivative[n] = un/(secDerivative[n-1] + 2.0);
342 
343  // The back-substitution loop for the triagonal algorithm of solving
344  // a linear system of equations.
345 
346  for(G4int k=n-1; k>0; --k)
347  {
348  secDerivative[k] *=
349  (secDerivative[k+1] -
350  u[k]*(binVector[k+1]-binVector[k-1])/(binVector[k+1]-binVector[k]));
351  }
352  secDerivative[0] = 0.5*(u[0] - secDerivative[1]);
353 
354  delete [] u;
355 }
G4PVDataVector dataVector
int G4int
Definition: G4Types.hh:78
G4PVDataVector binVector
Char_t n[5]
G4PVDataVector secDerivative
double G4double
Definition: G4Types.hh:76
void ComputeSecDerivatives()
Here is the call graph for this function:

◆ CopyData()

void G4PhysicsVector::CopyData ( const G4PhysicsVector vec)
protected

Definition at line 126 of file G4PhysicsVector.cc.

127 {
128  type = vec.type;
129  edgeMin = vec.edgeMin;
130  edgeMax = vec.edgeMax;
132  useSpline = vec.useSpline;
133 
134  size_t i;
135  dataVector.resize(numberOfNodes);
136  for(i=0; i<numberOfNodes; ++i) {
137  dataVector[i] = (vec.dataVector)[i];
138  }
139  binVector.resize(numberOfNodes);
140  for(i=0; i<numberOfNodes; ++i) {
141  binVector[i] = (vec.binVector)[i];
142  }
143  if(0 < (vec.secDerivative).size()) {
144  secDerivative.resize(numberOfNodes);
145  for(i=0; i<numberOfNodes; ++i){
146  secDerivative[i] = (vec.secDerivative)[i];
147  }
148  }
149 }
G4PVDataVector dataVector
G4PVDataVector binVector
G4PhysicsVectorType type
G4PVDataVector secDerivative
Here is the caller graph for this function:

◆ DeleteData()

void G4PhysicsVector::DeleteData ( )
protected

Definition at line 118 of file G4PhysicsVector.cc.

119 {
120  useSpline = false;
121  secDerivative.clear();
122 }
G4PVDataVector secDerivative
Here is the caller graph for this function:

◆ Energy()

G4double G4PhysicsVector::Energy ( size_t  index) const
inline

◆ FillSecondDerivatives()

void G4PhysicsVector::FillSecondDerivatives ( )

Definition at line 359 of file G4PhysicsVector.cc.

363 {
364  if(5 > numberOfNodes) // cannot compute derivatives for less than 4 points
365  {
367  return;
368  }
369 
370  if(!SplinePossible()) { return; }
371 
372  useSpline = true;
373 
374  G4int n = numberOfNodes-1;
375 
376  G4double* u = new G4double [n];
377 
378  G4double p, sig;
379 
380  u[1] = ((dataVector[2]-dataVector[1])/(binVector[2]-binVector[1]) -
381  (dataVector[1]-dataVector[0])/(binVector[1]-binVector[0]));
382  u[1] = 6.0*u[1]*(binVector[2]-binVector[1])
383  / ((binVector[2]-binVector[0])*(binVector[2]-binVector[0]));
384 
385  // Decomposition loop for tridiagonal algorithm. secDerivative[i]
386  // and u[i] are used for temporary storage of the decomposed factors.
387 
388  secDerivative[1] = (2.0*binVector[1]-binVector[0]-binVector[2])
389  / (2.0*binVector[2]-binVector[0]-binVector[1]);
390 
391  for(G4int i=2; i<n-1; ++i)
392  {
393  sig = (binVector[i]-binVector[i-1]) / (binVector[i+1]-binVector[i-1]);
394  p = sig*secDerivative[i-1] + 2.0;
395  secDerivative[i] = (sig - 1.0)/p;
396  u[i] = (dataVector[i+1]-dataVector[i])/(binVector[i+1]-binVector[i])
397  - (dataVector[i]-dataVector[i-1])/(binVector[i]-binVector[i-1]);
398  u[i] = (6.0*u[i]/(binVector[i+1]-binVector[i-1])) - sig*u[i-1]/p;
399  }
400 
401  sig = (binVector[n-1]-binVector[n-2]) / (binVector[n]-binVector[n-2]);
402  p = sig*secDerivative[n-3] + 2.0;
403  u[n-1] = (dataVector[n]-dataVector[n-1])/(binVector[n]-binVector[n-1])
404  - (dataVector[n-1]-dataVector[n-2])/(binVector[n-1]-binVector[n-2]);
405  u[n-1] = 6.0*sig*u[n-1]/(binVector[n]-binVector[n-2])
406  - (2.0*sig - 1.0)*u[n-2]/p;
407 
408  p = (1.0+sig) + (2.0*sig-1.0)*secDerivative[n-2];
409  secDerivative[n-1] = u[n-1]/p;
410 
411  // The back-substitution loop for the triagonal algorithm of solving
412  // a linear system of equations.
413 
414  for(G4int k=n-2; k>1; --k)
415  {
416  secDerivative[k] *=
417  (secDerivative[k+1] -
418  u[k]*(binVector[k+1]-binVector[k-1])/(binVector[k+1]-binVector[k]));
419  }
420  secDerivative[n] = (secDerivative[n-1] - (1.0-sig)*secDerivative[n-2])/sig;
421  sig = 1.0 - ((binVector[2]-binVector[1])/(binVector[2]-binVector[0]));
422  secDerivative[1] *= (secDerivative[2] - u[1]/(1.0-sig));
423  secDerivative[0] = (secDerivative[1] - sig*secDerivative[2])/(1.0-sig);
424 
425  delete [] u;
426 }
G4PVDataVector dataVector
int G4int
Definition: G4Types.hh:78
G4PVDataVector binVector
Char_t n[5]
G4PVDataVector secDerivative
double G4double
Definition: G4Types.hh:76
void ComputeSecDerivatives()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ FindBin()

size_t G4PhysicsVector::FindBin ( G4double  energy,
size_t  idx 
) const
inline
Here is the caller graph for this function:

◆ FindBinLocation()

size_t G4PhysicsVector::FindBinLocation ( G4double  theEnergy) const
inlineprivate

◆ FindLinearEnergy()

G4double G4PhysicsVector::FindLinearEnergy ( G4double  rand) const

Definition at line 516 of file G4PhysicsVector.cc.

517 {
518  if(1 >= numberOfNodes) { return 0.0; }
520  size_t bin = std::lower_bound(dataVector.begin(), dataVector.end(), y)
521  - dataVector.begin() - 1;
522  bin = std::min(bin, numberOfNodes-2);
523  G4double res = binVector[bin];
524  G4double del = dataVector[bin+1] - dataVector[bin];
525  if(del > 0.0) {
526  res += (y - dataVector[bin])*(binVector[bin+1] - res)/del;
527  }
528  return res;
529 }
G4PVDataVector dataVector
float bin[41]
Definition: plottest35.C:14
G4PVDataVector binVector
Double_t y
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetLowEdgeEnergy()

G4double G4PhysicsVector::GetLowEdgeEnergy ( size_t  binNumber) const

Definition at line 153 of file G4PhysicsVector.cc.

154 {
155  return binVector[binNumber];
156 }
G4PVDataVector binVector

◆ GetMaxEnergy()

G4double G4PhysicsVector::GetMaxEnergy ( ) const
inline
Here is the caller graph for this function:

◆ GetType()

G4PhysicsVectorType G4PhysicsVector::GetType ( ) const
inline

◆ GetValue()

G4double G4PhysicsVector::GetValue ( G4double  theEnergy,
G4bool isOutRange 
) const
inline
Here is the caller graph for this function:

◆ GetVectorLength()

size_t G4PhysicsVector::GetVectorLength ( ) const
inline

◆ GetVerboseLevel()

G4int G4PhysicsVector::GetVerboseLevel ( G4int  )
inline

◆ Interpolation()

G4double G4PhysicsVector::Interpolation ( size_t  idx,
G4double  energy 
) const
inlineprivate
Here is the caller graph for this function:

◆ IsFilledVectorExist()

G4bool G4PhysicsVector::IsFilledVectorExist ( ) const
inline
Here is the caller graph for this function:

◆ LinearInterpolation()

G4double G4PhysicsVector::LinearInterpolation ( size_t  idx,
G4double  energy 
) const
inlineprivate

◆ operator!=()

G4int G4PhysicsVector::operator!= ( const G4PhysicsVector right) const

Definition at line 111 of file G4PhysicsVector.cc.

112 {
113  return (this != &right);
114 }

◆ operator()()

G4double G4PhysicsVector::operator() ( const size_t  binNumber) const
inline

◆ operator=()

G4PhysicsVector & G4PhysicsVector::operator= ( const G4PhysicsVector right)

Definition at line 90 of file G4PhysicsVector.cc.

91 {
92  if (&right==this) { return *this; }
93  dBin = right.dBin;
94  baseBin = right.baseBin;
95  verboseLevel = right.verboseLevel;
96 
97  DeleteData();
98  CopyData(right);
99  return *this;
100 }
void CopyData(const G4PhysicsVector &vec)
Here is the call graph for this function:

◆ operator==()

G4int G4PhysicsVector::operator== ( const G4PhysicsVector right) const

Definition at line 104 of file G4PhysicsVector.cc.

105 {
106  return (this == &right);
107 }

◆ operator[]()

G4double G4PhysicsVector::operator[] ( const size_t  binNumber) const
inline

◆ PutValue()

void G4PhysicsVector::PutValue ( size_t  index,
G4double  theValue 
)
inline

◆ Retrieve()

G4bool G4PhysicsVector::Retrieve ( std::ifstream &  fIn,
G4bool  ascii = false 
)
virtual

Reimplemented in G4PhysicsLogVector, G4PhysicsLinearVector, and G4PhysicsLnVector.

Definition at line 193 of file G4PhysicsVector.cc.

194 {
195  // clear properties;
196  dataVector.clear();
197  binVector.clear();
198  secDerivative.clear();
199 
200  // retrieve in ascii mode
201  if (ascii){
202  // binning
203  fIn >> edgeMin >> edgeMax >> numberOfNodes;
204  if (fIn.fail()) { return false; }
205  // contents
206  G4int siz=0;
207  fIn >> siz;
208  if (fIn.fail()) { return false; }
209  if (siz<=0)
210  {
211 #ifdef G4VERBOSE
212  G4cerr << "G4PhysicsVector::Retrieve():";
213  G4cerr << " Invalid vector size: " << siz << G4endl;
214 #endif
215  return false;
216  }
217 
218  binVector.reserve(siz);
219  dataVector.reserve(siz);
220  G4double vBin, vData;
221 
222  for(G4int i = 0; i < siz ; i++)
223  {
224  vBin = 0.;
225  vData= 0.;
226  fIn >> vBin >> vData;
227  if (fIn.fail()) { return false; }
228  binVector.push_back(vBin);
229  dataVector.push_back(vData);
230  }
231 
232  // to remove any inconsistency
233  numberOfNodes = siz;
234  edgeMin = binVector[0];
235  edgeMax = binVector[numberOfNodes-1];
236  return true ;
237  }
238 
239  // retrieve in binary mode
240  // binning
241  fIn.read((char*)(&edgeMin), sizeof edgeMin);
242  fIn.read((char*)(&edgeMax), sizeof edgeMax);
243  fIn.read((char*)(&numberOfNodes), sizeof numberOfNodes );
244 
245  // contents
246  size_t size;
247  fIn.read((char*)(&size), sizeof size);
248 
249  G4double* value = new G4double[2*size];
250  fIn.read((char*)(value), 2*size*(sizeof(G4double)) );
251  if (G4int(fIn.gcount()) != G4int(2*size*(sizeof(G4double))) )
252  {
253  delete [] value;
254  return false;
255  }
256 
257  binVector.reserve(size);
258  dataVector.reserve(size);
259  for(size_t i = 0; i < size; ++i)
260  {
261  binVector.push_back(value[2*i]);
262  dataVector.push_back(value[2*i+1]);
263  }
264  delete [] value;
265 
266  // to remove any inconsistency
267  numberOfNodes = size;
268  edgeMin = binVector[0];
269  edgeMax = binVector[numberOfNodes-1];
270 
271  return true;
272 }
G4PVDataVector dataVector
int G4int
Definition: G4Types.hh:78
G4PVDataVector binVector
G4PVDataVector secDerivative
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr
Here is the caller graph for this function:

◆ ScaleVector()

void G4PhysicsVector::ScaleVector ( G4double  factorE,
G4double  factorV 
)
virtual

Reimplemented in G4PhysicsLogVector, G4PhysicsLinearVector, and G4PhysicsLnVector.

Definition at line 277 of file G4PhysicsVector.cc.

278 {
279  size_t n = dataVector.size();
280  size_t i;
281  for(i=0; i<n; ++i) {
282  binVector[i] *= factorE;
283  dataVector[i] *= factorV;
284  }
285  secDerivative.clear();
286 
287  edgeMin *= factorE;
288  edgeMax *= factorE;
289 }
G4PVDataVector dataVector
G4PVDataVector binVector
Char_t n[5]
G4PVDataVector secDerivative
Here is the caller graph for this function:

◆ SetSpline()

void G4PhysicsVector::SetSpline ( G4bool  )
inline
Here is the caller graph for this function:

◆ SetVerboseLevel()

void G4PhysicsVector::SetVerboseLevel ( G4int  value)
inline

◆ SplineInterpolation()

G4double G4PhysicsVector::SplineInterpolation ( size_t  idx,
G4double  energy 
) const
inlineprivate

◆ SplinePossible()

G4bool G4PhysicsVector::SplinePossible ( )
private

Definition at line 459 of file G4PhysicsVector.cc.

462 {
463  G4bool result = true;
464  for(size_t j=1; j<numberOfNodes; ++j)
465  {
466  if(binVector[j] <= binVector[j-1]) {
467  result = false;
468  useSpline = false;
469  secDerivative.clear();
470  break;
471  }
472  }
473  secDerivative.resize(numberOfNodes,0.0);
474  return result;
475 }
G4PVDataVector binVector
bool G4bool
Definition: G4Types.hh:79
G4PVDataVector secDerivative
Here is the caller graph for this function:

◆ Store()

G4bool G4PhysicsVector::Store ( std::ofstream &  fOut,
G4bool  ascii = false 
)
virtual

Definition at line 160 of file G4PhysicsVector.cc.

161 {
162  // Ascii mode
163  if (ascii)
164  {
165  fOut << *this;
166  return true;
167  }
168  // Binary Mode
169 
170  // binning
171  fOut.write((char*)(&edgeMin), sizeof edgeMin);
172  fOut.write((char*)(&edgeMax), sizeof edgeMax);
173  fOut.write((char*)(&numberOfNodes), sizeof numberOfNodes);
174 
175  // contents
176  size_t size = dataVector.size();
177  fOut.write((char*)(&size), sizeof size);
178 
179  G4double* value = new G4double[2*size];
180  for(size_t i = 0; i < size; ++i)
181  {
182  value[2*i] = binVector[i];
183  value[2*i+1]= dataVector[i];
184  }
185  fOut.write((char*)(value), 2*size*(sizeof (G4double)));
186  delete [] value;
187 
188  return true;
189 }
G4PVDataVector dataVector
G4PVDataVector binVector
double G4double
Definition: G4Types.hh:76

◆ Value() [1/2]

G4double G4PhysicsVector::Value ( G4double  theEnergy,
size_t &  lastidx 
) const

Definition at line 498 of file G4PhysicsVector.cc.

499 {
500  G4double y;
501  if(theEnergy <= edgeMin) {
502  lastIdx = 0;
503  y = dataVector[0];
504  } else if(theEnergy >= edgeMax) {
505  lastIdx = numberOfNodes-1;
506  y = dataVector[lastIdx];
507  } else {
508  lastIdx = FindBin(theEnergy, lastIdx);
509  y = Interpolation(lastIdx, theEnergy);
510  }
511  return y;
512 }
G4PVDataVector dataVector
G4double Interpolation(size_t idx, G4double energy) const
size_t FindBin(G4double energy, size_t idx) const
Double_t y
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ Value() [2/2]

G4double G4PhysicsVector::Value ( G4double  theEnergy) const
inline

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const G4PhysicsVector pv 
)
friend

Definition at line 479 of file G4PhysicsVector.cc.

480 {
481  // binning
482  out << std::setprecision(12) << pv.edgeMin << " "
483  << pv.edgeMax << " " << pv.numberOfNodes << G4endl;
484 
485  // contents
486  out << pv.dataVector.size() << G4endl;
487  for(size_t i = 0; i < pv.dataVector.size(); i++)
488  {
489  out << pv.binVector[i] << " " << pv.dataVector[i] << G4endl;
490  }
491  out << std::setprecision(6);
492 
493  return out;
494 }
G4PVDataVector dataVector
G4PVDataVector binVector
#define G4endl
Definition: G4ios.hh:61

Member Data Documentation

◆ baseBin

G4double G4PhysicsVector::baseBin
protected

Definition at line 237 of file G4PhysicsVector.hh.

◆ binVector

G4PVDataVector G4PhysicsVector::binVector
protected

Definition at line 215 of file G4PhysicsVector.hh.

◆ dataVector

G4PVDataVector G4PhysicsVector::dataVector
protected

Definition at line 214 of file G4PhysicsVector.hh.

◆ dBin

G4double G4PhysicsVector::dBin
protected

Definition at line 236 of file G4PhysicsVector.hh.

◆ edgeMax

G4double G4PhysicsVector::edgeMax
protected

Definition at line 210 of file G4PhysicsVector.hh.

◆ edgeMin

G4double G4PhysicsVector::edgeMin
protected

Definition at line 209 of file G4PhysicsVector.hh.

◆ numberOfNodes

size_t G4PhysicsVector::numberOfNodes
protected

Definition at line 212 of file G4PhysicsVector.hh.

◆ secDerivative

G4PVDataVector G4PhysicsVector::secDerivative
protected

Definition at line 216 of file G4PhysicsVector.hh.

◆ type

G4PhysicsVectorType G4PhysicsVector::type
protected

Definition at line 207 of file G4PhysicsVector.hh.

◆ useSpline

G4bool G4PhysicsVector::useSpline
private

Definition at line 232 of file G4PhysicsVector.hh.

◆ verboseLevel

G4int G4PhysicsVector::verboseLevel
protected

Definition at line 239 of file G4PhysicsVector.hh.


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