Geant4 11.1.1
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 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

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
 
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
 

Friends

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

Detailed Description

Definition at line 54 of file G4PhysicsVector.hh.

Constructor & Destructor Documentation

◆ G4PhysicsVector() [1/3]

G4PhysicsVector::G4PhysicsVector ( G4bool  spline = false)
explicit

Definition at line 39 of file G4PhysicsVector.cc.

40 : useSpline(val)
41{}

◆ G4PhysicsVector() [2/3]

G4PhysicsVector::G4PhysicsVector ( const G4PhysicsVector )
default

◆ G4PhysicsVector() [3/3]

G4PhysicsVector::G4PhysicsVector ( const G4PhysicsVector &&  )
delete

◆ ~G4PhysicsVector()

virtual G4PhysicsVector::~G4PhysicsVector ( )
virtualdefault

Member Function Documentation

◆ 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 163 of file G4PhysicsVector.cc.

164{
165 for(std::size_t i = 0; i < numberOfNodes; ++i)
166 {
167 G4cout << binVector[i] / unitE << " " << dataVector[i] / unitV
168 << G4endl;
169 }
170}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
std::size_t numberOfNodes
std::vector< G4double > dataVector
std::vector< G4double > binVector

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

◆ Energy()

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

Referenced by G4MaterialPropertiesTable::AddEntry(), G4MaterialPropertiesTable::AddProperty(), G4SPSEneDistribution::ArbInterpolate(), G4EmCorrections::BarkasCorrection(), G4DiffuseElasticV2::BuildAngleTable(), G4KokoulinMuonNuclearXS::BuildCrossSectionTable(), G4LossTableBuilder::BuildInverseRangeTable(), G4VLEPTSModel::BuildMeanFreePathTable(), G4Cerenkov::BuildPhysicsTable(), G4Scintillation::BuildPhysicsTable(), G4OpWLS::BuildPhysicsTable(), G4OpWLS2::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(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4PAIModelData::CrossSectionPerVolume(), G4PAIPhotData::CrossSectionPerVolume(), G4PAIModelData::DEDXPerVolume(), G4PAIPhotData::DEDXPerVolume(), G4EmModelManager::DumpModelList(), G4EmModelManager::FillDEDXVector(), G4EmModelManager::FillLambdaVector(), G4EmUtility::FillPeaksStructure(), G4EmUtility::FindCrossSectionMax(), G4HadElementSelector::G4HadElementSelector(), G4LivermorePhotoElectricModel::GetBindingEnergy(), G4PAIPhotData::GetPlasmonRatio(), G4PAIPhotData::Initialise(), G4PAIModelData::Initialise(), G4eeToHadronsModel::Initialise(), G4eplusTo2GammaOKVIModel::Initialise(), G4WentzelVIModel::Initialise(), G4DNABornExcitationModel2::Initialise(), G4HadronXSDataTable::Initialise(), G4NeutronCaptureXS::IsoCrossSection(), G4Cerenkov::PostStepDoIt(), G4GDMLWriteMaterials::PropertyVectorWrite(), G4PAIPhotData::SampleAlongStepPhotonTransfer(), G4PAIPhotData::SampleAlongStepPlasmonTransfer(), G4PAIPhotData::SampleAlongStepTransfer(), G4PAIModelData::SampleAlongStepTransfer(), G4PAIPhotData::SamplePostStepPhotonTransfer(), G4PAIPhotData::SamplePostStepPlasmonTransfer(), G4PAIPhotData::SamplePostStepTransfer(), G4PAIModelData::SamplePostStepTransfer(), and G4DiffuseElasticV2::SampleTableThetaCMS().

◆ FillSecondDerivatives()

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

Definition at line 205 of file G4PhysicsVector.cc.

208{
209 if(!useSpline) { return; }
210 // cannot compute derivatives for less than 5 points
211 const std::size_t nmin = (stype == G4SplineType::Base) ? 5 : 4;
212 if(nmin > numberOfNodes)
213 {
214 if(0 < verboseLevel)
215 {
216 G4cout << "### G4PhysicsVector: spline cannot be used for "
217 << numberOfNodes << " points - spline disabled"
218 << G4endl;
219 DumpValues();
220 }
221 useSpline = false;
222 return;
223 }
224 // check energies of free vector
226 {
227 for(std::size_t i=0; i<=idxmax; ++i)
228 {
229 if(binVector[i + 1] <= binVector[i])
230 {
231 if(0 < verboseLevel)
232 {
233 G4cout << "### G4PhysicsVector: spline cannot be used, because "
234 << " E[" << i << "]=" << binVector[i]
235 << " >= E[" << i+1 << "]=" << binVector[i + 1]
236 << G4endl;
237 DumpValues();
238 }
239 useSpline = false;
240 return;
241 }
242 }
243 }
244
245 // spline is possible
246 Initialise();
248
249 if(1 < verboseLevel)
250 {
251 G4cout << "### G4PhysicsVector:: FillSecondDerivatives N="
252 << numberOfNodes << G4endl;
253 DumpValues();
254 }
255
256 switch(stype)
257 {
258 case G4SplineType::Base:
259 ComputeSecDerivative1();
260 break;
261
262 case G4SplineType::FixedEdges:
263 ComputeSecDerivative2(dir1, dir2);
264 break;
265
266 default:
267 ComputeSecDerivative0();
268 }
269}
@ T_G4PhysicsFreeVector
std::size_t idxmax
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(), G4LossTableBuilder::BuildRangeTable(), G4LossTableBuilder::BuildTableForModel(), G4eplusTo2GammaOKVIModel::Initialise(), G4WentzelVIModel::Initialise(), G4HadronXSDataTable::Initialise(), and G4PenelopeBremsstrahlungAngular::PrepareTables().

◆ FindBin()

◆ FindLinearEnergy()

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

◆ GetEnergy()

G4double G4PhysicsVector::GetEnergy ( const G4double  value) const

Definition at line 431 of file G4PhysicsVector.cc.

432{
433 if(0 == numberOfNodes)
434 {
435 return 0.0;
436 }
437 if(1 == numberOfNodes || val <= dataVector[0])
438 {
439 return edgeMin;
440 }
441 if(val >= dataVector[numberOfNodes - 1])
442 {
443 return edgeMax;
444 }
445 std::size_t bin = std::lower_bound(dataVector.cbegin(), dataVector.cend(), val)
446 - dataVector.cbegin() - 1;
447 if(bin > idxmax) { bin = idxmax; }
448 G4double res = binVector[bin];
449 G4double del = dataVector[bin + 1] - dataVector[bin];
450 if(del > 0.0)
451 {
452 res += (val - dataVector[bin]) * (binVector[bin + 1] - res) / del;
453 }
454 return res;
455}
double G4double
Definition: G4Types.hh:83

Referenced by G4SPSRandomGenerator::GenRandEnergy(), G4SPSRandomGenerator::GenRandPhi(), G4SPSRandomGenerator::GenRandPosPhi(), G4SPSRandomGenerator::GenRandPosTheta(), G4SPSRandomGenerator::GenRandTheta(), G4SPSRandomGenerator::GenRandX(), G4SPSRandomGenerator::GenRandY(), G4SPSRandomGenerator::GenRandZ(), G4Scintillation::PostStepDoIt(), G4OpWLS::PostStepDoIt(), and G4OpWLS2::PostStepDoIt().

◆ GetLowEdgeEnergy()

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

◆ GetMaxEnergy()

◆ GetMaxValue()

◆ GetMinEnergy()

G4double G4PhysicsVector::GetMinEnergy ( ) const
inline

◆ GetMinValue()

G4double G4PhysicsVector::GetMinValue ( ) const
inline

◆ GetSpline()

G4bool G4PhysicsVector::GetSpline ( ) const
inline

◆ GetType()

G4PhysicsVectorType G4PhysicsVector::GetType ( ) const
inline

◆ GetValue()

◆ GetVectorLength()

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

Referenced by G4MaterialPropertiesTable::AddEntry(), G4MaterialPropertiesTable::AddProperty(), G4LossTableBuilder::BuildDEDXTable(), G4LossTableBuilder::BuildInverseRangeTable(), G4Cerenkov::BuildPhysicsTable(), G4Scintillation::BuildPhysicsTable(), G4OpWLS::BuildPhysicsTable(), G4OpWLS2::BuildPhysicsTable(), G4LossTableBuilder::BuildRangeTable(), G4PenelopeIonisationXSHandler::BuildXSTable(), G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom(), G4JAEAPolarizedElasticScatteringModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModel::ComputeCrossSectionPerAtom(), G4LivermoreNuclearGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(), G4LivermoreRayleighModel::ComputeCrossSectionPerAtom(), G4LowEPComptonModel::ComputeCrossSectionPerAtom(), G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(), G4PAIModelData::CrossSectionPerVolume(), G4PAIPhotData::CrossSectionPerVolume(), G4PAIModelData::DEDXPerVolume(), G4PAIPhotData::DEDXPerVolume(), G4PenelopeRayleighModel::DumpFormFactorTable(), G4PenelopeRayleighModelMI::DumpFormFactorTable(), G4EmModelManager::DumpModelList(), G4EmModelManager::FillDEDXVector(), G4EmModelManager::FillLambdaVector(), G4EmUtility::FillPeaksStructure(), G4EmUtility::FindCrossSectionMax(), G4SPSRandomGenerator::GenRandEnergy(), G4SPSRandomGenerator::GenRandPhi(), G4SPSRandomGenerator::GenRandPosPhi(), G4SPSRandomGenerator::GenRandPosTheta(), G4SPSRandomGenerator::GenRandTheta(), G4SPSRandomGenerator::GenRandX(), G4SPSRandomGenerator::GenRandY(), G4SPSRandomGenerator::GenRandZ(), G4SPSEneDistribution::GetArbEneWeight(), G4Cerenkov::GetAverageNumberOfPhotons(), G4PenelopeCrossSection::GetHardCrossSection(), G4PenelopeCrossSection::GetNormalizedShellCrossSection(), G4PAIPhotData::GetPlasmonRatio(), G4PenelopeCrossSection::GetShellCrossSection(), G4PenelopeCrossSection::GetSoftStoppingPower(), G4PenelopeCrossSection::GetTotalCrossSection(), G4eeToHadronsModel::Initialise(), G4GDMLWriteMaterials::PropertyVectorWrite(), G4PAIPhotData::SampleAlongStepPhotonTransfer(), G4PAIPhotData::SampleAlongStepPlasmonTransfer(), G4PAIPhotData::SampleAlongStepTransfer(), G4PAIModelData::SampleAlongStepTransfer(), G4PAIPhotData::SamplePostStepPhotonTransfer(), G4PAIPhotData::SamplePostStepPlasmonTransfer(), G4PAIPhotData::SamplePostStepTransfer(), and G4PAIModelData::SamplePostStepTransfer().

◆ Initialise()

void G4PhysicsVector::Initialise ( )
protectedvirtual

◆ LogVectorValue()

◆ 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 458 of file G4PhysicsVector.cc.

461{
463 ed << "Vector type: " << type << " length= " << numberOfNodes
464 << "; an attempt to put data at index= " << index
465 << " value= " << val << " in " << text;
466 G4Exception("G4PhysicsVector:", "gl0005",
467 FatalException, ed, "Wrong operation");
468}
@ 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

Referenced by G4PhysicsFreeVector::PutValues().

◆ PutValue()

◆ Retrieve()

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

Definition at line 87 of file G4PhysicsVector.cc.

88{
89 // clear properties;
90 dataVector.clear();
91 binVector.clear();
92 secDerivative.clear();
93
94 // retrieve in ascii mode
95 if(ascii)
96 {
97 // binning
98 fIn >> edgeMin >> edgeMax >> numberOfNodes;
99 if(fIn.fail() || numberOfNodes < 2)
100 {
101 return false;
102 }
103 // contents
104 G4int siz = 0;
105 fIn >> siz;
106 if(fIn.fail() || siz != G4int(numberOfNodes))
107 {
108 return false;
109 }
110
111 binVector.reserve(siz);
112 dataVector.reserve(siz);
113 G4double vBin, vData;
114
115 for(G4int i = 0; i < siz; ++i)
116 {
117 vBin = 0.;
118 vData = 0.;
119 fIn >> vBin >> vData;
120 if(fIn.fail())
121 {
122 return false;
123 }
124 binVector.push_back(vBin);
125 dataVector.push_back(vData);
126 }
127 Initialise();
128 return true;
129 }
130
131 // retrieve in binary mode
132 // binning
133 fIn.read((char*) (&edgeMin), sizeof edgeMin);
134 fIn.read((char*) (&edgeMax), sizeof edgeMax);
135 fIn.read((char*) (&numberOfNodes), sizeof numberOfNodes);
136
137 // contents
138 std::size_t size;
139 fIn.read((char*) (&size), sizeof size);
140
141 G4double* value = new G4double[2 * size];
142 fIn.read((char*) (value), 2 * size * (sizeof(G4double)));
143 if(G4int(fIn.gcount()) != G4int(2 * size * (sizeof(G4double))))
144 {
145 delete[] value;
146 return false;
147 }
148
149 binVector.reserve(size);
150 dataVector.reserve(size);
151 for(std::size_t i = 0; i < size; ++i)
152 {
153 binVector.push_back(value[2 * i]);
154 dataVector.push_back(value[2 * i + 1]);
155 }
156 delete[] value;
157
158 Initialise();
159 return true;
160}
int G4int
Definition: G4Types.hh:85

Referenced by G4PhysicsTable::RetrievePhysicsTable().

◆ ScaleVector()

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

Definition at line 193 of file G4PhysicsVector.cc.

195{
196 for(std::size_t i = 0; i < numberOfNodes; ++i)
197 {
198 binVector[i] *= factorE;
199 dataVector[i] *= factorV;
200 }
201 Initialise();
202}

◆ SetVerboseLevel()

void G4PhysicsVector::SetVerboseLevel ( G4int  value)
inline

◆ Store()

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

Definition at line 55 of file G4PhysicsVector.cc.

56{
57 // Ascii mode
58 if(ascii)
59 {
60 fOut << *this;
61 return true;
62 }
63 // Binary Mode
64
65 // binning
66 fOut.write((char*) (&edgeMin), sizeof edgeMin);
67 fOut.write((char*) (&edgeMax), sizeof edgeMax);
68 fOut.write((char*) (&numberOfNodes), sizeof numberOfNodes);
69
70 // contents
71 std::size_t size = dataVector.size();
72 fOut.write((char*) (&size), sizeof size);
73
74 G4double* value = new G4double[2 * size];
75 for(std::size_t i = 0; i < size; ++i)
76 {
77 value[2 * i] = binVector[i];
78 value[2 * i + 1] = dataVector[i];
79 }
80 fOut.write((char*) (value), 2 * size * (sizeof(G4double)));
81 delete[] value;
82
83 return true;
84}

◆ 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 G4EmCorrections::BarkasCorrection(), G4Track::CalculateVelocityForOpticalPhoton(), G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom(), G4JAEAPolarizedElasticScatteringModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModel::ComputeCrossSectionPerAtom(), G4LivermoreNuclearGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(), G4LivermoreRayleighModel::ComputeCrossSectionPerAtom(), G4LowEPComptonModel::ComputeCrossSectionPerAtom(), G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(), G4PenelopeGammaConversionModel::ComputeCrossSectionPerAtom(), G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom(), G4PenelopeRayleighModel::ComputeCrossSectionPerAtom(), G4PenelopeRayleighModelMI::ComputeCrossSectionPerAtom(), G4BoldyshevTripletModel::ComputeCrossSectionPerAtom(), G4eplusTo2GammaOKVIModel::ComputeCrossSectionPerAtom(), G4LivermoreGammaConversion5DModel::ComputeCrossSectionPerAtom(), G4LivermoreGammaConversionModel::ComputeCrossSectionPerAtom(), G4eeToHadronsModel::ComputeCrossSectionPerElectron(), G4eplusTo2GammaOKVIModel::CrossSectionPerVolume(), G4DNABornExcitationModel2::CrossSectionPerVolume(), G4EmCorrections::EffectiveChargeCorrection(), G4eeToHadronsModel::GenerateCMPhoton(), G4Cerenkov::GetAverageNumberOfPhotons(), G4ChannelingMaterialData::GetBR(), G4PenelopeIonisationXSHandler::GetDensityCorrection(), G4ChannelingECHARM::GetEC(), G4NeutronElectronElXsc::GetElementCrossSection(), G4KokoulinMuonNuclearXS::GetElementCrossSection(), G4PenelopeCrossSection::GetHardCrossSection(), G4OpAbsorption::GetMeanFreePath(), G4OpMieHG::GetMeanFreePath(), G4OpWLS::GetMeanFreePath(), G4OpWLS2::GetMeanFreePath(), G4PenelopeCrossSection::GetNormalizedShellCrossSection(), G4SPSEneDistribution::GetProbability(), G4Scintillation::GetScintillationYieldByParticleType(), G4PenelopePhotoElectricModel::GetShellCrossSection(), G4PenelopeCrossSection::GetShellCrossSection(), G4PenelopeCrossSection::GetSoftStoppingPower(), G4PenelopeCrossSection::GetTotalCrossSection(), G4ElementData::GetValueForElement(), G4NeutronElectronElXsc::Initialise(), G4eeToHadronsModel::Initialise(), G4EmCorrections::KShellCorrection(), G4EmCorrections::LShellCorrection(), G4Cerenkov::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4PenelopeBremsstrahlungAngular::PrepareTables(), G4PAIModelData::SampleAlongStepTransfer(), G4PenelopeBremsstrahlungAngular::SampleDirection(), G4PAIModelData::SamplePostStepTransfer(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4eplusTo2GammaOKVIModel::SampleSecondaries(), and G4EmCorrections::ShellCorrection().

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 412 of file G4PhysicsVector.cc.

413{
414 // binning
415 G4long prec = out.precision();
416 out << std::setprecision(12) << pv.edgeMin << " " << pv.edgeMax << " "
417 << pv.numberOfNodes << G4endl;
418
419 // contents
420 out << pv.dataVector.size() << G4endl;
421 for(std::size_t i = 0; i < pv.dataVector.size(); ++i)
422 {
423 out << pv.binVector[i] << " " << pv.dataVector[i] << G4endl;
424 }
425 out.precision(prec);
426
427 return out;
428}
long G4long
Definition: G4Types.hh:87

Member Data Documentation

◆ binVector

◆ dataVector

◆ edgeMax

◆ edgeMin

◆ idxmax

◆ invdBin

◆ logemin

G4double G4PhysicsVector::logemin = 0.0
protected

Definition at line 205 of file G4PhysicsVector.hh.

Referenced by G4PhysicsLogVector::Initialise().

◆ numberOfNodes

◆ secDerivative

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

Definition at line 216 of file G4PhysicsVector.hh.

Referenced by FillSecondDerivatives(), and Retrieve().

◆ type

◆ verboseLevel

G4int G4PhysicsVector::verboseLevel = 0
protected

Definition at line 207 of file G4PhysicsVector.hh.

Referenced by FillSecondDerivatives().


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