Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4PhysicsLinearVector Class Reference

#include <G4PhysicsLinearVector.hh>

Inheritance diagram for G4PhysicsLinearVector:
G4PhysicsVector

Public Member Functions

 G4PhysicsLinearVector ()
 
 G4PhysicsLinearVector (size_t theNbin)
 
 G4PhysicsLinearVector (G4double theEmin, G4double theEmax, size_t theNbin)
 
virtual ~G4PhysicsLinearVector ()
 
virtual G4bool Retrieve (std::ifstream &fIn, G4bool ascii)
 
virtual void ScaleVector (G4double factorE, G4double factorV)
 
- Public Member Functions inherited from G4PhysicsVector
 G4PhysicsVector (G4bool spline=false)
 
 G4PhysicsVector (const G4PhysicsVector &)
 
G4PhysicsVectoroperator= (const G4PhysicsVector &)
 
virtual ~G4PhysicsVector ()
 
voidoperator new (size_t)
 
void operator delete (void *)
 
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)
 
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)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel (G4int)
 

Additional Inherited Members

- Protected Member Functions inherited from G4PhysicsVector
void DeleteData ()
 
void CopyData (const G4PhysicsVector &vec)
 
- Protected Attributes inherited from G4PhysicsVector
G4PhysicsVectorType type
 
G4double edgeMin
 
G4double edgeMax
 
size_t numberOfNodes
 
G4PVDataVector dataVector
 
G4PVDataVector binVector
 
G4PVDataVector secDerivative
 
G4double dBin
 
G4double baseBin
 
G4int verboseLevel
 

Detailed Description

Definition at line 60 of file G4PhysicsLinearVector.hh.

Constructor & Destructor Documentation

G4PhysicsLinearVector::G4PhysicsLinearVector ( )

Definition at line 43 of file G4PhysicsLinearVector.cc.

References T_G4PhysicsLinearVector, and G4PhysicsVector::type.

44  : G4PhysicsVector()
45 {
47 }
G4PhysicsVector(G4bool spline=false)
G4PhysicsVectorType type
G4PhysicsLinearVector::G4PhysicsLinearVector ( size_t  theNbin)
explicit

Definition at line 49 of file G4PhysicsLinearVector.cc.

References G4PhysicsVector::binVector, G4PhysicsVector::dataVector, G4PhysicsVector::numberOfNodes, T_G4PhysicsLinearVector, and G4PhysicsVector::type.

50  : G4PhysicsVector()
51 {
53 
54  numberOfNodes = theNbin + 1;
55  dataVector.reserve(numberOfNodes);
56  binVector.reserve(numberOfNodes);
57 
58  for (size_t i=0; i<numberOfNodes; i++)
59  {
60  binVector.push_back(0.0);
61  dataVector.push_back(0.0);
62  }
63 }
G4PVDataVector dataVector
G4PhysicsVector(G4bool spline=false)
G4PVDataVector binVector
G4PhysicsVectorType type
G4PhysicsLinearVector::G4PhysicsLinearVector ( G4double  theEmin,
G4double  theEmax,
size_t  theNbin 
)

Definition at line 65 of file G4PhysicsLinearVector.cc.

References G4PhysicsVector::baseBin, G4PhysicsVector::binVector, G4PhysicsVector::dataVector, G4PhysicsVector::dBin, G4PhysicsVector::edgeMax, G4PhysicsVector::edgeMin, G4PhysicsVector::numberOfNodes, T_G4PhysicsLinearVector, and G4PhysicsVector::type.

67  : G4PhysicsVector()
68 {
70 
71  dBin = (theEmax-theEmin)/theNbin;
72  baseBin = theEmin/dBin;
73 
74  numberOfNodes = theNbin + 1;
75  dataVector.reserve(numberOfNodes);
76  binVector.reserve(numberOfNodes);
77 
78  binVector.push_back(theEmin);
79  dataVector.push_back(0.0);
80 
81  for (size_t i=1; i<numberOfNodes-1; i++)
82  {
83  binVector.push_back( theEmin + i*dBin );
84  dataVector.push_back(0.0);
85  }
86  binVector.push_back(theEmax);
87  dataVector.push_back(0.0);
88 
89  edgeMin = binVector[0];
90  edgeMax = binVector[numberOfNodes-1];
91 
92 }
G4PVDataVector dataVector
G4PhysicsVector(G4bool spline=false)
G4PVDataVector binVector
G4PhysicsVectorType type
G4PhysicsLinearVector::~G4PhysicsLinearVector ( )
virtual

Definition at line 94 of file G4PhysicsLinearVector.cc.

94 {}

Member Function Documentation

G4bool G4PhysicsLinearVector::Retrieve ( std::ifstream &  fIn,
G4bool  ascii 
)
virtual

Reimplemented from G4PhysicsVector.

Definition at line 96 of file G4PhysicsLinearVector.cc.

References G4PhysicsVector::baseBin, G4PhysicsVector::binVector, G4PhysicsVector::dBin, and G4PhysicsVector::Retrieve().

97 {
98  G4bool success = G4PhysicsVector::Retrieve(fIn, ascii);
99  if (success)
100  {
101  G4double theEmin = binVector[0];
102  dBin = binVector[1]-theEmin;
103  baseBin = theEmin/dBin;
104  }
105  return success;
106 }
G4PVDataVector binVector
bool G4bool
Definition: G4Types.hh:79
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
double G4double
Definition: G4Types.hh:76
void G4PhysicsLinearVector::ScaleVector ( G4double  factorE,
G4double  factorV 
)
virtual

Reimplemented from G4PhysicsVector.

Definition at line 108 of file G4PhysicsLinearVector.cc.

References G4PhysicsVector::baseBin, G4PhysicsVector::binVector, G4PhysicsVector::dBin, and G4PhysicsVector::ScaleVector().

109 {
110  G4PhysicsVector::ScaleVector(factorE, factorV);
111  G4double theEmin = binVector[0];
112  dBin = binVector[1]-theEmin;
113  baseBin = theEmin/dBin;
114 }
G4PVDataVector binVector
virtual void ScaleVector(G4double factorE, G4double factorV)
double G4double
Definition: G4Types.hh:76

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