Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G3G4Interface.hh File Reference
#include "globals.hh"
Include dependency graph for G3G4Interface.hh:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

void G4gsvolu (G4String name, G4String shape, G4int nmed, G4double *par, G4int npar)
 
void G4gspos (G4String name, G4int num, G4String moth, G4double x, G4double y, G4double z, G4int irot, G4String only)
 
void G4gsposp (G4String name, G4int num, G4String moth, G4double x, G4double y, G4double z, G4int irot, G4String only, G4double Rpar[], G4int npar)
 
void G4gsbool (G4String volName, G4String manyVolName)
 
void G4gsrotm (G4int irot, G4double theta1, G4double phi1, G4double theta2, G4double phi2, G4double theta3, G4double phi3)
 
void G4gsatt (G4String name, G4String attr, G4int ival)
 
void G4gsdvn (G4String vname, G4String vmoth, G4int ndiv, G4int iaxis)
 
void G4gsdvt (G4String name, G4String moth, G4double Step, G4int iaxis, G4int numed, G4int ndvmx)
 
void G4gsdvx (G4String name, G4String moth, G4int ndiv, G4int iaxis, G4double Step, G4double c0, G4int numed, G4int ndvmx)
 
void G4gsdvn2 (G4String name, G4String moth, G4int ndiv, G4int iaxis, G4double c0, G4int numed)
 
void G4gsdvt2 (G4String name, G4String moth, G4double Step, G4int iaxis, G4double c0, G4int numed, G4int ndvmx)
 
void G4gsmate (G4int imate, G4String name, G4double a, G4double z, G4double dens, G4double radl, G4int nwbf, G4double *ubuf)
 
void G4gsmixt (G4int imate, G4String name, G4double a[], G4double *z, G4double dens, G4int nlmat, G4double *wmat)
 
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 npar)
 
void G4gstpar (G4int itmed, G4String chpar, G4double parval)
 
void G4gspart (G4int ipart, G4String chnpar, G4int itrtyp, G4double amass, G4double charge, G4double tlife, G4double *ubuf, G4int nwb)
 
void G4gsdk (G4int ipart, G4double *bratio, G4int *mode)
 
void G4gsdet (G4String chset, G4String chdet, G4int nv, G4String *chnmsv, G4int *nbitsv, G4int idtyp, G4int nwhi, G4int nwdi)
 
void G4gsdetv (G4String chset, G4String chdet, G4int idtyp, G4int nwhi, G4int nwdi)
 
void G4gsdeta (G4String chset, G4String chdet, G4String chali, G4int nwhi, G4int nwdi)
 
void G4gsdeth (G4String chset, G4String chdet, G4int nh, G4String *chnamh, G4int *nbitsh, G4double *orig, G4double *fact)
 
void G4gsdetd (G4String chset, G4String chdet, G4int nd, G4String *chnmsd, G4int *nbitsd)
 
void G4gsdetu (G4String chset, G4String chdet, G4int nupar, G4double *upar)
 
void G4ggclos ()
 
G4LogicalVolumeG4BuildGeom (G4String &inFile)
 

Function Documentation

G4LogicalVolume* G4BuildGeom ( G4String inFile)

Definition at line 55 of file G4BuildGeom.cc.

