Geant4 11.2.2
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
 
G4DensityEffectDataoperator= (const G4DensityEffectData &right)=delete
 
 G4DensityEffectData (const G4DensityEffectData &)=delete
 
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
 

Detailed Description

Definition at line 52 of file G4DensityEffectData.hh.

Constructor & Destructor Documentation

◆ G4DensityEffectData() [1/2]

G4DensityEffectData::G4DensityEffectData ( )

Definition at line 42 of file G4DensityEffectData.cc.

42{ Initialize(); }

◆ ~G4DensityEffectData()

G4DensityEffectData::~G4DensityEffectData ( )
default

◆ G4DensityEffectData() [2/2]

G4DensityEffectData::G4DensityEffectData ( const G4DensityEffectData & )
delete

Member Function Documentation

◆ DumpData()

void G4DensityEffectData::DumpData ( ) const

Definition at line 1600 of file G4DensityEffectData.cc.

1601{
1602 G4cout << "======================================================================" << G4endl;
1603 G4cout << " Material Eplasma(eV) rho -C x0 x1 a m d0 err" << G4endl;
1604 G4cout << "======================================================================" << G4endl;
1605 for (G4int i = 0; i < NDENSDATA; ++i) {
1606 G4cout << std::setw(3) << i << ". " << std::setw(25) << names[i] << std::setw(8)
1607 << data[i][0] / eV;
1608 for (G4int j = 1; j < NDENSARRAY; ++j) {
1609 G4cout << std::setw(8) << data[i][j];
1610 }
1611 G4cout << G4endl;
1612 }
1613 G4cout << "======================================================================" << G4endl;
1614}
const G4int NDENSARRAY
const G4int NDENSDATA
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

Referenced by PrintData().

◆ GetAdensity()

G4double G4DensityEffectData::GetAdensity ( G4int idx) const
inline

Definition at line 127 of file G4DensityEffectData.hh.

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

◆ GetAdjustmentFactor()

G4double G4DensityEffectData::GetAdjustmentFactor ( G4int idx) const
inline

Definition at line 107 of file G4DensityEffectData.hh.

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

◆ GetCdensity()

G4double G4DensityEffectData::GetCdensity ( G4int idx) const
inline

Definition at line 112 of file G4DensityEffectData.hh.

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

◆ GetDelta0density()

G4double G4DensityEffectData::GetDelta0density ( G4int idx) const
inline

Definition at line 137 of file G4DensityEffectData.hh.

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

◆ GetElementIndex()

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

Definition at line 1552 of file G4DensityEffectData.cc.

1553{
1554 return (Z >= 0 && Z < NDENSELEM) ? indexZ[Z] : -1;
1555}
const G4int NDENSELEM

◆ GetErrorDensity()

G4double G4DensityEffectData::GetErrorDensity ( G4int idx) const
inline

Definition at line 142 of file G4DensityEffectData.hh.

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

◆ GetIndex()

G4int G4DensityEffectData::GetIndex ( const G4String & matName) const

Definition at line 1557 of file G4DensityEffectData.cc.

1558{
1559 G4int idx = -1;
1560
1561 for (G4int i = 0; i < NDENSDATA; ++i) {
1562 if (names[i] == matName) {
1563 idx = i;
1564 break;
1565 }
1566 }
1567 return idx;
1568}

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

◆ GetMdensity()

G4double G4DensityEffectData::GetMdensity ( G4int idx) const
inline

Definition at line 132 of file G4DensityEffectData.hh.

133{
134 return (idx >= 0 && idx < NDENSDATA) ? data[idx][6] : DBL_MAX;
135}

◆ GetMeanIonisationPotential()

G4double G4DensityEffectData::GetMeanIonisationPotential ( G4int idx) const
inline

Definition at line 147 of file G4DensityEffectData.hh.

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

Referenced by G4IonisParamMat::FindMeanExcitationEnergy().

◆ GetPlasmaEnergy()

G4double G4DensityEffectData::GetPlasmaEnergy ( G4int idx) const
inline

Definition at line 102 of file G4DensityEffectData.hh.

103{
104 return (idx >= 0 && idx < NDENSDATA) ? data[idx][0] : DBL_MAX;
105}

◆ GetX0density()

G4double G4DensityEffectData::GetX0density ( G4int idx) const
inline

Definition at line 117 of file G4DensityEffectData.hh.

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

◆ GetX1density()

G4double G4DensityEffectData::GetX1density ( G4int idx) const
inline

Definition at line 122 of file G4DensityEffectData.hh.

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

◆ operator=()

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

◆ PrintData()

void G4DensityEffectData::PrintData ( const G4String & matName) const

Definition at line 1581 of file G4DensityEffectData.cc.

1582{
1583 if (matName.empty() || "all" == matName) {
1584 DumpData();
1585 return;
1586 }
1587 G4int idx = GetIndex(matName);
1588 if (idx >= 0) {
1589 G4cout << "G4DensityEffectData for <" << matName << "> index= " << idx << G4endl;
1590 G4cout << "I(eV)= " << data[idx][9] / CLHEP::eV << "Eplasma(eV)= " << data[idx][0] / CLHEP::eV
1591 << " rho= " << data[idx][1] << " -C= " << data[idx][2] << " x0= " << data[idx][3]
1592 << " x1= " << data[idx][4] << " a= " << data[idx][5] << " m= " << data[idx][6]
1593 << " d0= " << data[idx][7] << " err= " << data[idx][8] << G4endl;
1594 }
1595 else {
1596 G4cout << "G4DensityEffectData does not have <" << matName << ">" << G4endl;
1597 }
1598}
G4int GetIndex(const G4String &matName) const

Referenced by G4NistMessenger::SetNewValue().


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