Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions
G4gsmixt.cc File Reference
#include <iomanip>
#include "globals.hh"
#include "G4SystemOfUnits.hh"
#include "G3toG4.hh"
#include "G3EleTable.hh"
#include "G3MatTable.hh"
#include "G4Material.hh"
#include "G4Isotope.hh"

Go to the source code of this file.

Functions

void PG4gsmixt (G4String *tokens)
 
void G4gsmixt (G4int imate, G4String name, G4double *a, G4double *z, G4double dens, G4int nlmat, G4double *wmat)
 

Function Documentation

void G4gsmixt ( G4int  imate,
G4String  name,
G4double a,
G4double z,
G4double  dens,
G4int  nlmat,
G4double wmat 
)

Definition at line 73 of file G4gsmixt.cc.

References G4Material::AddElement(), python.hepunit::cm3, FatalException, g(), G3Ele, G3Mat, G4Exception(), G3EleTable::GetEle(), eplot::material, and G3MatTable::put().

Referenced by PG4gsmixt().

75 {
76  // in Geant3:
77  // After a call with ratios by number (negative number of elements),
78  // the ratio array is changed to the ratio by weight, so all successive
79  // calls with the same array must specify the number of elements as
80  // positive
81  G4int i=0;
82  if (nlmat<0) {
83  // in case of proportions given in atom counts (nlmat<0),
84  // the wmat[i] are converted to weight fractions
85  G4double aMol = 0.;
86  for (i=0; i<std::abs(nlmat); i++) {
87  // total molecular weight
88  aMol += wmat[i]*a[i];
89  }
90  if (aMol == 0.) {
91  G4String text = "G4mixt: Total molecular weight in " + name + " = 0.";
92  G4Exception("G4gsmixt()", "G3toG40016", FatalException, text);
93  return;
94  }
95  for (i=0; i<std::abs(nlmat); i++) {
96  // weight fractions
97  wmat[i] = wmat[i]*a[i]/aMol;
98  }
99  }
100 
101  // create material with given number of components
102  // (elements)
103 
105  = new G4Material(name, dens*g/cm3, std::abs(nlmat));
106  for (i=0; i<std::abs(nlmat); i++) {
107  // add units
108  // G4Element* element = G4Element(z[i], a[i]*g/mole, name);
109  G4Element* element = G3Ele.GetEle(z[i]);
110  material->AddElement(element, wmat[i]);
111  }
112 
113  // add the material to the List
114  G3Mat.put(imate, material);
115 }
G4Element * GetEle(G4double Z)
Definition: G3EleTable.cc:52
G4double z
Definition: TRTMaterials.hh:39
void put(G4int id, G4Material *material)
Definition: G3MatTable.cc:53
G3G4DLL_API G3MatTable G3Mat
Definition: clparse.cc:55
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G3G4DLL_API G3EleTable G3Ele
Definition: clparse.cc:60
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:345
double G4double
Definition: G4Types.hh:76
void PG4gsmixt ( G4String tokens)

Definition at line 43 of file G4gsmixt.cc.

References test::a, G4String::data(), G3fillParams(), G4gsmixt(), Ipar, PTgsmixt, Rpar, Spar, and z.

Referenced by G3CLEval().

44 {
45  // fill the parameter containers
46  G3fillParams(tokens,PTgsmixt);
47 
48  // interpret the parameters
49  G4String name = Spar[0].data();
50  G4int imate = Ipar[0];
51  G4int nlmat = Ipar[1];
52  //G4double dens = Rpar[0]*g/cm3;
53  G4double dens = Rpar[0];
54  G4double *a = Rpar + 1;
55  G4double *z = Rpar + 1+std::abs(nlmat);
56  G4double *wmat = Rpar + 1 + 2*std::abs(nlmat);
57 /*
58  for (int i=0; i<std::abs(nlmat); i++){
59  //Rpar[i]=Rpar[i]*g/mole;
60  Rpar[i]=Rpar[i];
61  };
62 */
63  G4gsmixt(imate,name,a,z,dens,nlmat,wmat);
64 }
G3G4DLL_API G4double Rpar[1000]
Definition: clparse.cc:67
G4double z
Definition: TRTMaterials.hh:39
const XML_Char * name
void G3fillParams(G4String *tokens, const char *ptypes)
Definition: clparse.cc:219
int G4int
Definition: G4Types.hh:78
void G4gsmixt(G4int imate, G4String name, G4double *a, G4double *z, G4double dens, G4int nlmat, G4double *wmat)
Definition: G4gsmixt.cc:73
G3G4DLL_API G4String Spar[1000]
Definition: clparse.cc:68
G3G4DLL_API G4int Ipar[1000]
Definition: clparse.cc:66
const char * data() const
double G4double
Definition: G4Types.hh:76
#define PTgsmixt
Definition: G3toG4.hh:64