55  {
56 
57  G4int irot=0;
58  G4gsrotm(0, 90, 0, 90, 90, 0, 0);
59 
60  G4cout << "Instantiated unit rotation matrix irot=" << irot << G4endl;
61 
62  // Read the call List and interpret to Generate Geant4 geometry
63 
64  G4cout << "Reading the call List file " << inFile << "..." << G4endl;
65 
66  G3CLRead(inFile, 0);
67 
68  G3Part.PrintAll();
69 
70  G3Det.PrintAll();
71 
72  G3Vol.PrintAll();
73 
74  G4cout << "Call List file read completed. Build geometry" << G4endl;
75 
76  // Build the geometry
77 
78  G3VolTableEntry* topVTE = G3Vol.GetFirstVTE();
79  G4cout << "G3toG4 top level volume is " << topVTE->GetName() << G4endl;
80 
81  // modified
82  G3toG4BuildTree(topVTE, 0);
83 
84  // Retrieve the top-level G3toG4 logical mother volume pointer
85 
86  G4LogicalVolume* topLV = topVTE->GetLV();
87 
88  // position the top logical volume
89  // (in Geant3 the top volume is not positioned)
90  //
91  new G4PVPlacement(0, G4ThreeVector(), topLV->GetName(), topLV, 0, false, 0);
92 
93  // mark as invisible
94 
96 
97  G4cout << "Top-level G3toG4 logical volume " << topLV->GetName() << " "
98  << *(topLV->GetVisAttributes()) << G4endl;
99 
100  // check the geometry here
101 
102  #ifdef G3G4DEBUG
103  G4cout << "scan through G4LogicalVolumeStore:" << G4endl;
104  checkVol();
105  #endif
106 
107  return topLV;
108 }
G3G4DLL_API G3PartTable G3Part
Definition: clparse.cc:58
void PrintAll()
Definition: G3DetTable.cc:96
void G4gsrotm(G4int irot, G4double theta1, G4double phi1, G4double theta2, G4double phi2, G4double theta3, G4double phi3)
Definition: G4gsrotm.cc:54
CLHEP::Hep3Vector G4ThreeVector
G3VolTableEntry * GetFirstVTE()
Definition: G3VolTable.cc:107
G4LogicalVolume * GetLV()
void PrintAll()
Definition: G3PartTable.cc:79
G3G4DLL_API G3DetTable G3Det
Definition: clparse.cc:59
int G4int
Definition: G4Types.hh:78
void PrintAll()
Definition: G3VolTable.cc:61
void checkVol(G4LogicalVolume *, G4int)
Definition: G4BuildGeom.cc:118
G4GLOB_DLL std::ostream G4cout
const G4VisAttributes * GetVisAttributes() const
void G3CLRead(G4String &fname, char *select=0)
Definition: clparse.cc:99
void G3toG4BuildTree(G3VolTableEntry *curVTE, G3VolTableEntry *motherVTE)
G4String GetName()
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
static const G4VisAttributes & GetInvisible()
void SetVisAttributes(const G4VisAttributes *pVA)
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the call graph for this function:

void G4ggclos ( )

Definition at line 36 of file G4ggclos.cc.

36  {
37  G4cout << "G4ggclos: setting top-level VolTableEntry" << G4endl;
39 }
G4GLOB_DLL std::ostream G4cout
void SetFirstVTE()
Definition: G3VolTable.cc:97
#define G4endl
Definition: G4ios.hh:61
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsatt ( G4String  name,
G4String  attr,
G4int  ival 
)

Definition at line 46 of file G4gsatt.cc.

47 {
48  // get logical volume pointer
49  // G4LogicalVolume *lvol = G3Vol.GetVTE(name)->GetLV();
50  G4cerr << "G4gsatt not implemented" << G4endl;
51 }
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr

Here is the caller graph for this function:

void G4gsbool ( G4String  volName,
G4String  manyVolName 
)

Definition at line 35 of file G4gsbool.cc.

36 {
37  // find VTEs
38  G3VolTableEntry* vte = G3Vol.GetVTE(volName);
39  G3VolTableEntry* manyVTE = G3Vol.GetVTE(manyVolName);
40 
41  if (vte == 0) {
42  G4String text = "G4gsbool: '" + volName + "' has no VolTableEntry";
43  G4Exception("G4gsbool()", "G3toG40012", FatalException, text);
44  return;
45  }
46  else if (manyVTE == 0) {
47  // warning
48  G4cerr << "G4gsbool: '" << manyVolName << "' has no VolTableEntry."
49  << G4endl
50  << " Specified overlap will be ignored."
51  << G4endl;
52  return;
53  }
54  else {
55  manyVTE->AddOverlap(vte);
56  }
57 }
void AddOverlap(G3VolTableEntry *aOverlap)
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54
G4GLOB_DLL std::ostream G4cerr

Here is the call graph for this function:

void G4gsdet ( G4String  chset,
G4String  chdet,
G4int  nv,
G4String chnmsv,
G4int nbitsv,
G4int  idtyp,
G4int  nwhi,
G4int  nwdi 
)

