VGM Version 5.3
Loading...
Searching...
No Matches
Element.cxx
Go to the documentation of this file.
1// $Id$
2
3// -----------------------------------------------------------------------
4// The Geant4GM package of the Virtual Geometry Model
5// Copyright (C) 2007, Ivana Hrivnacova
6// All rights reserved.
7//
8// For the licensing terms see vgm/LICENSE.
9// Contact: ivana@ipno.in2p3.fr
10// -----------------------------------------------------------------------
11
12//
13// Class Element
14// ---------------
15// VGM implementation for Geant4 element.
16//
17// Author: Ivana Hrivnacova; IPN Orsay
18
19#include "ClhepVGM/Units.h"
20
25
26#include "G4Element.hh"
27
28//_____________________________________________________________________________
30 const std::string& name, const std::string& symbol, double z, double a)
31 : VGM::IElement(),
32 fElement(
33 new G4Element(name, symbol, z, a / ClhepVGM::Units::AtomicWeight()))
34{
41
42 // Register element in the map
43 ElementMap::Instance()->AddElement(this, fElement);
44}
45
46//_____________________________________________________________________________
47Geant4GM::Element::Element(const std::string& name, const std::string& symbol,
48 const VGM::IsotopeVector& isotopes,
49 const VGM::RelAbundanceVector& relAbundances)
50 : VGM::IElement(), fElement(0)
51{
58
59 if (!isotopes.size()) {
60 std::cerr << " Geant4GM::Element::Element: " << std::endl;
61 std::cerr << " No isotopes defined.";
62 std::cerr << "*** Error: Aborting execution ***" << std::endl;
63 exit(1);
64 }
65
66 // Check coherence
67 if (isotopes.size() != relAbundances.size()) {
68 std::cerr << " Geant4GM::Element::Element: " << std::endl;
69 std::cerr << " Isotopes size and relAbundances size differ.";
70 std::cerr << "*** Error: Aborting execution ***" << std::endl;
71 exit(1);
72 }
73
74 // Create element
75 fElement = new G4Element(name, symbol, isotopes.size());
76
77 // Add isotopes
78 for (unsigned int i = 0; i < isotopes.size(); i++) {
79 G4Isotope* g4Isotope = IsotopeMap::Instance()->GetIsotope(isotopes[i]);
80 fElement->AddIsotope(g4Isotope, relAbundances[i]);
81 }
82
83 // Register element in the map
84 ElementMap::Instance()->AddElement(this, fElement);
85}
86
87//_____________________________________________________________________________
88Geant4GM::Element::Element(G4Element* element)
89 : VGM::IElement(), fElement(element)
90{
92
93 // Register element in the map
94 ElementMap::Instance()->AddElement(this, fElement);
95}
96
97//_____________________________________________________________________________
99{
101}
102
103//_____________________________________________________________________________
104Geant4GM::Element::Element(const Element& rhs) : VGM::IElement(rhs)
105{
107}
108
109//_____________________________________________________________________________
114
115//
116// private functions
117//
118
119//_____________________________________________________________________________
120void Geant4GM::Element::CheckIndex(int i) const
121{
122 if (i < 0 || i >= NofIsotopes()) {
123 std::cerr << " Geant4GM::Element::CheckIndex: " << std::endl;
124 std::cerr << " Index of isotope outside limits." << std::endl;
125 std::cerr << "*** Error: Aborting execution ***" << std::endl;
126 exit(1);
127 }
128}
129
130//
131// public functions
132//
133
134//_____________________________________________________________________________
135std::string Geant4GM::Element::Name() const { return fElement->GetName(); }
136
137//_____________________________________________________________________________
138std::string Geant4GM::Element::Symbol() const { return fElement->GetSymbol(); }
139
140//_____________________________________________________________________________
141double Geant4GM::Element::Z() const { return fElement->GetZ(); }
142
143//_____________________________________________________________________________
144double Geant4GM::Element::N() const { return fElement->GetN(); }
145
146//_____________________________________________________________________________
148{
149 return fElement->GetA() * ClhepVGM::Units::AtomicWeight();
150}
151
152//_____________________________________________________________________________
154{
155 return fElement->GetNumberOfIsotopes();
156}
157
158//_____________________________________________________________________________
160{
161 CheckIndex(i);
162
163 const G4Isotope* g4Isotope = fElement->GetIsotope(i);
164 return IsotopeMap::Instance()->GetIsotope(const_cast<G4Isotope*>(g4Isotope));
165}
166
167//_____________________________________________________________________________
169{
170 CheckIndex(i);
171
172 return fElement->GetRelativeAbundanceVector()[i];
173}
static double AtomicWeight()
Return CLHEP default atomic weight unit in VGM units.
Definition Units.cxx:107
void AddElement(VGM::IElement *, G4Element *)
static ElementMap * Instance()
VGM implementation for Geant4 element.
Definition Element.h:32
virtual double N() const
Return the effective number of nucleons.
Definition Element.cxx:144
virtual int NofIsotopes() const
Return the number of isotopes constituing this element.
Definition Element.cxx:153
virtual double Z() const
Return the effective atomic number.
Definition Element.cxx:141
virtual VGM::IIsotope * Isotope(int i) const
Return the i-th isotope constituing this element.
Definition Element.cxx:159
virtual std::string Symbol() const
Return the symbol of this element.
Definition Element.cxx:138
virtual std::string Name() const
Return the name of this element.
Definition Element.cxx:135
virtual ~Element()
Definition Element.cxx:110
virtual double A() const
Return the effective effective mass of a mole in g/mole.
Definition Element.cxx:147
virtual double RelAbundance(int i) const
Return the relative abundance (the fraction of nb of atomes per volume) of the i-th isotope constitui...
Definition Element.cxx:168
static IsotopeMap * Instance()
G4Isotope * GetIsotope(VGM::IIsotope *iIsotope) const
The VGM interface to elements.
Definition IIsotope.h:28
ClhepVGM utilities.
Definition transform.h:29
VGM interfaces.
Definition VMedium.h:28
std::vector< double > RelAbundanceVector
Definition IElement.h:31
std::vector< IIsotope * > IsotopeVector
Definition IElement.h:30