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

#include <G4PhysicsLinearVector.hh>

+ Inheritance diagram for G4PhysicsLinearVector:

Public Member Functions

 G4PhysicsLinearVector (G4bool spline=false)
 
 G4PhysicsLinearVector (G4double Emin, G4double Emax, std::size_t Nbin, G4bool spline=false)
 
 ~G4PhysicsLinearVector () override=default
 
- Public Member Functions inherited from G4PhysicsVector
 G4PhysicsVector (G4bool spline=false)
 
 G4PhysicsVector (const G4PhysicsVector &)=default
 
G4PhysicsVectoroperator= (const G4PhysicsVector &)=default
 
 G4PhysicsVector (const G4PhysicsVector &&)=delete
 
G4PhysicsVectoroperator= (const G4PhysicsVector &&)=delete
 
G4bool operator== (const G4PhysicsVector &right) const =delete
 
G4bool operator!= (const G4PhysicsVector &right) const =delete
 
virtual ~G4PhysicsVector ()=default
 
G4double Value (const G4double energy, std::size_t &lastidx) const
 
G4double Value (const G4double energy) const
 
G4double GetValue (const G4double energy, G4bool &isOutRange) const
 
G4double LogVectorValue (const G4double energy, const G4double theLogEnergy) const
 
G4double operator[] (const std::size_t index) const
 
G4double operator() (const std::size_t index) const
 
void PutValue (const std::size_t index, const G4double value)
 
G4double Energy (const std::size_t index) const
 
G4double GetLowEdgeEnergy (const std::size_t index) const
 
G4double GetMinEnergy () const
 
G4double GetMaxEnergy () const
 
G4double GetMinValue () const
 
G4double GetMaxValue () const
 
std::size_t GetVectorLength () const
 
std::size_t ComputeLogVectorBin (const G4double logenergy) const
 
G4PhysicsVectorType GetType () const
 
G4bool GetSpline () const
 
void SetVerboseLevel (G4int value)
 
G4double FindLinearEnergy (const G4double rand) const
 
std::size_t FindBin (const G4double energy, std::size_t idx) const
 
void ScaleVector (const G4double factorE, const G4double factorV)
 
void FillSecondDerivatives (const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
 
G4double GetEnergy (const G4double value) const
 
G4bool Store (std::ofstream &fOut, G4bool ascii=false) const
 
G4bool Retrieve (std::ifstream &fIn, G4bool ascii=false)
 
void DumpValues (G4double unitE=1.0, G4double unitV=1.0) const
 

Protected Member Functions

void Initialise () final
 
- Protected Member Functions inherited from G4PhysicsVector
virtual void Initialise ()
 
void PrintPutValueError (std::size_t index, G4double value, const G4String &text)
 

Additional Inherited Members

- Protected Attributes inherited from G4PhysicsVector
G4double edgeMin = 0.0
 
G4double edgeMax = 0.0
 
G4double invdBin = 0.0
 
G4double logemin = 0.0
 
G4int verboseLevel = 0
 
std::size_t idxmax = 0
 
std::size_t numberOfNodes = 0
 
G4PhysicsVectorType type = T_G4PhysicsFreeVector
 
std::vector< G4doublebinVector
 
std::vector< G4doubledataVector
 
std::vector< G4doublesecDerivative
 

Detailed Description

Definition at line 47 of file G4PhysicsLinearVector.hh.

Constructor & Destructor Documentation

◆ G4PhysicsLinearVector() [1/2]

G4PhysicsLinearVector::G4PhysicsLinearVector ( G4bool  spline = false)
explicit

Definition at line 38 of file G4PhysicsLinearVector.cc.

39 : G4PhysicsVector(spline)
40{
42}
@ T_G4PhysicsLinearVector
G4PhysicsVectorType type

◆ G4PhysicsLinearVector() [2/2]

G4PhysicsLinearVector::G4PhysicsLinearVector ( G4double  Emin,
G4double  Emax,
std::size_t  Nbin,
G4bool  spline = false 
)
explicit

Definition at line 45 of file G4PhysicsLinearVector.cc.

47 : G4PhysicsVector(spline)
48{
49 numberOfNodes = Nbin + 1;
50 if(Nbin < 1 || Emin >= Emax)
51 {
53 ed << "G4PhysicsLinearVector with wrong parameters: theNbin= " << Nbin
54 << " Emin= " << Emin << " Emax= " << Emax;
55 G4Exception("G4PhysicsLinearVector::G4PhysicsLinearVector()", "glob03",
56 FatalException, ed, "theNbins should be > 0 and Emax > Emin");
57 }
58 if(numberOfNodes < 2)
59 {
60 numberOfNodes = 2;
61 }
63
65 dataVector.resize(numberOfNodes, 0.0);
66 binVector[0] = Emin;
67 binVector[numberOfNodes - 1] = Emax;
68 Initialise();
69
70 for(std::size_t i = 1; i <= idxmax; ++i)
71 {
72 binVector[i] = edgeMin + i / invdBin;
73 }
74}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::size_t idxmax
std::size_t numberOfNodes
std::vector< G4double > dataVector
std::vector< G4double > binVector

◆ ~G4PhysicsLinearVector()

G4PhysicsLinearVector::~G4PhysicsLinearVector ( )
overridedefault

Member Function Documentation

◆ Initialise()

void G4PhysicsLinearVector::Initialise ( )
finalprotectedvirtual

Reimplemented from G4PhysicsVector.

Definition at line 77 of file G4PhysicsLinearVector.cc.

78{
80 edgeMin = binVector[0];
82 invdBin = (idxmax + 1) / (edgeMax - edgeMin);
83}

Referenced by G4PhysicsLinearVector().


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