Definition at line 51 of file G4gsdet.cc.

53 {
54  G4gsdetv(chset, chdet, idtyp, nwhi, nwdi);
55 }
void G4gsdetv(G4String chset, G4String chdet, G4int idtyp, G4int nwhi, G4int nwdi)
Definition: G4gsdetv.cc:51

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsdeta ( G4String  chset,
G4String  chdet,
G4String  chali,
G4int  nwhi,
G4int  nwdi 
)

Definition at line 53 of file G4gsdeta.cc.

55 {
56  G4int idtyp = G3Det.GetID(chset, chdet);
57  // just associate another sensitive detector structure with
58  // the volume chdet
59  G4gsdetv(chset, chdet, idtyp, nwhi, nwdi);
60 }
void G4gsdetv(G4String chset, G4String chdet, G4int idtyp, G4int nwhi, G4int nwdi)
Definition: G4gsdetv.cc:51
G3G4DLL_API G3DetTable G3Det
Definition: clparse.cc:59
G4int GetID(G4String &set, G4String &det)
Definition: G3DetTable.cc:70
int G4int
Definition: G4Types.hh:78

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsdetd ( G4String  chset,
G4String  chdet,
G4int  nd,
G4String chnmsd,
G4int nbitsd 
)

Definition at line 50 of file G4gsdetd.cc.

51 {
52  // Get pointer to detector chset
53  // G4VSensitiveDetector* sdet = G3Det.GetSD(chset, chdet);
54  // Add hits to sensitive detector
55  // for (G4int i=0; i<nd; i++) {
56  // $$$ sdet->AddDigi(chnmsd[i],nbitsd[i]);
57  // }
58 }

Here is the caller graph for this function:

void G4gsdeth ( G4String  chset,
G4String  chdet,
G4int  nh,
G4String chnamh,
G4int nbitsh,
G4double orig,
G4double fact 
)

Definition at line 52 of file G4gsdeth.cc.

54 {
55  // Get pointer to sensitive detector chset
56  // G4VSensitiveDetector* sdet = G3Det.GetSD(chset, chdet);
57  // Add hits to sensitive detector
58  // for (G4int i=0; i<nh; i++) {
59  // $$$ sdet->AddHit(chnamh[i],nbitsh[i],orig[i],fact[i]);
60  // }
61 }

Here is the caller graph for this function:

void G4gsdetu ( G4String  chset,
G4String  chdet,
G4int  nupar,
G4double upar 
)

Definition at line 45 of file G4gsdetu.cc.

46 {
47  // $$$ nothing right now
48 }

Here is the caller graph for this function:

void G4gsdetv ( G4String  chset,
G4String  chdet,
G4int  idtyp,
G4int  nwhi,
G4int  nwdi 
)

Definition at line 51 of file G4gsdetv.cc.

52 {
53  G4cout << "G4gsdetv not currently implemented." << G4endl;
54  /*
55  // get lvol for detector chdet
56  G4LogicalVolume *lvol = G3Vol.GetLV(chdet);
57  if (lvol == 0) {
58  G4cout << "G4gsdetv: Logical volume " << chdet << " not available. Skip." << G4endl;
59  return;
60  }
61  // Generate a sensitive detector structure
62  // G4VSensitiveDetector *sdet;
63  // $$$ G4VSensitiveDetector *sdet = new G4VSensitiveDetector(chset);
64  // inform the logical volume of its sensitive detector
65  // lvol->SetSensitiveDetector(sdet);
66  // $$$ sdet->SetID(idtyp);
67  // Add the sensitive detector to the table
68  // G3Det.put(chset,idtyp,sdet);
69  */
70 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the caller graph for this function:

void G4gsdk ( G4int  ipart,
G4double bratio,
G4int mode 
)

Definition at line 46 of file G4gsdk.cc.

47 {
48 /*
49  // create decay object for the particle
50  G4Decay *decay = new G4Decay();
51  // add decay modes
52  for (G4int i=0; i<6; i++) {
53  if (mode[i] != 0) {
54 // $$$ decay->AddMode(mode[i],bratio[i]);
55  }
56  }
57  // associate decay object with particle ipart
58  G4ParticleDefinition *part = G3Part.Get(ipart);
59 // $$$ part->SetDecay(decay);
60 */
61 }

Here is the caller graph for this function:

void G4gsdvn ( G4String  vname,
G4String  vmoth,
G4int  ndiv,
G4int  iaxis 
)

Definition at line 103 of file G4gsdvn.cc.

104 {
105  // find mother VTE
106  G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
107 
108  if (mvte == 0) {
109  G4String text = "G4gsdvn:'" + vmoth + "' has no VolTableEntry";
110  G4Exception("G4gsdvn()", "G3toG40013", FatalException, text);
111  return;
112  }
113  else {
114  // a new vte clone copy with division is created
115  // for each mother (clone copy)
116 
117  G4CreateCloneVTEWithDivision(vname, mvte, kDvn, ndiv, iaxis, 0, 0., 0.);
118  }
119 }
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int, G4double c0, G4double step)
Definition: G4gsdvn.cc:52
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsdvn2 ( G4String  name,
G4String  moth,
G4int  ndiv,
G4int  iaxis,
G4double  c0,
G4int  numed 
)

