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

#include <G4DensityEffectData.hh>

Public Member Functions

 G4DensityEffectData ()
 
 ~G4DensityEffectData ()=default
 
G4int GetElementIndex (G4int Z, G4State st=kStateUndefined) const
 
G4int GetIndex (const G4String &matName) const
 
void PrintData (const G4String &matName) const
 
void DumpData () const
 
G4double GetPlasmaEnergy (G4int idx) const
 
G4double GetAdjustmentFactor (G4int idx) const
 
G4double GetCdensity (G4int idx) const
 
G4double GetX0density (G4int idx) const
 
G4double GetX1density (G4int idx) const
 
G4double GetAdensity (G4int idx) const
 
G4double GetMdensity (G4int idx) const
 
G4double GetDelta0density (G4int idx) const
 
G4double GetErrorDensity (G4int idx) const
 
G4double GetMeanIonisationPotential (G4int idx) const
 
G4DensityEffectDataoperator= (const G4DensityEffectData &right)=delete
 
 G4DensityEffectData (const G4DensityEffectData &)=delete
 

Detailed Description

Definition at line 55 of file G4DensityEffectData.hh.

Constructor & Destructor Documentation

◆ G4DensityEffectData() [1/2]

G4DensityEffectData::G4DensityEffectData ( )
explicit

Definition at line 46 of file G4DensityEffectData.cc.

46 :index(0)
47{
48 Initialize();
49}

◆ ~G4DensityEffectData()

G4DensityEffectData::~G4DensityEffectData ( )
default

◆ G4DensityEffectData() [2/2]

G4DensityEffectData::G4DensityEffectData ( const G4DensityEffectData )
delete

Member Function Documentation

◆ DumpData()

void G4DensityEffectData::DumpData ( ) const

Definition at line 1349 of file G4DensityEffectData.cc.

1350{
1351 G4cout << "======================================================================"
1352 << G4endl;
1353 G4cout << " Material Eplasma(eV) rho -C x0 x1 a m d0 err"
1354 << G4endl;
1355 G4cout << "======================================================================"
1356 << G4endl;
1357 for(G4int i=0; i<NDENSDATA; ++i) {
1358 G4cout << std::setw(3) << i << ". " << std::setw(25) << names[i]
1359 << std::setw(8) << data[i][0]/eV;
1360 for(G4int j=1; j<NDENSARRAY; ++j) { G4cout << std::setw(8) << data[i][j]; }
1361 G4cout << G4endl;
1362 }
1363 G4cout << "======================================================================"
1364 << G4endl;
1365}
const G4int NDENSARRAY
const G4int NDENSDATA
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

Referenced by PrintData().

◆ GetAdensity()

G4double G4DensityEffectData::GetAdensity ( G4int  idx) const
inline

Definition at line 132 of file G4DensityEffectData.hh.

133{
134 return (idx >= 0 && idx < NDENSDATA) ? data[idx][5] : DBL_MAX;
135}
#define DBL_MAX
Definition: templates.hh:62

◆ GetAdjustmentFactor()

G4double G4DensityEffectData::GetAdjustmentFactor ( G4int  idx) const
inline

Definition at line 112 of file G4DensityEffectData.hh.

113{
114 return (idx >= 0 && idx < NDENSDATA) ? data[idx][1] : DBL_MAX;
115}

◆ GetCdensity()

G4double G4DensityEffectData::GetCdensity ( G4int  idx) const
inline

Definition at line 117 of file G4DensityEffectData.hh.

118{
119 return (idx >= 0 && idx < NDENSDATA) ? data[idx][2] : DBL_MAX;
120}

◆ GetDelta0density()

G4double G4DensityEffectData::GetDelta0density ( G4int  idx) const
inline

Definition at line 142 of file G4DensityEffectData.hh.

143{
144 return (idx >= 0 && idx < NDENSDATA) ? data[idx][7] : DBL_MAX;
145}

◆ GetElementIndex()

G4int G4DensityEffectData::GetElementIndex ( G4int  Z,
G4State  st = kStateUndefined 
) const

