Geant4 11.2.2
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 EnableLogBinSearch (const G4int n=1)
 
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 LogFreeVectorValue (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
 
G4double iBin1 = 0.0
 
G4double lmin1 = 0.0
 
G4int verboseLevel = 0
 
std::size_t idxmax = 0
 
std::size_t imax1 = 0
 
std::size_t numberOfNodes = 0
 
std::size_t nLogNodes = 0
 
G4PhysicsVectorType type = T_G4PhysicsFreeVector
 
std::vector< G4doublebinVector
 
std::vector< G4doubledataVector
 
std::vector< G4doublesecDerivative
 
std::vector< std::size_t > scale
 

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 40 of file G4PhysicsFreeVector.cc.

41 : G4PhysicsVector(spline)
42{}
G4PhysicsVector(G4bool spline=false)

◆ G4PhysicsFreeVector() [2/6]

G4PhysicsFreeVector::G4PhysicsFreeVector ( G4int length)
explicit

Definition at line 45 of file G4PhysicsFreeVector.cc.

46 : G4PhysicsFreeVector(static_cast<std::size_t>(length), false)
47{}
G4PhysicsFreeVector(G4bool spline=false)

◆ G4PhysicsFreeVector() [3/6]

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

Definition at line 50 of file G4PhysicsFreeVector.cc.

51 : G4PhysicsVector(spline)
52{
53 numberOfNodes = length;
54
55 if (0 < length) {
56 binVector.resize(numberOfNodes, 0.0);
57 dataVector.resize(numberOfNodes, 0.0);
58 }
59 Initialise();
60}
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 63 of file G4PhysicsFreeVector.cc.

65 : G4PhysicsFreeVector(length, spline)
66{}

◆ G4PhysicsFreeVector() [5/6]

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

Definition at line 69 of file G4PhysicsFreeVector.cc.

72 : G4PhysicsVector(spline)
73{
74 numberOfNodes = energies.size();
75
76 if (numberOfNodes != values.size())
77 {
79 ed << "The size of energy vector " << numberOfNodes << " != " << values.size();
80 G4Exception("G4PhysicsFreeVector constructor: ","glob04", FatalException, ed);
81 }
82
83 binVector = energies;
84 dataVector = values;
85 Initialise();
86}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription

◆ G4PhysicsFreeVector() [6/6]

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

Definition at line 89 of file G4PhysicsFreeVector.cc.

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

◆ ~G4PhysicsFreeVector()

G4PhysicsFreeVector::~G4PhysicsFreeVector ( )
overridedefault

Member Function Documentation

◆ EnableLogBinSearch()

void G4PhysicsFreeVector::EnableLogBinSearch ( const G4int n = 1)

Definition at line 149 of file G4PhysicsFreeVector.cc.

150{
151 // check if log search is applicable
152 if (0 >= n || edgeMin <= 0.0 || edgeMin == edgeMax || numberOfNodes < 3)
153 {
154 return;
155 }
156 nLogNodes = static_cast<std::size_t>(static_cast<G4int>(numberOfNodes)/n);
157 if (nLogNodes < 3) { nLogNodes = 3; }
158 scale.resize(nLogNodes, 0);
159 imax1 = nLogNodes - 2;
160 iBin1 = (imax1 + 1) / G4Log(edgeMax/edgeMin);
162 scale[0] = 0;
163 scale[imax1 + 1] = idxmax;
164 std::size_t j = 0;
165 for (std::size_t i = 1; i <= imax1; ++i)
166 {
168 for (; j <= idxmax; ++j)
169 {
170 if (binVector[j] <= e && e < binVector[j+1])
171 {
172 scale[i] = j;
173 break;
174 }
175 }
176 }
177}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
std::size_t nLogNodes
std::vector< std::size_t > scale

◆ 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 112 of file G4PhysicsFreeVector.cc.

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

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


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