Definition at line 56 of file G4gsdvn2.cc.

58 {
59  // find mother VTE
60  G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
61  if (mvte == 0) {
62  G4String text = "G4gsdvn2:'" + vmoth + "' has no VolTableEntry";
63  G4Exception("G4gsdvn2()", "G3toG40025", FatalException, text);
64  return;
65  }
66  else {
67  // a new vte clone copy with division is created
68  // for each mother (clone copy)
69 
70  G4CreateCloneVTEWithDivision(vname, mvte,
71  kDvn2, ndiv, iaxis, numed, c0, 0.);
72  }
73 }
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int, G4double c0, G4double step)
Definition: G4gsdvn.cc:52
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the caller graph for this function:

void G4gsdvt ( G4String  name,
G4String  moth,
G4double  Step,
G4int  iaxis,
G4int  numed,
G4int  ndvmx 
)

Definition at line 57 of file G4gsdvt.cc.

59 {
60  // find mother VTE
61  G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
62  if (mvte == 0) {
63  G4String text = "G4gsdvt:'" + vmoth + "' has no VolTableEntry";
64  G4Exception("G4gsdvt()", "G3toG40014", FatalException, text);
65  return;
66  }
67  else {
68  // a new vte clone copy with division is created
69  // for each mother (clone copy)
70 
71  G4CreateCloneVTEWithDivision(vname, mvte,
72  kDvt, ndvmx, iaxis, numed, 0., step);
73  }
74 }
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int, G4double c0, G4double step)
Definition: G4gsdvn.cc:52
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsdvt2 ( G4String  name,
G4String  moth,
G4double  Step,
G4int  iaxis,
G4double  c0,
G4int  numed,
G4int  ndvmx 
)

Definition at line 58 of file G4gsdvt2.cc.

60 {
61  // find mother VTE
62  G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
63  if (mvte == 0) {
64  G4String text = "G4gsdvt2:'" + vmoth + "' has no VolTableEntry";
65  G4Exception("G4gsdvt2()", "G3toG40015", FatalException, text);
66  return;
67  }
68  else {
69  // a new vte clone copy with division is created
70  // for each mother (clone copy)
71 
72  G4CreateCloneVTEWithDivision(vname, mvte,
73  kDvt2, ndvmx, iaxis, numed, c0, step);
74  }
75 }
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int, G4double c0, G4double step)
Definition: G4gsdvn.cc:52
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the caller graph for this function:

void G4gsdvx ( G4String  name,
G4String  moth,
G4int  ndiv,
G4int  iaxis,
G4double  Step,
G4double  c0,
G4int  numed,
G4int  ndvmx 
)

Definition at line 58 of file G4gsdvx.cc.

