File indexing completed on 2025-01-18 09:58:06
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040 #ifndef DensityEffectData_h
0041 #define DensityEffectData_h 1
0042
0043 #include "G4Material.hh"
0044 #include "globals.hh"
0045
0046 #include <vector>
0047
0048 const G4int NDENSDATA = 278;
0049 const G4int NDENSARRAY = 10;
0050 const G4int NDENSELEM = 98;
0051
0052 class G4DensityEffectData
0053 {
0054 public:
0055 G4DensityEffectData();
0056
0057 ~G4DensityEffectData() = default;
0058
0059
0060 G4DensityEffectData& operator=(const G4DensityEffectData& right) = delete;
0061 G4DensityEffectData(const G4DensityEffectData&) = delete;
0062
0063
0064 G4int GetElementIndex(G4int Z, G4State st = kStateUndefined) const;
0065
0066
0067 G4int GetIndex(const G4String& matName) const;
0068
0069
0070 void PrintData(const G4String& matName) const;
0071
0072
0073 void DumpData() const;
0074
0075
0076 inline G4double GetPlasmaEnergy(G4int idx) const;
0077 inline G4double GetAdjustmentFactor(G4int idx) const;
0078 inline G4double GetCdensity(G4int idx) const;
0079 inline G4double GetX0density(G4int idx) const;
0080 inline G4double GetX1density(G4int idx) const;
0081 inline G4double GetAdensity(G4int idx) const;
0082 inline G4double GetMdensity(G4int idx) const;
0083 inline G4double GetDelta0density(G4int idx) const;
0084 inline G4double GetErrorDensity(G4int idx) const;
0085 inline G4double GetMeanIonisationPotential(G4int idx) const;
0086
0087 private:
0088 void Initialize();
0089
0090 void AddMaterial(G4double* val, const G4String& matName);
0091
0092 G4double data[NDENSDATA][NDENSARRAY];
0093 std::vector<G4String> names;
0094
0095
0096 G4int indexZ[NDENSELEM];
0097 G4State state[NDENSELEM];
0098
0099 G4int index{0};
0100 };
0101
0102 inline G4double G4DensityEffectData::GetPlasmaEnergy(G4int idx) const
0103 {
0104 return (idx >= 0 && idx < NDENSDATA) ? data[idx][0] : DBL_MAX;
0105 }
0106
0107 inline G4double G4DensityEffectData::GetAdjustmentFactor(G4int idx) const
0108 {
0109 return (idx >= 0 && idx < NDENSDATA) ? data[idx][1] : DBL_MAX;
0110 }
0111
0112 inline G4double G4DensityEffectData::GetCdensity(G4int idx) const
0113 {
0114 return (idx >= 0 && idx < NDENSDATA) ? data[idx][2] : DBL_MAX;
0115 }
0116
0117 inline G4double G4DensityEffectData::GetX0density(G4int idx) const
0118 {
0119 return (idx >= 0 && idx < NDENSDATA) ? data[idx][3] : DBL_MAX;
0120 }
0121
0122 inline G4double G4DensityEffectData::GetX1density(G4int idx) const
0123 {
0124 return (idx >= 0 && idx < NDENSDATA) ? data[idx][4] : DBL_MAX;
0125 }
0126
0127 inline G4double G4DensityEffectData::GetAdensity(G4int idx) const
0128 {
0129 return (idx >= 0 && idx < NDENSDATA) ? data[idx][5] : DBL_MAX;
0130 }
0131
0132 inline G4double G4DensityEffectData::GetMdensity(G4int idx) const
0133 {
0134 return (idx >= 0 && idx < NDENSDATA) ? data[idx][6] : DBL_MAX;
0135 }
0136
0137 inline G4double G4DensityEffectData::GetDelta0density(G4int idx) const
0138 {
0139 return (idx >= 0 && idx < NDENSDATA) ? data[idx][7] : DBL_MAX;
0140 }
0141
0142 inline G4double G4DensityEffectData::GetErrorDensity(G4int idx) const
0143 {
0144 return (idx >= 0 && idx < NDENSDATA) ? data[idx][8] : DBL_MAX;
0145 }
0146
0147 inline G4double G4DensityEffectData::GetMeanIonisationPotential(G4int idx) const
0148 {
0149 return (idx >= 0 && idx < NDENSDATA) ? data[idx][9] : DBL_MAX;
0150 }
0151
0152 #endif