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

#include <G4PhysicsVector.hh>

Inheritance diagram for G4PhysicsVector:

Public Member Functions

 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
G4bool CheckIndex (const G4double energy, std::size_t &lastidx) 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)
void SetDataLength (G4int dlength)
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

virtual void Initialise ()
void PrintPutValueError (std::size_t index, G4double value, const G4String &text)

Protected Attributes

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

Friends

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

Detailed Description

Definition at line 54 of file G4PhysicsVector.hh.

Constructor & Destructor Documentation

◆ G4PhysicsVector() [1/3]

◆ G4PhysicsVector() [2/3]

G4PhysicsVector::G4PhysicsVector ( const G4PhysicsVector & )
default

◆ G4PhysicsVector() [3/3]

G4PhysicsVector::G4PhysicsVector ( const G4PhysicsVector && )
delete

◆ ~G4PhysicsVector()

virtual G4PhysicsVector::~G4PhysicsVector ( )
virtualdefault

Member Function Documentation

◆ CheckIndex()

G4bool G4PhysicsVector::CheckIndex ( const G4double energy,
std::size_t & lastidx ) const
inline

◆ ComputeLogVectorBin()

std::size_t G4PhysicsVector::ComputeLogVectorBin ( const G4double logenergy) const
inline

◆ DumpValues()

void G4PhysicsVector::DumpValues ( G4double unitE = 1.0,
G4double unitV = 1.0 ) const

Definition at line 186 of file G4PhysicsVector.cc.

187{
188 for (std::size_t i = 0; i < numberOfNodes; ++i)
189 {
190 G4cout << binVector[i] / unitE << " " << dataVector[i] / unitV
191 << G4endl;
192 }
193}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::size_t numberOfNodes
std::vector< G4double > dataVector
std::vector< G4double > binVector

Referenced by G4OpWLS2::DumpPhysicsTable(), G4OpWLS::DumpPhysicsTable(), and FillSecondDerivatives().

◆ Energy()

G4double G4PhysicsVector::Energy ( const std::size_t index) const
inline

Referenced by G4MaterialPropertiesTable::AddEntry(), G4MaterialPropertiesTable::AddProperty(), G4LossTableBuilder::BuildInverseRangeTable(), G4Cerenkov::BuildPhysicsTable(), G4OpWLS2::BuildPhysicsTable(), G4OpWLS::BuildPhysicsTable(), G4QuasiCerenkov::BuildPhysicsTable(), G4LossTableBuilder::BuildRangeTable(), G4PenelopeBremsstrahlungFS::BuildScaledXSTable(), G4LossTableBuilder::BuildTableForModel(), G4AdjointCSManager::BuildTotalSigmaTables(), G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom(), G4JAEAPolarizedElasticScatteringModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(), G4LivermoreRayleighModel::ComputeCrossSectionPerAtom(), G4LowEPComptonModel::ComputeCrossSectionPerAtom(), G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(), G4EmModelManager::DumpModelList(), G4EmModelManager::FillDEDXVector(), G4EmModelManager::FillLambdaVector(), G4EmUtility::FillPeaksStructure(), G4EmUtility::FindCrossSectionMax(), G4WentzelVIModel::Initialise(), G4StandardCerenkovModel::InitialiseModel(), G4Cerenkov::PostStepDoIt(), G4QuasiCerenkov::PostStepDoIt(), G4GDMLWriteMaterials::PropertyVectorWrite(), G4PAIPhotData::SampleAlongStepPhotonTransfer(), G4PAIPhotData::SampleAlongStepPlasmonTransfer(), G4PAIModelData::SampleAlongStepTransfer(), G4PAIPhotData::SampleAlongStepTransfer(), G4PAIModelData::SamplePostStepTransfer(), and G4VChannelingFastSimCrystalData::SetCrystallineUndulatorParameters().

◆ FillSecondDerivatives()

void G4PhysicsVector::FillSecondDerivatives ( const G4SplineType stype = G4SplineType::Base,
const G4double dir1 = 0.0,
const G4double dir2 = 0.0 )

Definition at line 228 of file G4PhysicsVector.cc.