60 {
61  // pass to gsdvn2 or gsdvt2
62  if (Step > 0.) {
63  G4gsdvt2(name,moth,Step,iaxis,c0,numed,ndvmx);
64  } else if (ndiv > 0) {
65  G4gsdvn2(name,moth,ndiv,iaxis,c0,numed);
66  }
67 }
void G4gsdvn2(G4String name, G4String moth, G4int ndiv, G4int iaxis, G4double c0, G4int numed)
Definition: G4gsdvn2.cc:56
void G4gsdvt2(G4String name, G4String moth, G4double Step, G4int iaxis, G4double c0, G4int numed, G4int ndvmx)
Definition: G4gsdvt2.cc:58

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsmate ( G4int  imate,
G4String  name,
G4double  a,
G4double  z,
G4double  dens,
G4double  radl,
G4int  nwbf,
G4double ubuf 
)

Definition at line 105 of file G4gsmate.cc.

107 {
108  G4double G3_minimum_density = 1.e-10*g/cm3;
109 
110  // add units
111  G4double z = zin;
112  G4double a = ain*g/mole;
113  G4double dens = densin*g/cm3;
114 
115  G4Material* material=0;
116 
117  G4String sname = name.strip(G4String::both);
118  if (sname == "AIR") {
119  // handle the built in AIR mixture
120  G4double aa[2], zz[2], wmat[2];
121  aa[0] = 14.01*g/mole;
122  aa[1] = 16.00*g/mole;
123  zz[0] = 7;
124  zz[1] = 8;
125  wmat[0] = 0.7;
126  wmat[1] = 0.3;
127  // G4double theDensity = 1.2931*mg/cm3;
128  G4double theDensity = 0.0012931;
129  G4int n=2;
130  G4gsmixt(imate, sname, aa, zz, theDensity, n, wmat);
131  }
132  else if ( z<1 || dens < G3_minimum_density ) {
133  // define vacuum according to definition from N03 example
134  G4double density = universe_mean_density; //from PhysicalConstants.h
135  G4double pressure = 3.e-18*pascal;
136  G4double temperature = 2.73*kelvin;
137  material = new G4Material(name, z=1., a=1.01*g/mole, density,
138  kStateGas,temperature,pressure);
139  }
140  else {
141  //G4Element* element = CreateElement(z, a, name);
142  G4Element* element = G3Ele.GetEle(z);
143  material = new G4Material(name, dens, 1);
144  material->AddElement(element, 1.);
145  }
146 
147  // add the material to the List
148  G3Mat.put(imate, material);
149 }
G4Element * GetEle(G4double Z)
Definition: G3EleTable.cc:52
G4String strip(G4int strip_Type=trailing, char c=' ')
void put(G4int id, G4Material *material)
Definition: G3MatTable.cc:53
G3G4DLL_API G3MatTable G3Mat
Definition: clparse.cc:55
static constexpr double g
Definition: G4SIunits.hh:183
int G4int
Definition: G4Types.hh:78
void G4gsmixt(G4int imate, G4String name, G4double a[], G4double *z, G4double dens, G4int nlmat, G4double *wmat)
G3G4DLL_API G3EleTable G3Ele
Definition: clparse.cc:60
#define pascal
static constexpr double kelvin
Definition: G4SIunits.hh:281
static constexpr double cm3
Definition: G4SIunits.hh:121
static constexpr double universe_mean_density
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:362
double G4double
Definition: G4Types.hh:76
static constexpr double mole
Definition: G4SIunits.hh:286

Here is the call graph for this function:

Here is the caller graph for this function:

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

Here is the caller graph for this function:

void G4gspart ( G4int  ipart,
G4String  chnpar,
G4int  itrtyp,
G4double  amass,
G4double  charge,
G4double  tlife,
G4double ubuf,
G4int  nwb 
)

Definition at line 51 of file G4gspart.cc.

53 {
54 }

Here is the caller graph for this function:

void G4gspos ( G4String  name,
G4int  num,
G4String  moth,
G4double  x,
G4double  y,
G4double  z,
G4int  irot,
G4String  only 
)

Definition at line 65 of file G4gspos.cc.

