Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G3G4Interface.hh File Reference
#include "globals.hh"

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

◆ G4BuildGeom()

G4LogicalVolume * G4BuildGeom ( G4String & inFile)

Definition at line 54 of file G4BuildGeom.cc.

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

◆ G4ggclos()

void G4ggclos ( )

Definition at line 35 of file G4ggclos.cc.

35 {
36 G4cout << "G4ggclos: setting top-level VolTableEntry" << G4endl;
38}
void SetFirstVTE()
Definition G3VolTable.cc:96

Referenced by PG4ggclos().

◆ G4gsatt()

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

Definition at line 45 of file G4gsatt.cc.

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

Referenced by PG4gsatt().

◆ G4gsbool()

void G4gsbool ( G4String volName,
G4String manyVolName )

Definition at line 34 of file G4gsbool.cc.

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

◆ G4gsdet()

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

Definition at line 50 of file G4gsdet.cc.

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

Referenced by PG4gsdet().

◆ G4gsdeta()

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

Definition at line 52 of file G4gsdeta.cc.

54{
55 G4int idtyp = G3Det.GetID(chset, chdet);
56 // just associate another sensitive detector structure with
57 // the volume chdet
58 G4gsdetv(chset, chdet, idtyp, nwhi, nwdi);
59}
G4int GetID(G4String &set, G4String &det)
Definition G3DetTable.cc:69

Referenced by PG4gsdeta().

◆ G4gsdetd()

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

Definition at line 49 of file G4gsdetd.cc.

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

Referenced by PG4gsdetd().

◆ G4gsdeth()

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

Definition at line 51 of file G4gsdeth.cc.

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

Referenced by PG4gsdeth().

◆ G4gsdetu()

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

Definition at line 44 of file G4gsdetu.cc.

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

Referenced by PG4gsdetu().

◆ G4gsdetv()

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

Definition at line 50 of file G4gsdetv.cc.

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

Referenced by G4gsdet(), G4gsdeta(), and PG4gsdetv().

◆ G4gsdk()

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

Definition at line 45 of file G4gsdk.cc.

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

Referenced by PG4gsdk().

◆ G4gsdvn()

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

Definition at line 102 of file G4gsdvn.cc.

103{
104 // find mother VTE
105 G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
106
107 if (mvte == 0) {
108 G4String text = "G4gsdvn:'" + vmoth + "' has no VolTableEntry";
109 G4Exception("G4gsdvn()", "G3toG40013", FatalException, text);
110 return;
111 }
112 else {
113 // a new vte clone copy with division is created
114 // for each mother (clone copy)
115
116 G4CreateCloneVTEWithDivision(vname, mvte, kDvn, ndiv, iaxis, 0, 0., 0.);
117 }
118}
@ kDvn
Definition G3Division.hh:52
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int, G4double c0, G4double step)
Definition G4gsdvn.cc:51

Referenced by PG4gsdvn().

◆ G4gsdvn2()

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

Definition at line 55 of file G4gsdvn2.cc.

57{
58 // find mother VTE
59 G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
60 if (mvte == 0) {
61 G4String text = "G4gsdvn2:'" + vmoth + "' has no VolTableEntry";
62 G4Exception("G4gsdvn2()", "G3toG40025", FatalException, text);
63 return;
64 }
65 else {
66 // a new vte clone copy with division is created
67 // for each mother (clone copy)
68
70 kDvn2, ndiv, iaxis, numed, c0, 0.);
71 }
72}
@ kDvn2
Definition G3Division.hh:52
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int nmed, G4double c0, G4double step)
Definition G4gsdvn.cc:51

◆ G4gsdvt()

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

Definition at line 56 of file G4gsdvt.cc.

58{
59 // find mother VTE
60 G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
61 if (mvte == 0) {
62 G4String text = "G4gsdvt:'" + vmoth + "' has no VolTableEntry";
63 G4Exception("G4gsdvt()", "G3toG40014", FatalException, text);
64 return;
65 }
66 else {
67 // a new vte clone copy with division is created
68 // for each mother (clone copy)
69
71 kDvt, ndvmx, iaxis, numed, 0., step);
72 }
73}
@ kDvt
Definition G3Division.hh:52
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int nmed, G4double c0, G4double step)
Definition G4gsdvn.cc:51

Referenced by PG4gsdvt().

◆ G4gsdvt2()

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

Definition at line 57 of file G4gsdvt2.cc.

59{
60 // find mother VTE
61 G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
62 if (mvte == 0) {
63 G4String text = "G4gsdvt2:'" + vmoth + "' has no VolTableEntry";
64 G4Exception("G4gsdvt2()", "G3toG40015", FatalException, text);
65 return;
66 }
67 else {
68 // a new vte clone copy with division is created
69 // for each mother (clone copy)
70
72 kDvt2, ndvmx, iaxis, numed, c0, step);
73 }
74}
@ kDvt2
Definition G3Division.hh:52
void G4CreateCloneVTEWithDivision(G4String vname, G3VolTableEntry *mvte, G3DivType divType, G4int nofDivisions, G4int iaxis, G4int nmed, G4double c0, G4double step)
Definition G4gsdvn.cc:51

◆ G4gsdvx()

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

Definition at line 57 of file G4gsdvx.cc.

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

Referenced by PG4gsdvx().

◆ G4gsmate()

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

Definition at line 58 of file G4gsmate.cc.

