Geant4-11
Functions
G4gsmate.cc File Reference
#include <cmath>
#include "G4PhysicalConstants.hh"
#include "G4SystemOfUnits.hh"
#include "G3toG4.hh"
#include "G3MatTable.hh"
#include "G3EleTable.hh"
#include "G4Material.hh"
#include "G4Isotope.hh"
#include "G4UnitsTable.hh"

Go to the source code of this file.

Functions

void G4gsmate (G4int imate, G4String name, G4double ain, G4double zin, G4double densin, G4double, G4int, G4double *)
 
void PG4gsmate (G4String *tokens)
 

Function Documentation

◆ G4gsmate()

void G4gsmate ( G4int  imate,
G4String  name,
G4double  ain,
G4double  zin,
G4double  densin,
G4double  radl,
G4int  nwbf,
G4double ubuf 
)

Definition at line 104 of file G4gsmate.cc.

106{
107 G4double G3_minimum_density = 1.e-10*g/cm3;
108
109 // add units
110 G4double z = zin;
111 G4double a = ain*g/mole;
112 G4double dens = densin*g/cm3;
113
115
117 if (sname == "AIR") {
118 // handle the built in AIR mixture
119 G4double aa[2], zz[2], wmat[2];
120 aa[0] = 14.01*g/mole;
121 aa[1] = 16.00*g/mole;
122 zz[0] = 7;
123 zz[1] = 8;
124 wmat[0] = 0.7;
125 wmat[1] = 0.3;
126 // G4double theDensity = 1.2931*mg/cm3;
127 G4double theDensity = 0.0012931;
128 G4int n=2;
129 G4gsmixt(imate, sname, aa, zz, theDensity, n, wmat);
130 }
131 else if ( z<1 || dens < G3_minimum_density ) {
132 // define vacuum according to definition from N03 example
133 G4double density = universe_mean_density; //from PhysicalConstants.h
134 G4double pressure = 3.e-18*pascal;
135 G4double temperature = 2.73*kelvin;
136 material = new G4Material(name, z=1., a=1.01*g/mole, density,
137 kStateGas,temperature,pressure);
138 }
139 else {
140 //G4Element* element = CreateElement(z, a, name);
141 G4Element* element = G3Ele.GetEle(z);
142 material = new G4Material(name, dens, 1);
143 material->AddElement(element, 1.);
144 }
145
146 // add the material to the List
147 G3Mat.put(imate, material);
148}
G3G4DLL_API G3EleTable G3Ele
Definition: clparse.cc:59
void G4gsmixt(G4int imate, G4String name, G4double a[], G4double *z, G4double dens, G4int nlmat, G4double *wmat)
G3G4DLL_API G3MatTable G3Mat
Definition: clparse.cc:54
@ kStateGas
Definition: G4Material.hh:111
static constexpr double kelvin
Definition: G4SIunits.hh:274
static constexpr double mole
Definition: G4SIunits.hh:279
static constexpr double cm3
Definition: G4SIunits.hh:101
static constexpr double g
Definition: G4SIunits.hh:168
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define pascal
G4Element * GetEle(G4double Z)
Definition: G3EleTable.cc:51
void put(G4int id, G4Material *material)
Definition: G3MatTable.cc:52
const char * name(G4int ptype)
G4String strip_copy(G4String str, char c=' ')
Return copy of string with leading and trailing characters removed.
string material
Definition: eplot.py:19
int universe_mean_density
Definition: hepunit.py:306

References cm3, g, G3Ele, G3Mat, G4gsmixt(), G3EleTable::GetEle(), kelvin, kStateGas, eplot::material, mole, CLHEP::detail::n, G4InuclParticleNames::name(), pascal, G3MatTable::put(), G4StrUtil::strip_copy(), and source.hepunit::universe_mean_density.

Referenced by PG4gsmate().

◆ PG4gsmate()

void PG4gsmate ( G4String tokens)

Definition at line 41 of file G4gsmate.cc.

42{
43 // fill the parameter containers
44 G3fillParams(tokens,PTgsmate);
45 G4String name = Spar[0];
46 G4int imate = Ipar[0];
47 G4int nwbf = Ipar[1];
48 G4double a = Rpar[0];
49 G4double z = Rpar[1];
50 G4double dens = Rpar[2];
51 G4double radl = Rpar[3];
52 // G4double absl = Rpar[4];
53 G4double *ubuf = &Rpar[5];
54
55 G4gsmate(imate, name, a, z, dens, radl, nwbf, ubuf);
56}
G3G4DLL_API G4int Ipar[1000]
Definition: clparse.cc:65
void G3fillParams(G4String *tokens, const char *ptypes)
Definition: clparse.cc:218
G3G4DLL_API G4double Rpar[1000]
Definition: clparse.cc:66
G3G4DLL_API G4String Spar[1000]
Definition: clparse.cc:67
#define PTgsmate
Definition: G3toG4.hh:61
void G4gsmate(G4int imate, G4String name, G4double ain, G4double zin, G4double densin, G4double, G4int, G4double *)
Definition: G4gsmate.cc:104

References G3fillParams(), G4gsmate(), Ipar, G4InuclParticleNames::name(), PTgsmate, Rpar, and Spar.

Referenced by G3CLEval().