67 {
68  // find VTEs
69  G3VolTableEntry* vte = G3Vol.GetVTE(vname);
70  G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
71 
72  if (vte == 0) {
73  G4String text = "G4gspos: '" + vname + "' has no VolTableEntry";
74  G4Exception("G4gspos()", "G3toG40017", FatalException, text);
75  return;
76  }
77  else if (mvte == 0) {
78  G4String text = "G4gspos: '" + vmoth + "' has no VolTableEntry";
79  G4Exception("G4gspos()", "G3toG40018", FatalException, text);
80  return;
81  }
82  else {
83  if (!vte->HasNegPars()) {
84  // position vector
85  G4ThreeVector* offset = new G4ThreeVector(x*cm, y*cm, z*cm);
86 
87  // create a G3Pos object and add it to the vte
88  G3Pos* aG3Pos = new G3Pos(vmoth, num, offset, irot, vonly);
89  vte->AddG3Pos(aG3Pos);
90 
91  // loop over all mothers
92  for (G4int i=0; i<mvte->GetNoClones(); i++) {
93  // (mvte is retrieved from its "master" name
94  // -> there is no need to call GetMasterClone()
95  G3VolTableEntry* mvteClone = mvte->GetClone(i);
96  vte->AddMother(mvteClone);
97  mvteClone->AddDaughter(vte);
98  }
99  }
100  else {
101  // if vte has neg parameters
102  // a new vte clone copy is created for each mother (clone copy)
103  // and its parameters are derived from it if possible
104 
105  G4CreateCloneVTE(vte, mvte, vte->GetRpar(), vte->GetNpar(), num,
106  x, y, z, irot, vonly);
107  }
108  }
109 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G3Pos.hh:43
G4double * GetRpar()
G3VolTableEntry * GetClone(G4int i)
int G4int
Definition: G4Types.hh:78
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
G4int GetNoClones()
G4int GetNpar()
static constexpr double cm
Definition: G4SIunits.hh:119
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void AddG3Pos(G3Pos *aG3Pos)
void AddMother(G3VolTableEntry *aDaughter)
G4bool HasNegPars()
void G4CreateCloneVTE(G3VolTableEntry *vte, G3VolTableEntry *mvte, G4double pars[], G4int npar, G4int num, G4double x, G4double y, G4double z, G4int irot, G4String vonly)
Definition: G4gsposp.cc:193
unsigned offset
Definition: inflate.h:102
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54
void AddDaughter(G3VolTableEntry *aDaughter)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsposp ( G4String  name,
G4int  num,
G4String  moth,
G4double  x,
G4double  y,
G4double  z,
G4int  irot,
G4String  only,
G4double  Rpar[],
G4int  npar 
)

Definition at line 297 of file G4gsposp.cc.

300 {
301  // find VTEs
302  G3VolTableEntry* vte = G3Vol.GetVTE(vname);
303  G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
304 
305  if (vte == 0) {
306  G4String err_message1 = "G4gsposp: '" + vname + "' has no VolTableEntry";
307  G4Exception("G4psposp()", "G3toG40021", FatalException, err_message1);
308  return;
309  }
310  if (mvte == 0) {
311  G4String err_message2 = "G4gsposp: '" + vmoth + "' has no VolTableEntry";
312  G4Exception("G4psposp()", "G3toG40022", FatalException, err_message2);
313  return;
314  }
315  else {
316  // a new vte clone copy is created for each mother (clone copy)
317  // and its parameters are derived from it if possible
318 
319  G4CreateCloneVTE(vte, mvte, pars, npar, num, x, y, z, irot, vonly);
320  }
321 }
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void G4CreateCloneVTE(G3VolTableEntry *vte, G3VolTableEntry *mvte, G4double pars[], G4int npar, G4int num, G4double x, G4double y, G4double z, G4int irot, G4String vonly)
Definition: G4gsposp.cc:193
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gsrotm ( G4int  irot,
G4double  theta1,
G4double  phi1,
G4double  theta2,
G4double  phi2,
G4double  theta3,
G4double  phi3 
)

Definition at line 54 of file G4gsrotm.cc.

56 {
57  G4double degrad = pi/180;
58 
59  G4double th1r = theta1*degrad;
60  G4double th2r = theta2*degrad;
61  G4double th3r = theta3*degrad;
62 
63  G4double phi1r = phi1*degrad;
64  G4double phi2r = phi2*degrad;
65  G4double phi3r = phi3*degrad;
66 
67  // Construct unit vectors
68 
69  G4ThreeVector x(std::sin(th1r)*std::cos(phi1r), std::sin(th1r)*std::sin(phi1r), std::cos(th1r));
70  G4ThreeVector y(std::sin(th2r)*std::cos(phi2r), std::sin(th2r)*std::sin(phi2r), std::cos(th2r));
71  G4ThreeVector z(std::sin(th3r)*std::cos(phi3r), std::sin(th3r)*std::sin(phi3r), std::cos(th3r));
72 
73  // check for orthonormality and left-handedness
74 
75  G4double check = (x.cross(y))*z;
76  G4double tol = 1.0e-3;
77 
78  if (1-std::abs(check)>tol) {
79  G4cerr << "Coordinate axes forming rotation matrix "
80  << irot << " are not orthonormal.(" << 1-std::abs(check) << ")"
81  << G4endl;
82  G4cerr << " theta1=" << theta1;
83  G4cerr << " phi1=" << phi1;
84  G4cerr << " theta2=" << theta2;
85  G4cerr << " phi2=" << phi2;
86  G4cerr << " theta3=" << theta3;
87  G4cerr << " phi3=" << phi3;
88  G4cerr << G4endl;
89  G4Exception("G4gsrotm()", "G3toG40023", FatalException,
90  "Non orthogonal axes!");
91  return;
92  }
93  //else if (1+check<=tol) {
94  // G4cerr << "G4gsrotm warning: coordinate axes forming rotation "
95  // << "matrix " << irot << " are left-handed" << G4endl;
96  //}
97 
99 
100  rotp->SetRotationMatrixByRow(x, y, z);
101 
102  // add it to the List
103 
104  G3Rot.Put(irot, rotp);
105 }
G3G4DLL_API G3RotTable G3Rot
Definition: clparse.cc:57
unsigned long check
Definition: inflate.h:88
void SetRotationMatrixByRow(const G4ThreeVector &Row1, const G4ThreeVector &Row2, const G4ThreeVector &Row3)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
void Put(G4int id, G4RotationMatrix *matrix)
Definition: G3RotTable.cc:53
G4GLOB_DLL std::ostream G4cerr

Here is the call graph for this function:

Here is the caller graph for this function:

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  npar 
)