231{
232 if (!useSpline) { return; }
233 // cannot compute derivatives for less than 5 points
234 const std::size_t nmin = (stype == G4SplineType::Base) ? 5 : 4;
235 if (nmin > numberOfNodes)
236 {
237 if (0 < verboseLevel)
238 {
239 G4cout << "### G4PhysicsVector: spline cannot be used for "
240 << numberOfNodes << " points - spline disabled"
241 << G4endl;
242 DumpValues();
243 }
244 useSpline = false;
245 return;
246 }
247 // check energies of free vector
249 {
250 for (std::size_t i=0; i<=idxmax; ++i)
251 {
252 if (binVector[i + 1] <= binVector[i])
253 {
254 if (0 < verboseLevel)
255 {
256 G4cout << "### G4PhysicsVector: spline cannot be used, because "
257 << " E[" << i << "]=" << binVector[i]
258 << " >= E[" << i+1 << "]=" << binVector[i + 1] << G4endl;
259 DumpValues();
260 }
261 useSpline = false;
262 return;
263 }
264 }
265 }
266
267 // spline is possible
268 Initialise();
270
271 if (1 < verboseLevel)
272 {
273 G4cout << "### G4PhysicsVector:: FillSecondDerivatives N="
274 << numberOfNodes << G4endl;
275 DumpValues();
276 }
277
278 switch(stype)
279 {
281 ComputeSecDerivative1();
282 break;
283
285 ComputeSecDerivative2(dir1, dir2);
286 break;
287
288 default:
289 ComputeSecDerivative0();
290 }
291}
@ T_G4PhysicsFreeVector
G4PhysicsVectorType type
std::vector< G4double > secDerivative
virtual void Initialise()
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const

Referenced by G4EmTableUtil::BuildDEDXTable(), G4EmTableUtil::BuildLambdaTable(), G4EmTableUtil::BuildLambdaTable(), G4LossTableBuilder::BuildRangeTable(), G4LossTableBuilder::BuildTableForModel(), G4WentzelVIModel::Initialise(), G4PenelopeBremsstrahlungAngular::PrepareTables(), G4ExtDEDXTable::RetrievePhysicsTable(), and G4VChannelingFastSimCrystalData::SetCrystallineUndulatorParameters().

◆ FindBin()

std::size_t G4PhysicsVector::FindBin ( const G4double energy,
std::size_t idx ) const

Definition at line 196 of file G4PhysicsVector.cc.

198{
199 if (idx + 1 < numberOfNodes &&
200 energy >= binVector[idx] && energy <= binVector[idx])
201 {
202 return idx;
203 }
204 if (energy <= binVector[1])
205 {
206 return 0;
207 }
208 if (energy >= binVector[idxmax])
209 {
210 return idxmax;
211 }
212 return GetBin(energy);
213}

◆ FindLinearEnergy()

G4double G4PhysicsVector::FindLinearEnergy ( const G4double rand) const
inline

◆ GetEnergy()

G4double G4PhysicsVector::GetEnergy ( const G4double value) const

Definition at line 453 of file G4PhysicsVector.cc.

454{
455 if (0 == numberOfNodes)
456 {
457 return 0.0;
458 }
459 if (1 == numberOfNodes || val <= dataVector[0])
460 {
461 return edgeMin;
462 }
463 if (val >= dataVector[numberOfNodes - 1])
464 {
465 return edgeMax;
466 }
467 std::size_t bin = std::lower_bound(dataVector.cbegin(), dataVector.cend(), val)
468 - dataVector.cbegin() - 1;
469 if (bin > idxmax) { bin = idxmax; }
470 G4double res = binVector[bin];
471 G4double del = dataVector[bin + 1] - dataVector[bin];
472 if (del > 0.0)
473 {
474 res += (val - dataVector[bin]) * (binVector[bin + 1] - res) / del;
475 }
476 return res;
477}
double G4double
Definition G4Types.hh:83

Referenced by G4OpWLS2::PostStepDoIt(), and G4OpWLS::PostStepDoIt().

◆ GetLowEdgeEnergy()

◆ GetMaxEnergy()

◆ GetMaxValue()

◆ GetMinEnergy()

◆ GetMinValue()

G4double G4PhysicsVector::GetMinValue ( ) const
inline

◆ GetSpline()

G4bool G4PhysicsVector::GetSpline ( ) const
inline

◆ GetType()

G4PhysicsVectorType G4PhysicsVector::GetType ( ) const
inline

◆ GetValue()

