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

#include <G4PhysicsFreeVector.hh>

+ Inheritance diagram for G4PhysicsFreeVector:

Public Member Functions

 G4PhysicsFreeVector (G4bool spline=false)
 
 G4PhysicsFreeVector (G4int length)
 
 G4PhysicsFreeVector (std::size_t length, G4bool spline=false)
 
 G4PhysicsFreeVector (std::size_t length, G4double emin, G4double emax, G4bool spline=false)
 
 G4PhysicsFreeVector (const std::vector< G4double > &energies, const std::vector< G4double > &values, G4bool spline=false)
 
 G4PhysicsFreeVector (const G4double *energies, const G4double *values, std::size_t length, G4bool spline=false)
 
 ~G4PhysicsFreeVector () override=default
 
void PutValues (const std::size_t index, const G4double energy, const G4double value)
 
void InsertValues (const G4double energy, const G4double value)
 
void PutValue (const std::size_t index, const G4double e, const G4double value)
 
- 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4PhysicsVector
virtual void Initialise ()
 
void PrintPutValueError (std::size_t index, G4double value, const G4String &text)
 
- 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 53 of file G4PhysicsFreeVector.hh.

Constructor & Destructor Documentation

◆ G4PhysicsFreeVector() [1/6]

G4PhysicsFreeVector::G4PhysicsFreeVector ( G4bool  spline = false)
explicit

Definition at line 39 of file G4PhysicsFreeVector.cc.

◆ G4PhysicsFreeVector() [2/6]

G4PhysicsFreeVector::G4PhysicsFreeVector ( G4int  length)
explicit

Definition at line 44 of file G4PhysicsFreeVector.cc.

45 : G4PhysicsFreeVector((std::size_t)length, false)
46{}

◆ G4PhysicsFreeVector() [3/6]

G4PhysicsFreeVector::G4PhysicsFreeVector ( std::size_t  length,
G4bool  spline = false 
)
explicit

Definition at line 49 of file G4PhysicsFreeVector.cc.

50 : G4PhysicsVector(spline)
51{
52 numberOfNodes = length;
53
54 if(0 < length) {
55 binVector.resize(numberOfNodes, 0.0);
56 dataVector.resize(numberOfNodes, 0.0);
57 }
58 Initialise();
59}
std::size_t numberOfNodes
std::vector< G4double > dataVector
std::vector< G4double > binVector
virtual void Initialise()

◆ G4PhysicsFreeVector() [4/6]

G4PhysicsFreeVector::G4PhysicsFreeVector ( std::size_t  length,
G4double  emin,
G4double  emax,
G4bool  spline = false 
)
explicit

Definition at line 62 of file G4PhysicsFreeVector.cc.

64 : G4PhysicsFreeVector(length, spline)
65{}

◆ G4PhysicsFreeVector() [5/6]

G4PhysicsFreeVector::G4PhysicsFreeVector ( const std::vector< G4double > &  energies,
const std::vector< G4double > &  values,
G4bool  spline = false 
)
explicit

Definition at line 68 of file G4PhysicsFreeVector.cc.

71 : G4PhysicsVector(spline)
72{
73 numberOfNodes = energies.size();
74
75 if(numberOfNodes != values.size())
76 {
78 ed << "The size of energy vector " << numberOfNodes << " != " << values.size();
79 G4Exception("G4PhysicsFreeVector constructor: ","glob04", FatalException, ed);
80 }
81
82 binVector = energies;
83 dataVector = values;
84 Initialise();
85}
@ 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

◆ G4PhysicsFreeVector() [6/6]

G4PhysicsFreeVector::G4PhysicsFreeVector ( const G4double energies,
const G4double values,
std::size_t  length,
G4bool  spline = false 
)
explicit

Definition at line 88 of file G4PhysicsFreeVector.cc.

92 : G4PhysicsVector(spline)
93{
94 numberOfNodes = length;
95
96 if(0 < numberOfNodes)
97 {
100
101 for(std::size_t i = 0; i < numberOfNodes; ++i)
102 {
103 binVector[i] = energies[i];
104 dataVector[i] = values[i];
105 }
106 }
107 Initialise();
108}

◆ ~G4PhysicsFreeVector()

G4PhysicsFreeVector::~G4PhysicsFreeVector ( )
overridedefault

Member Function Documentation

◆ InsertValues()

void G4PhysicsFreeVector::InsertValues ( const G4double  energy,
const G4double  value 
)

◆ PutValue()

◆ PutValues()

void G4PhysicsFreeVector::PutValues ( const std::size_t  index,
const G4double  energy,
const G4double  value 
)

Definition at line 111 of file G4PhysicsFreeVector.cc.

114{
115 if(index >= numberOfNodes)
116 {
117 PrintPutValueError(index, value, "G4PhysicsFreeVector::PutValues ");
118 return;
119 }
120 binVector[index] = e;
121 dataVector[index] = value;
122 if(index == 0)
123 {
124 edgeMin = e;
125 }
126 else if(numberOfNodes == index + 1)
127 {
128 edgeMax = e;
129 }
130}
void PrintPutValueError(std::size_t index, G4double value, const G4String &text)

Referenced by G4PenelopeCrossSection::AddCrossSectionPoint(), G4PenelopeCrossSection::AddShellCrossSectionPoint(), G4PenelopeBremsstrahlungFS::BuildScaledXSTable(), G4ParticleHPData::DoPhysicsVector(), G4PenelopeCrossSection::NormalizeShellCrossSections(), G4PenelopeBremsstrahlungAngular::PrepareTables(), PutValue(), and G4PenelopeBremsstrahlungFS::SampleGammaEnergy().


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