Definition at line 68 of file G4gstmed.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4gstpar ( G4int  itmed,
G4String  chpar,
G4double  parval 
)

Definition at line 45 of file G4gstpar.cc.

46 {
47  // set special tracking medium parameter. Apply to all logical
48  // volumes making use of the specified tracking medium.
49  G4cerr << "G4gstpar: not implemented." << G4endl;
50 }
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr

Here is the caller graph for this function:

void G4gsvolu ( G4String  name,
G4String  shape,
G4int  nmed,
G4double par,
G4int  npar 
)

Definition at line 73 of file G4gsvolu.cc.

75 {
76  /*
77  G4cout << "Creating logical volume " << vname << " shape " << shape
78  << " nmed " << nmed << " #pars "<< npar << " parameters (cm): ";
79  for (int ipar=0; ipar< npar; ipar++) G4cout << std::setw(8) << rpar[ipar];
80  G4cout << G4endl;
81  */
82  if (G3Vol.GetVTE(vname)) {
83  // abort if VTE with given name exists
84  G4String text = "G4gsvolu: Attempt to create volume " + vname + " twice.";
85  G4Exception("G4gsvolu()", "G3toG40024", FatalException, text);
86  return;
87  }
88  else {
89  G4CreateVTE(vname, shape, nmed, rpar, npar);
90  }
91 }
G3VolTableEntry * GetVTE(const G4String &Vname)
Definition: G3VolTable.cc:54
G3VolTableEntry * G4CreateVTE(G4String vname, G4String shape, G4int nmed, G4double Rpar[], G4int npar)
Definition: G4gsvolu.cc:51
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:54

Here is the call graph for this function:

Here is the caller graph for this function: