34#include "G4Element.hh"
35#include "G4Material.hh"
36#include "G4NistManager.hh"
37#include "G4SystemOfUnits.hh"
41const double Geant4GM::MaterialFactory::fgkTolerance = 1e-09;
45 :
VGM::IMaterialFactory(),
46 BaseVGM::VMaterialFactory(
"Geant4_GM_Material_Factory"),
54 :
VGM::IMaterialFactory(rhs),
56 fVacuumElements(rhs.fVacuumElements)
81void Geant4GM::MaterialFactory::ImportIsotope(G4Isotope* isotope)
87 std::cout <<
"Importing isotope: ";
89 void* isotopePtr = (
void*)isotope;
90 std::cout << isotopePtr;
92 std::cout << std::endl;
94 std::cout << *isotope << std::endl;
102void Geant4GM::MaterialFactory::ImportElement(G4Element* element)
108 std::cout <<
"Importing element: ";
110 void* elementPtr = (
void*)element;
111 std::cout << elementPtr;
113 std::cout << std::endl;
115 std::cout << *element << std::endl;
123void Geant4GM::MaterialFactory::ImportMaterial(G4Material* material)
129 std::cout <<
"Importing material: ";
131 void* materialPtr = (
void*)material;
132 std::cout << materialPtr;
134 std::cout << std::endl;
136 std::cout << *material << std::endl;
144bool Geant4GM::MaterialFactory::CompareIsotopes(
const G4Element* g4Element,
152 if (isotopes.size() != g4Element->GetNumberOfIsotopes())
return false;
154 for (G4int i = 0; i < G4int(g4Element->GetNumberOfIsotopes()); ++i) {
156 const G4Isotope* g4Isotope = g4Element->GetIsotope(i);
157 double g4RelAbundance = g4Element->GetRelativeAbundanceVector()[i];
158 G4bool match =
false;
160 for (G4int j = 0; j < G4int(isotopes.size()); ++j) {
162 if (std::abs(vgmIsotope->
Z() - g4Isotope->GetZ()) < fgkTolerance &&
163 std::abs(vgmIsotope->
N() - g4Isotope->GetN()) < fgkTolerance &&
164 std::abs(vgmIsotope->
A() - g4Isotope->GetA() / (g / mole)) <
166 std::abs(relAbundances[j] - g4RelAbundance) < fgkTolerance) {
171 if (!match)
return false;
184 const std::string& name,
int z,
int n,
double a)
188 G4Isotope* g4Isotope = G4Isotope::GetIsotope(name,
false);
192 (std::abs(g4Isotope->GetZ() - z) >= fgkTolerance ||
193 std::abs(g4Isotope->GetN() - n) >= fgkTolerance ||
194 std::abs(g4Isotope->GetA() / (g / mole) - a) >= fgkTolerance)) {
204 IsotopeStore().push_back(vgmIsotope);
212 const std::string& name,
const std::string& symbol,
double z,
double a)
217 G4Element* g4Element = G4Element::GetElement(name,
false);
221 (std::abs(g4Element->GetZ() - z) >= fgkTolerance ||
222 std::abs(g4Element->GetA() / (g / mole) - a) >= fgkTolerance)) {
231 ElementStore().push_back(vgmElement);
240 fVacuumElements.insert(vgmElement);
246 ElementStore().push_back(vgmElement);
250 g4Element = G4Element::GetElement(name,
false);
253 for (
unsigned int i = 0; i < g4Element->GetNumberOfIsotopes(); ++i) {
254 G4Isotope* g4Isotope =
const_cast<G4Isotope*
>(g4Element->GetIsotope(i));
257 if (!vgmIsotope) ImportIsotope(g4Isotope);
271 G4Element* g4Element = G4Element::GetElement(name,
false);
274 if (g4Element && !CompareIsotopes(g4Element, isotopes, relAbundances)) {
283 ElementStore().push_back(vgmElement);
295 G4NistManager* nistManager = G4NistManager::Instance();
296 G4Element* g4Element = nistManager->FindOrBuildElement(z, isotopes);
298 std::cerr <<
"No element with z=" << z <<
" defined." << std::endl;
307 for (
unsigned i = 0; i < g4Element->GetNumberOfIsotopes(); i++)
308 ImportIsotope(
const_cast<G4Isotope*
>(g4Element->GetIsotope(i)));
311 ElementStore().push_back(vgmElement);
318 const std::string& name,
double density,
VGM::IElement* element,
323 bool isVacuum =
false;
324 if (fVacuumElements.find(element) != fVacuumElements.end()) isVacuum =
true;
329 MaterialStore().push_back(vgmMaterial);
335 const std::string& name,
double density,
VGM::IElement* element,
337 double temperature,
double pressure)
341 bool isVacuum =
false;
342 if (fVacuumElements.find(element) != fVacuumElements.end()) isVacuum =
true;
345 name, density, element, state, temperature, pressure, isVacuum);
347 MaterialStore().push_back(vgmMaterial);
361 MaterialStore().push_back(vgmMaterial);
369 double temperature,
double pressure)
374 name, density, elements, fractions, state, temperature, pressure);
376 MaterialStore().push_back(vgmMaterial);
390 MaterialStore().push_back(vgmMaterial);
398 double temperature,
double pressure)
403 name, density, elements, atomCounts, state, temperature, pressure);
405 MaterialStore().push_back(vgmMaterial);
411 int mediumId,
VGM::IMaterial* material,
int nofParameters,
double* parameters)
418 MediumStore().push_back(vgmMedium);
427 const G4IsotopeTable* isotopeTable = G4Isotope::GetIsotopeTable();
429 for (G4int i = 0; i < G4int(isotopeTable->size()); i++) {
430 G4Isotope* isotope = (*isotopeTable)[i];
431 ImportIsotope(isotope);
434 const G4ElementTable* elementTable = G4Element::GetElementTable();
436 for (G4int i = 0; i < G4int(elementTable->size()); i++) {
437 G4Element* element = (*elementTable)[i];
438 ImportElement(element);
441 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
443 for (G4int j = 0; j < G4int(materialTable->size()); j++) {
444 G4Material* material = (*materialTable)[j];
445 ImportMaterial(material);
static ElementMap * Instance()
G4Element * GetElement(VGM::IElement *iElement) const
VGM implementation for Geant4 element.
static IsotopeMap * Instance()
G4Isotope * GetIsotope(VGM::IIsotope *iIsotope) const
VGM implementation for Geant4 Isotope.
VGM material factory for Geant4.
virtual VGM::IIsotope * CreateIsotope(const std::string &name, int z, int n, double a)
Create a chemical isotope.
virtual ~MaterialFactory()
virtual bool Import()
Import native materials.
virtual VGM::IMaterial * CreateMaterial(const std::string &name, double density, VGM::IElement *element, double radlen, double intlen)
Create a material.
virtual VGM::IMedium * CreateMedium(const std::string &name, int mediumId, VGM::IMaterial *material, int nofParameters, double *parameters)
Create a tracking medium.
virtual VGM::IElement * CreateElement(const std::string &name, const std::string &symbol, double z, double a)
Create a chemical element.
static MaterialMap * Instance()
VGM implementation for Geant4 material.
The VGM implementation of interface to tracking medium.
The VGM interface to elements.
The VGM interface to elements.
virtual int N() const =0
Return the effective number of nucleons.
virtual double A() const =0
Return the effective effective mass of a mole in g/mole.
virtual int Z() const =0
Return the effective atomic number.
The VGM interface to materials.
The VGM interface to tracking medium.
void DebugInfo()
Debug printing.
std::vector< IMaterial * > MaterialStore
std::vector< double > RelAbundanceVector
std::vector< int > AtomCountVector
std::vector< double > MassFractionVector
std::vector< IElement * > ElementVector
std::vector< IIsotope * > IsotopeStore
std::vector< IIsotope * > IsotopeVector
std::vector< IElement * > ElementStore