Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4gstmed.cc File Reference
#include "G4SystemOfUnits.hh"
#include "G4LogicalVolume.hh"
#include "G3toG4.hh"
#include "G3MatTable.hh"
#include "G3MedTable.hh"
#include "G4UserLimits.hh"
#include "G4MagneticField.hh"
#include "G4Material.hh"

Go to the source code of this file.

Functions

void PG4gstmed (G4String *tokens)
 
void G4gstmed (G4int itmed, G4String, G4int nmat, G4int isvol, G4int, G4double, G4double, G4double stemax, G4double, G4double, G4double, G4double *, G4int useG3TMLimits)
 

Function Documentation

◆ G4gstmed()

void G4gstmed ( G4int  itmed,
G4String  name,
G4int  nmat,
G4int  isvol,
G4int  ifield,
G4double  fieldm,
G4double  tmaxfd,
G4double  stemax,
G4double  deemax,
G4double  epsil,
G4double  stmin,
G4double par,
G4int  useG3TMLimits 
)

Definition at line 67 of file G4gstmed.cc.

71{
72 // get the pointer to material nmat
73 G4Material* material = G3Mat.get(nmat);
74
75 // NB. there is the possibility for redundancy in the mag field
76 // and user limits objects. Who cares.
77 // Generate the mag field object
78 // $$$ G4MagneticField* field = new G4MagneticField(ifield, fieldm, tmaxfd);
79 G4MagneticField* field = 0;
80
81 // Generate the user limits object
82 // !!! only "stemax" has its equivalent in G4
83
84 G4UserLimits* limits = 0;
85 if (useG3TMLimits) {
86 limits = new G4UserLimits();
87 limits->SetMaxAllowedStep(stemax*cm);
88 // limits->SetG3DefaultCuts();
89 // this is arranged globally by physics manager
90 }
91
92 // Store this medium in the G3Med structure
93 G3Med.put(itmed, material, field, limits, isvol);
94}
G3G4DLL_API G3MatTable G3Mat
Definition: clparse.cc:54
G3G4DLL_API G3MedTable G3Med
Definition: clparse.cc:55
G4Material * get(G4int id) const
Definition: G3MatTable.cc:43
void put(G4int id, G4Material *material, G4MagneticField *field, G4UserLimits *limits, G4int isvol)
Definition: G3MedTable.cc:52
virtual void SetMaxAllowedStep(G4double ustepMax)

Referenced by PG4gstmed().

◆ PG4gstmed()

void PG4gstmed ( G4String tokens)

Definition at line 43 of file G4gstmed.cc.

44{
45 // fill the parameter containers
46 G3fillParams(tokens,PTgstmed);
47
48 // interpret the parameters
49 G4String name = Spar[0];
50 G4int itmed = Ipar[0];
51 G4int nmat = Ipar[1];
52 G4int isvol = Ipar[2];
53 G4int ifield = Ipar[3];
54 G4int nwbuf = Ipar[4];
55 G4double fieldm = Rpar[0];
56 G4double tmaxfd = Rpar[1];
57 G4double stemax = Rpar[2];
58 G4double deemax = Rpar[3];
59 G4double epsil = Rpar[4];
60 G4double stmin = Rpar[5];
61 G4double *ubuf = &Rpar[6];
62
63 G4gstmed(itmed,name,nmat,isvol,ifield,fieldm,tmaxfd,stemax,
64 deemax,epsil,stmin,ubuf,nwbuf);
65}
#define PTgstmed
Definition: G3toG4.hh:64
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
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void G4gstmed(G4int itmed, G4String, G4int nmat, G4int isvol, G4int, G4double, G4double, G4double stemax, G4double, G4double, G4double, G4double *, G4int useG3TMLimits)
Definition: G4gstmed.cc:67
const char * name(G4int ptype)

Referenced by G3CLEval().