Warning, file /include/Geant4/G4NistMaterialBuilder.hh was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
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 #ifndef G4NistMaterialBuilder_h
0028 #define G4NistMaterialBuilder_h 1
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056 #include "G4Material.hh"
0057 #include "globals.hh"
0058
0059 #include <CLHEP/Units/PhysicalConstants.h>
0060
0061 #include <vector>
0062
0063 class G4NistElementBuilder;
0064
0065 class G4NistMaterialBuilder
0066 {
0067 public:
0068 G4NistMaterialBuilder(G4NistElementBuilder*, G4int verb = 0);
0069
0070 ~G4NistMaterialBuilder() = default;
0071
0072
0073 inline G4Material* FindMaterial(const G4String& name) const;
0074 G4Material* FindOrBuildMaterial(const G4String& name, G4bool warning = true);
0075
0076
0077 inline G4Material* FindSimpleMaterial(G4int Z) const;
0078 G4Material* FindOrBuildSimpleMaterial(G4int Z, G4bool warning);
0079
0080
0081 G4Material* ConstructNewMaterial(const G4String& name, const std::vector<G4String>& elm,
0082 const std::vector<G4int>& nbAtoms, G4double dens, G4State state = kStateSolid,
0083 G4double temp = NTP_Temperature, G4double pressure = CLHEP::STP_Pressure);
0084
0085
0086 G4Material* ConstructNewMaterial(const G4String& name, const std::vector<G4String>& elm,
0087 const std::vector<G4double>& weight, G4double dens, G4State state = kStateSolid,
0088 G4double temp = NTP_Temperature, G4double pressure = CLHEP::STP_Pressure);
0089
0090
0091 G4Material* ConstructNewGasMaterial(
0092 const G4String& name, const G4String& nameDB, G4double temp, G4double pres);
0093
0094
0095 G4Material* ConstructNewIdealGasMaterial(const G4String& name, const std::vector<G4String>& elm,
0096 const std::vector<G4int>& nbAtoms, G4double temp = NTP_Temperature,
0097 G4double pressure = CLHEP::STP_Pressure);
0098
0099
0100 void SetVerbose(G4int val);
0101
0102
0103
0104
0105
0106
0107
0108 void ListMaterials(const G4String&) const;
0109
0110
0111 void ListNistSimpleMaterials() const;
0112 void ListNistCompoundMaterials() const;
0113 void ListHepMaterials() const;
0114 void ListSpaceMaterials() const;
0115 void ListBioChemicalMaterials() const;
0116
0117
0118 const std::vector<G4String>& GetMaterialNames() const;
0119
0120
0121 inline G4double GetMeanIonisationEnergy(G4int index) const;
0122 inline G4double GetNominalDensity(G4int index) const;
0123
0124 G4bool operator==(const G4NistMaterialBuilder&) const = delete;
0125 G4bool operator!=(const G4NistMaterialBuilder&) const = delete;
0126 G4NistMaterialBuilder(const G4NistMaterialBuilder&) = delete;
0127 const G4NistMaterialBuilder& operator=(const G4NistMaterialBuilder&) = delete;
0128
0129 private:
0130 void Initialise();
0131 void NistSimpleMaterials();
0132 void NistCompoundMaterials();
0133 void NistCompoundMaterials2();
0134 void HepAndNuclearMaterials();
0135 void SpaceMaterials();
0136 void BioChemicalMaterials();
0137
0138
0139
0140
0141 void AddMaterial(const G4String& nameMat, G4double dens, G4int Z = 0, G4double pot = 0.0,
0142 G4int ncomp = 1, G4State = kStateSolid, G4bool stp = true);
0143
0144 void AddGas(const G4String& nameMat, G4double T, G4double P);
0145
0146 void AddElementByWeightFraction(G4int Z, G4double);
0147 void AddElementByAtomCount(G4int Z, G4int);
0148
0149 void AddElementByWeightFraction(const G4String& name, G4double);
0150 void AddElementByAtomCount(const G4String& name, G4int);
0151
0152
0153 G4Material* BuildNistMaterial(const G4String& matname, G4bool warning);
0154 G4Material* BuildMaterial(G4int idx);
0155
0156 void DumpElm(G4int) const;
0157 void DumpMix(G4int) const;
0158
0159 private:
0160 G4NistElementBuilder* elmBuilder;
0161
0162 G4int verbose;
0163 G4int nMaterials{0};
0164 G4int nComponents{0};
0165 G4int nCurrent{0};
0166 G4int nElementary;
0167 G4int nNIST;
0168 G4int nHEP;
0169 G4int nSpace;
0170
0171 std::vector<G4String> names;
0172 std::vector<G4String> chFormulas;
0173
0174 std::vector<G4double> densities;
0175 std::vector<G4double> ionPotentials;
0176 std::vector<G4State> states;
0177 std::vector<G4double> fractions;
0178 std::vector<G4bool> atomCount;
0179 std::vector<G4int> components;
0180 std::vector<G4int> indexes;
0181 std::vector<G4int> elements;
0182 std::vector<G4int> matIndex;
0183 std::vector<G4bool> STP;
0184
0185 std::vector<G4int> idxGas;
0186 std::vector<G4double> gasTemperature;
0187 std::vector<G4double> gasPressure;
0188 };
0189
0190 inline const std::vector<G4String>& G4NistMaterialBuilder::GetMaterialNames() const
0191 {
0192 return names;
0193 }
0194
0195 inline G4double G4NistMaterialBuilder::GetMeanIonisationEnergy(G4int index) const
0196 {
0197 return (index >= 0 && index < nMaterials) ? ionPotentials[index] : 10.0 * index;
0198 }
0199
0200 inline G4double G4NistMaterialBuilder::GetNominalDensity(G4int index) const
0201 {
0202 return (index >= 0 && index < nMaterials) ? densities[index] : 0.0;
0203 }
0204
0205 inline G4Material* G4NistMaterialBuilder::FindMaterial(const G4String& name) const
0206 {
0207 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
0208 G4Material* ptr = nullptr;
0209 for (auto& mat : *theMaterialTable) {
0210 if (name == mat->GetName()) {
0211 ptr = mat;
0212 break;
0213 }
0214 }
0215 return ptr;
0216 }
0217
0218 inline G4Material* G4NistMaterialBuilder::FindSimpleMaterial(G4int Z) const
0219 {
0220 return (Z > 0 && Z < nElementary) ? FindMaterial(names[Z]) : nullptr;
0221 }
0222
0223 #endif