Definition at line 1293 of file G4DensityEffectData.cc.

1294{
1295 return (Z >= 0 && Z < NDENSELEM) ? indexZ[Z] : -1;
1296}
const G4int NDENSELEM
const G4int Z[17]

◆ GetErrorDensity()

G4double G4DensityEffectData::GetErrorDensity ( G4int  idx) const
inline

Definition at line 147 of file G4DensityEffectData.hh.

148{
149 return (idx >= 0 && idx < NDENSDATA) ? data[idx][8] : DBL_MAX;
150}

◆ GetIndex()

G4int G4DensityEffectData::GetIndex ( const G4String matName) const

Definition at line 1298 of file G4DensityEffectData.cc.

1299{
1300 G4int idx = -1;
1301
1302 for (G4int i=0; i<NDENSDATA; ++i) {
1303 if ( names[i] == matName ) {
1304 idx = i;
1305 break;
1306 }
1307 }
1308 return idx;
1309}

Referenced by G4IonisParamMat::FindMeanExcitationEnergy(), and PrintData().

◆ GetMdensity()

G4double G4DensityEffectData::GetMdensity ( G4int  idx) const
inline

Definition at line 137 of file G4DensityEffectData.hh.

138{
139 return (idx >= 0 && idx < NDENSDATA) ? data[idx][6] : DBL_MAX;
140}

◆ GetMeanIonisationPotential()

G4double G4DensityEffectData::GetMeanIonisationPotential ( G4int  idx) const
inline

Definition at line 152 of file G4DensityEffectData.hh.

153{
154 return (idx >= 0 && idx < NDENSDATA) ? data[idx][9] : DBL_MAX;
155}

Referenced by G4IonisParamMat::FindMeanExcitationEnergy().

◆ GetPlasmaEnergy()

G4double G4DensityEffectData::GetPlasmaEnergy ( G4int  idx) const
inline

Definition at line 107 of file G4DensityEffectData.hh.

108{
109 return (idx >= 0 && idx < NDENSDATA) ? data[idx][0] : DBL_MAX;
110}

◆ GetX0density()

G4double G4DensityEffectData::GetX0density ( G4int  idx) const
inline

Definition at line 122 of file G4DensityEffectData.hh.

123{
124 return (idx >= 0 && idx < NDENSDATA) ? data[idx][3] : DBL_MAX;
125}

◆ GetX1density()

G4double G4DensityEffectData::GetX1density ( G4int  idx) const
inline

Definition at line 127 of file G4DensityEffectData.hh.

128{
129 return (idx >= 0 && idx < NDENSDATA) ? data[idx][4] : DBL_MAX;
130}

◆ operator=()

G4DensityEffectData & G4DensityEffectData::operator= ( const G4DensityEffectData right)
delete

◆ PrintData()

void G4DensityEffectData::PrintData ( const G4String matName) const

Definition at line 1322 of file G4DensityEffectData.cc.

1323{
1324 if(matName.empty() || "all" == matName)
1325 {
1326 DumpData();
1327 return;
1328 }
1329 G4int idx = GetIndex(matName);
1330 if(idx >= 0) {
1331 G4cout << "G4DensityEffectData for <" << matName
1332 << "> index= " << idx << G4endl;
1333 G4cout << "I(eV)= " << data[idx][9]/CLHEP::eV
1334 << "Eplasma(eV)= " << data[idx][0]/CLHEP::eV
1335 << " rho= " << data[idx][1]
1336 << " -C= " << data[idx][2]
1337 << " x0= " << data[idx][3]
1338 << " x1= " << data[idx][4]
1339 << " a= " << data[idx][5]
1340 << " m= " << data[idx][6]
1341 << " d0= " << data[idx][7]
1342 << " err= " << data[idx][8]
1343 << G4endl;
1344 } else {
1345 G4cout << "G4DensityEffectData does not have <" << matName << ">" << G4endl;
1346 }
1347}
G4int GetIndex(const G4String &matName) const

Referenced by G4NistMessenger::SetNewValue().


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