60{
61 G4double G3_minimum_density = 1.e-10*g/cm3;
62
63 // add units
64 G4double z = zin;
65 G4double a = ain*g/mole;
66 G4double dens = densin*g/cm3;
67
68 G4Material* material=0;
69
70 G4String sname = G4StrUtil::strip_copy(name);
71 if (sname == "AIR") {
72 // handle the built in AIR mixture
73 G4double aa[2], zz[2], wmat[2];
74 aa[0] = 14.01*g/mole;
75 aa[1] = 16.00*g/mole;
76 zz[0] = 7;
77 zz[1] = 8;
78 wmat[0] = 0.7;
79 wmat[1] = 0.3;
80 // G4double theDensity = 1.2931*mg/cm3;
81 G4double theDensity = 0.0012931;
82 G4int n=2;
83 G4gsmixt(imate, sname, aa, zz, theDensity, n, wmat);
84 }
85 else if ( z<1 || dens < G3_minimum_density ) {
86 // define vacuum according to definition from N03 example
87 G4double density = universe_mean_density; //from PhysicalConstants.h
88 G4double pressure = 3.e-18*pascal;
89 G4double temperature = 2.73*kelvin;
90 material = new G4Material(name, z=1., a=1.01*g/mole, density,
91 kStateGas,temperature,pressure);
92 }
93 else {
94 //G4Element* element = CreateElement(z, a, name);
95 G4Element* element = G3Ele.GetEle(z);
96 material = new G4Material(name, dens, 1);
97 material->AddElement(element, 1.);
98 }
99
100 // add the material to the List
101 G3Mat.put(imate, material);
102}
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
double G4double
Definition G4Types.hh:83
#define pascal
G4Element * GetEle(G4double Z)
Definition G3EleTable.cc:51
void put(G4int id, G4Material *material)
Definition G3MatTable.cc:52
void AddElement(G4Element *elm, G4int nAtoms)

Referenced by PG4gsmate().

◆ G4gsmixt()

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

Referenced by G4gsmate().

◆ G4gspart()

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

Definition at line 50 of file G4gspart.cc.

52{
53}

Referenced by PG4gspart().

◆ G4gspos()

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

Definition at line 64 of file G4gspos.cc.

66{
67 // find VTEs
68 G3VolTableEntry* vte = G3Vol.GetVTE(vname);
69 G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
70
71 if (vte == 0) {
72 G4String text = "G4gspos: '" + vname + "' has no VolTableEntry";
73 G4Exception("G4gspos()", "G3toG40017", FatalException, text);
74 return;
75 }
76 else if (mvte == 0) {
77 G4String text = "G4gspos: '" + vmoth + "' has no VolTableEntry";
78 G4Exception("G4gspos()", "G3toG40018", FatalException, text);
79 return;
80 }
81 else {
82 if (!vte->HasNegPars()) {
83 // position vector
84 G4ThreeVector* offset = new G4ThreeVector(x*cm, y*cm, z*cm);
85
86 // create a G3Pos object and add it to the vte
87 G3Pos* aG3Pos = new G3Pos(vmoth, num, offset, irot, vonly);
88 vte->AddG3Pos(aG3Pos);
89
90 // loop over all mothers
91 for (G4int i=0; i<mvte->GetNoClones(); i++) {
92 // (mvte is retrieved from its "master" name
93 // -> there is no need to call GetMasterClone()
94 G3VolTableEntry* mvteClone = mvte->GetClone(i);
95 vte->AddMother(mvteClone);
96 mvteClone->AddDaughter(vte);
97 }
98 }
99 else {
100 // if vte has neg parameters
101 // a new vte clone copy is created for each mother (clone copy)
102 // and its parameters are derived from it if possible
103
104 G4CreateCloneVTE(vte, mvte, vte->GetRpar(), vte->GetNpar(), num,
105 x, y, z, irot, vonly);
106 }
107 }
108}
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:192
Definition G3Pos.hh:43
void AddMother(G3VolTableEntry *aDaughter)
G3VolTableEntry * GetClone(G4int i)
void AddG3Pos(G3Pos *aG3Pos)
G4double * GetRpar()
void AddDaughter(G3VolTableEntry *aDaughter)

Referenced by PG4gspos().

◆ G4gsposp()

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 296 of file G4gsposp.cc.

299{
300 // find VTEs
301 G3VolTableEntry* vte = G3Vol.GetVTE(vname);
302 G3VolTableEntry* mvte = G3Vol.GetVTE(vmoth);
303
304 if (vte == 0) {
305 G4String err_message1 = "G4gsposp: '" + vname + "' has no VolTableEntry";
306 G4Exception("G4psposp()", "G3toG40021", FatalException, err_message1);
307 return;
308 }
309 if (mvte == 0) {
310 G4String err_message2 = "G4gsposp: '" + vmoth + "' has no VolTableEntry";
311 G4Exception("G4psposp()", "G3toG40022", FatalException, err_message2);
312 return;
313 }
314 else {
315 // a new vte clone copy is created for each mother (clone copy)
316 // and its parameters are derived from it if possible
317
318 G4CreateCloneVTE(vte, mvte, pars, npar, num, x, y, z, irot, vonly);
319 }
320}
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:192

Referenced by PG4gsposp().

◆ G4gsrotm()

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

Definition at line 53 of file G4gsrotm.cc.

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

Referenced by G4BuildGeom(), and PG4gsrotm().

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

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 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().

◆ G4gstpar()

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

Definition at line 44 of file G4gstpar.cc.

45{
46 // set special tracking medium parameter. Apply to all logical
47 // volumes making use of the specified tracking medium.
48 G4cerr << "G4gstpar: not implemented." << G4endl;
49}

Referenced by PG4gstpar().

◆ G4gsvolu()

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

Definition at line 72 of file G4gsvolu.cc.

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

Referenced by PG4gsvolu().