Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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

◆ G4gsmixt()

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

Definition at line 72 of file G4gsmixt.cc.

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

Referenced by PG4gsmixt().

◆ PG4gsmixt()

void PG4gsmixt ( G4String * tokens)

Definition at line 42 of file G4gsmixt.cc.

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

Referenced by G3CLEval().