G4double G4PhysicsVector::GetValue ( const G4double energy,
G4bool & isOutRange ) const
inline

◆ GetVectorLength()

std::size_t G4PhysicsVector::GetVectorLength ( ) const
inline

Referenced by G4MaterialPropertiesTable::AddEntry(), G4MaterialPropertiesTable::AddProperty(), G4LossTableBuilder::BuildDEDXTable(), G4LossTableBuilder::BuildInverseRangeTable(), G4Cerenkov::BuildPhysicsTable(), G4OpWLS2::BuildPhysicsTable(), G4OpWLS::BuildPhysicsTable(), G4QuasiCerenkov::BuildPhysicsTable(), G4LossTableBuilder::BuildRangeTable(), G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom(), G4JAEAPolarizedElasticScatteringModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModel::ComputeCrossSectionPerAtom(), G4LivermoreNuclearGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(), G4LivermoreRayleighModel::ComputeCrossSectionPerAtom(), G4LowEPComptonModel::ComputeCrossSectionPerAtom(), G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(), G4PenelopeRayleighModel::DumpFormFactorTable(), G4PenelopeRayleighModelMI::DumpFormFactorTable(), G4EmModelManager::DumpModelList(), G4EmModelManager::FillDEDXVector(), G4EmModelManager::FillLambdaVector(), G4EmUtility::FillPeaksStructure(), G4EmUtility::FindCrossSectionMax(), G4Cerenkov::GetAverageNumberOfPhotons(), G4QuasiCerenkov::GetAverageNumberOfPhotons(), G4PenelopeCrossSection::GetHardCrossSection(), G4PenelopeCrossSection::GetNormalizedShellCrossSection(), G4PenelopeCrossSection::GetShellCrossSection(), G4PenelopeCrossSection::GetSoftStoppingPower(), G4PenelopeCrossSection::GetTotalCrossSection(), G4StandardCerenkovModel::InitialiseModel(), G4GDMLWriteMaterials::PropertyVectorWrite(), and G4VChannelingFastSimCrystalData::SetCrystallineUndulatorParameters().

◆ Initialise()

◆ LogFreeVectorValue()

G4double G4PhysicsVector::LogFreeVectorValue ( const G4double energy,
const G4double theLogEnergy ) const
inline

◆ LogVectorValue()

G4double G4PhysicsVector::LogVectorValue ( const G4double energy,
const G4double theLogEnergy ) const
inline

◆ operator!=()

G4bool G4PhysicsVector::operator!= ( const G4PhysicsVector & right) const
delete

◆ operator()()

G4double G4PhysicsVector::operator() ( const std::size_t index) const
inline

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ operator==()

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

◆ operator[]()

G4double G4PhysicsVector::operator[] ( const std::size_t index) const
inline

◆ PrintPutValueError()

void G4PhysicsVector::PrintPutValueError ( std::size_t index,
G4double value,
const G4String & text )
protected

Definition at line 480 of file G4PhysicsVector.cc.

483{
485 ed << "Vector type: " << type << " length= " << numberOfNodes
486 << "; an attempt to put data at index= " << index
487 << " value= " << val << " in " << text;
488 G4Exception("G4PhysicsVector:", "gl0005",
489 FatalException, ed, "Wrong operation");
490}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription

Referenced by G4PhysicsFreeVector::PutValues().

◆ PutValue()

◆ Retrieve()

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

Definition at line 108 of file G4PhysicsVector.cc.

109{
110 // clear properties;
111 dataVector.clear();
112 binVector.clear();
113 secDerivative.clear();
114
115 // retrieve in ascii mode
116 if (ascii)
117 {
118 // binning
119 fIn >> edgeMin >> edgeMax >> numberOfNodes;
120 if (fIn.fail() || numberOfNodes < 2)
121 {
122 return false;
123 }
124 // contents
125 G4int siz0 = 0;
126 fIn >> siz0;
127 if (siz0 < 2) { return false; }
128 auto siz = static_cast<std::size_t>(siz0);
129 if (fIn.fail() || siz != numberOfNodes)
130 {
131 return false;
132 }
133
134 binVector.reserve(siz);
135 dataVector.reserve(siz);
136 G4double vBin, vData;
137
138 for (std::size_t i = 0; i < siz; ++i)
139 {
140 vBin = 0.;
141 vData = 0.;
142 fIn >> vBin >> vData;
143 if (fIn.fail())
144 {
145 return false;
146 }
147 binVector.push_back(vBin);
148 dataVector.push_back(vData);
149 }
150 Initialise();
151 return true;
152 }
153
154 // retrieve in binary mode
155 // binning
156 fIn.read((char*) (&edgeMin), sizeof edgeMin);
157 fIn.read((char*) (&edgeMax), sizeof edgeMax);
158 fIn.read((char*) (&numberOfNodes), sizeof numberOfNodes);
159
160 // contents
161 std::size_t size;
162 fIn.read((char*) (&size), sizeof size);
163
164 auto value = new G4double[2 * size];
165 fIn.read((char*) (value), 2 * size * (sizeof(G4double)));
166 if (static_cast<G4int>(fIn.gcount()) != static_cast<G4int>(2 * size * (sizeof(G4double))))
167 {
168 delete[] value;
169 return false;
170 }
171
172 binVector.reserve(size);
173 dataVector.reserve(size);
174 for (std::size_t i = 0; i < size; ++i)
175 {
176 binVector.push_back(value[2 * i]);
177 dataVector.push_back(value[2 * i + 1]);
178 }
179 delete[] value;
180
181 Initialise();
182 return true;
183}
int G4int
Definition G4Types.hh:85

Referenced by G4ExtDEDXTable::RetrievePhysicsTable(), G4PhysicsTable::RetrievePhysicsTable(), and G4VChannelingFastSimCrystalData::SetCrystallineUndulatorParameters().

◆ ScaleVector()

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

Definition at line 216 of file G4PhysicsVector.cc.

218{
219 for (std::size_t i = 0; i < numberOfNodes; ++i)
220 {
221 binVector[i] *= factorE;
222 dataVector[i] *= factorV;
223 }
224 Initialise();
225}

◆ SetDataLength()

void G4PhysicsVector::SetDataLength ( G4int dlength)

Definition at line 55 of file G4PhysicsVector.cc.

56{
57 // this method may be applied for empty vector only
58 if (numberOfNodes > 0) { return; }
59
60 //
61 if (dlength < 1)
62 {
63 if (0 < verboseLevel)
64 {
65 G4cout << "### G4PhysicsVector::SetDataLength length="
66 << dlength << " is ignored." << G4endl;
67 }
68 return;
69 }
70 numberOfNodes = dlength;
71 binVector.resize(numberOfNodes, 0.0);
72 dataVector.resize(numberOfNodes, 0.0);
73}

◆ SetVerboseLevel()

void G4PhysicsVector::SetVerboseLevel ( G4int value)
inline

◆ Store()

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

Definition at line 76 of file G4PhysicsVector.cc.

77{
78 // Ascii mode
79 if (ascii)
80 {
81 fOut << *this;
82 return true;
83 }
84 // Binary Mode
85
86 // binning
87 fOut.write((char*) (&edgeMin), sizeof edgeMin);
88 fOut.write((char*) (&edgeMax), sizeof edgeMax);
89 fOut.write((char*) (&numberOfNodes), sizeof numberOfNodes);
90
91 // contents
92 std::size_t size = dataVector.size();
93 fOut.write((char*) (&size), sizeof size);
94
95 auto value = new G4double[2 * size];
96 for (std::size_t i = 0; i < size; ++i)
97 {
98 value[2 * i] = binVector[i];
99 value[2 * i + 1] = dataVector[i];
100 }
101 fOut.write((char*) (value), 2 * size * (sizeof(G4double)));
102 delete[] value;
103
104 return true;
105}

Referenced by G4ExtDEDXTable::StorePhysicsTable().

◆ Value() [1/2]

G4double G4PhysicsVector::Value ( const G4double energy) const
inline

◆ Value() [2/2]

G4double G4PhysicsVector::Value ( const G4double energy,
std::size_t & lastidx ) const
inline

Referenced by G4BoldyshevTripletModel::ComputeCrossSectionPerAtom(), G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom(), G4JAEAPolarizedElasticScatteringModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModel::ComputeCrossSectionPerAtom(), G4LivermoreGammaConversion5DModel::ComputeCrossSectionPerAtom(), G4LivermoreGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermoreNuclearGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(), G4LivermoreRayleighModel::ComputeCrossSectionPerAtom(), G4LowEPComptonModel::ComputeCrossSectionPerAtom(), G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(), G4PenelopeGammaConversionModel::ComputeCrossSectionPerAtom(), G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom(), G4PenelopeRayleighModel::ComputeCrossSectionPerAtom(), G4PenelopeRayleighModelMI::ComputeCrossSectionPerAtom(), G4Cerenkov::GetAverageNumberOfPhotons(), G4QuasiCerenkov::GetAverageNumberOfPhotons(), G4PenelopeIonisationXSHandler::GetDensityCorrection(), G4PenelopeCrossSection::GetHardCrossSection(), G4OpAbsorption::GetMeanFreePath(), G4OpMieHG::GetMeanFreePath(), G4OpWLS2::GetMeanFreePath(), G4OpWLS::GetMeanFreePath(), G4PenelopeCrossSection::GetNormalizedShellCrossSection(), G4QuasiScintillation::GetScintillationYieldByParticleType(), G4Scintillation::GetScintillationYieldByParticleType(), G4PenelopeCrossSection::GetShellCrossSection(), G4PenelopePhotoElectricModel::GetShellCrossSection(), G4PenelopeCrossSection::GetSoftStoppingPower(), G4PenelopeCrossSection::GetTotalCrossSection(), G4Cerenkov::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4QuasiCerenkov::PostStepDoIt(), G4QuasiScintillation::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4PenelopeBremsstrahlungAngular::PrepareTables(), G4XrayReflection::Reflectivity(), G4PAIModelData::SampleAlongStepTransfer(), G4PenelopeBremsstrahlungAngular::SampleDirection(), G4PAIModelData::SamplePostStepTransfer(), G4PenelopeRayleighModel::SampleSecondaries(), and G4PenelopeRayleighModelMI::SampleSecondaries().

◆ operator<<

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

Definition at line 434 of file G4PhysicsVector.cc.

435{
436 // binning
437 G4long prec = out.precision();
438 out << std::setprecision(12) << pv.edgeMin << " " << pv.edgeMax << " "
439 << pv.numberOfNodes << G4endl;
440
441 // contents
442 out << pv.dataVector.size() << G4endl;
443 for (std::size_t i = 0; i < pv.dataVector.size(); ++i)
444 {
445 out << pv.binVector[i] << " " << pv.dataVector[i] << G4endl;
446 }
447 out.precision(prec);
448
449 return out;
450}
long G4long
Definition G4Types.hh:87

Member Data Documentation

◆ binVector

◆ dataVector

◆ edgeMax

◆ edgeMin

◆ iBin1

G4double G4PhysicsVector::iBin1 = 0.0
protected

Definition at line 224 of file G4PhysicsVector.hh.

Referenced by G4PhysicsFreeVector::EnableLogBinSearch().

◆ idxmax

◆ imax1

std::size_t G4PhysicsVector::imax1 = 0
protected

Definition at line 229 of file G4PhysicsVector.hh.

Referenced by G4PhysicsFreeVector::EnableLogBinSearch().

◆ invdBin

◆ lmin1

G4double G4PhysicsVector::lmin1 = 0.0
protected

Definition at line 225 of file G4PhysicsVector.hh.

Referenced by G4PhysicsFreeVector::EnableLogBinSearch().

◆ logemin

G4double G4PhysicsVector::logemin = 0.0
protected

Definition at line 222 of file G4PhysicsVector.hh.

Referenced by G4PhysicsLogVector::Initialise().

◆ nLogNodes

std::size_t G4PhysicsVector::nLogNodes = 0
protected

Definition at line 231 of file G4PhysicsVector.hh.

Referenced by G4PhysicsFreeVector::EnableLogBinSearch().

◆ numberOfNodes

◆ scale

std::vector<std::size_t> G4PhysicsVector::scale
protected

Definition at line 239 of file G4PhysicsVector.hh.

Referenced by G4PhysicsFreeVector::EnableLogBinSearch().

◆ secDerivative

std::vector<G4double> G4PhysicsVector::secDerivative
protected

Definition at line 238 of file G4PhysicsVector.hh.

Referenced by FillSecondDerivatives(), and Retrieve().

◆ type

◆ verboseLevel

G4int G4PhysicsVector::verboseLevel = 0
protected

Definition at line 227 of file G4PhysicsVector.hh.

Referenced by FillSecondDerivatives(), and SetDataLength().


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