FairShip
Loading...
Searching...
No Matches
geomGeant4 Namespace Reference

Functions

None _collectVolumeNames (volume, listOfVolumes)
 
None check4OrphanVolumes (fGeo)
 
None setMagnetField (flag=None)
 
def printWF (vl, alreadyPrinted, bool onlyWithField=True)
 
def nextLevel (lv, int magnetMass, onlyWithField, exclude, alreadyPrinted)
 
None printWeightsandFields (bool onlyWithField=True, exclude=None)
 
def addVMCFields (shipGeo, str controlFile="", bool verbose=False, bool withVirtualMC=True)
 
None printVMCFields ()
 
def getRunManager ()
 
None startUI ()
 
None debug ()
 

Function Documentation

◆ _collectVolumeNames()

None geomGeant4._collectVolumeNames (   volume,
  listOfVolumes 
)
protected
Recursively walk the geometry tree collecting volume names from nodes.

Definition at line 18 of file geomGeant4.py.

18def _collectVolumeNames(volume, listOfVolumes) -> None:
19 """Recursively walk the geometry tree collecting volume names from nodes."""
20 nodes = volume.GetNodes()
21 if not nodes:
22 return
23 for i in range(nodes.GetEntries()):
24 node = nodes.At(i)
25 name = node.GetVolume().GetName()
26 if name not in listOfVolumes:
27 listOfVolumes.append(name)
28 _collectVolumeNames(node.GetVolume(), listOfVolumes)
29
30

◆ addVMCFields()

def geomGeant4.addVMCFields (   shipGeo,
str   controlFile = "",
bool   verbose = False,
bool   withVirtualMC = True 
)
Define VMC B fields, e.g. global field, field maps, local or local+global fields

Definition at line 165 of file geomGeant4.py.

165def addVMCFields(shipGeo, controlFile: str = "", verbose: bool = False, withVirtualMC: bool = True):
166 """
167 Define VMC B fields, e.g. global field, field maps, local or local+global fields
168 """
169 print("Calling addVMCFields")
170
171 fieldMaker = ROOT.ShipFieldMaker(verbose)
172
173 # Read the input control file. If this is empty then the only fields that are
174 # defined (so far) are those within the C++ geometry classes
175 if controlFile != "":
176 fieldMaker.readInputFile(controlFile)
177
178 # Set the main spectrometer field map as a global field
179 if hasattr(shipGeo, "Bfield"):
180 fieldsList = []
181 fieldMaker.defineFieldMap("MainSpecMap", shipGeo.Bfield.fieldMap, ROOT.TVector3(0.0, 0.0, shipGeo.Bfield.z))
182 fieldsList.append("MainSpecMap")
183
184 if not shipGeo.muShield.WithConstField:
185 offset = shipGeo.muShield.Entrance[0]
186 quadSymm = True
187 file_name = f"files/{shipGeo.shieldName}.root"
188 fieldMaker.defineFieldMap(
189 "muonShieldField", file_name, ROOT.TVector3(0.0, 0.0, offset), ROOT.TVector3(0.0, 0.0, 0.0), quadSymm
190 )
191 fieldsList.append("muonShieldField")
192
193 elif not shipGeo.hadronAbsorber.WithConstField:
194 fieldMaker.defineFieldMap(
195 "HadronAbsorberMap",
196 "files/FieldHadronStopper_raised_20190411.root",
197 ROOT.TVector3(0.0, 0.0, shipGeo.hadronAbsorber.z),
198 )
199 fieldsList.append("HadronAbsorberMap")
200
201 # Combine the fields to obtain the global field
202 if len(fieldsList) > 1:
203 fieldMaker.defineComposite("TotalField", *fieldsList) # fieldsList MUST have length <=4
204 fieldMaker.defineGlobalField("TotalField")
205 else:
206 fieldMaker.defineGlobalField("MainSpecMap")
207 if withVirtualMC:
208 # Force the VMC to update/reset the fields defined by the fieldMaker object.
209 # Get the ROOT/Geant4 geometry manager
210 geom = ROOT.TG4GeometryManager.Instance()
211 # Let the geometry know about the fieldMaker object
212 geom.SetUserPostDetConstruction(fieldMaker)
213 # Update the fields via the overridden ShipFieldMaker::Construct() function
214 geom.ConstructSDandField()
215
216 # Return the fieldMaker object, otherwise it will "go out of scope" and its
217 # content will be deleted
218 return fieldMaker
219
220

◆ check4OrphanVolumes()

None geomGeant4.check4OrphanVolumes (   fGeo)

Definition at line 31 of file geomGeant4.py.

31def check4OrphanVolumes(fGeo) -> None:
32 # fill list with volumes from nodes and compare with list of volumes
33 top = fGeo.GetTopVolume()
34 listOfVolumes = [top.GetName()]
35 _collectVolumeNames(top, listOfVolumes)
36 orphan = []
37 gIndex: dict[str, list] = {}
38 for v in fGeo.GetListOfVolumes():
39 name = v.GetName()
40 if name not in listOfVolumes:
41 orphan.append(name)
42 if name not in gIndex:
43 gIndex[name] = []
44 gIndex[name].append(v.GetNumber())
45 print("list of orphan volumes:", orphan)
46 vSame = {}
47 for x in gIndex:
48 if len(gIndex[x]) > 1:
49 vSame[x] = len(gIndex[x])
50 print("list of volumes with same name", vSame)
51
52

◆ debug()

None geomGeant4.debug ( )

Definition at line 255 of file geomGeant4.py.

255def debug() -> None:
256 gt = ROOT.G4TransportationManager.GetTransportationManager()
257 gn = gt.GetNavigatorForTracking()
258 world = gn.GetWorldVolume().GetLogicalVolume()
259 vmap = {}
260 for da in range(world.GetNoDaughters()):
261 vl = world.GetDaughter(da)
262 vmap[vl.GetName().c_str()] = vl
263 print(da, vl.GetName())
264 lvl = vmap["MagB"].GetLogicalVolume()
265 print(lvl.GetMass() / G4Unit.kg, lvl.GetMaterial().GetName())
266 print(lvl.GetFieldManager())
267 #
268 for da in range(world.GetNoDaughters()):
269 vl = world.GetDaughter(da)
270 vln = vl.GetName().c_str()
271 lvl = vl.GetLogicalVolume()
272 fm = lvl.GetFieldManager()
273 if fm:
274 v = fm.GetDetectorField().GetConstantFieldValue()
275 print(vln, fm, v.getX(), v.getY())
276 # FairROOT view
277 fgeom = ROOT.gGeoManager
278 magB = fgeom.GetVolume("MagB")
279 fl = magB.GetField()
280 print(fl.GetFieldValue()[0], fl.GetFieldValue()[1], fl.GetFieldValue()[2])
const Int_t debug

◆ getRunManager()

def geomGeant4.getRunManager ( )

Definition at line 246 of file geomGeant4.py.

246def getRunManager():
247 return ROOT.G4RunManager.GetRunManager()
248
249

◆ nextLevel()

def geomGeant4.nextLevel (   lv,
int  magnetMass,
  onlyWithField,
  exclude,
  alreadyPrinted 
)

Definition at line 134 of file geomGeant4.py.

134def nextLevel(lv, magnetMass: int, onlyWithField, exclude, alreadyPrinted):
135 tmp = 0
136 for da in range(lv.GetNoDaughters()):
137 lvn = lv.GetDaughter(da)
138 name = lvn.GetName().c_str()
139 if name in exclude:
140 continue
141 lvln = lvn.GetLogicalVolume()
142 if lvln.GetNoDaughters() > 0:
143 xtmp, _dummy = nextLevel(lvln, magnetMass, onlyWithField, exclude, alreadyPrinted)
144 magnetMass += xtmp
145 else:
146 tmp += printWF(lvn, alreadyPrinted, onlyWithField)
147 return tmp, magnetMass
148
149

◆ printVMCFields()

None geomGeant4.printVMCFields ( )
Method to print out information about VMC fields

Definition at line 221 of file geomGeant4.py.

221def printVMCFields() -> None:
222 """
223 Method to print out information about VMC fields
224 """
225 print("Printing VMC fields and associated volumes")
226
227 fGeo = ROOT.gGeoManager
228 vols = fGeo.GetListOfVolumes()
229
230 for v in vols:
231 field = v.GetField()
232 if field:
233 print(f"Vol is {v.GetName()}, field is {field}")
234 else:
235 print(f"Vol is {v.GetName()}")
236
237 if field:
238 # Get the field value assuming the global coordinate origin.
239 # This needs to be modified to use the local volume centre
240 centre = array("d", [0.0, 0.0, 0.0])
241 B = array("d", [0.0, 0.0, 0.0])
242 field.Field(centre, B)
243 print(f"Volume {v.GetName()} has B = ({B[0] / u.tesla}, {B[1] / u.tesla}, {B[2] / u.tesla}) T")
244
245

◆ printWeightsandFields()

None geomGeant4.printWeightsandFields ( bool   onlyWithField = True,
  exclude = None 
)

Definition at line 150 of file geomGeant4.py.

150def printWeightsandFields(onlyWithField: bool = True, exclude=None) -> None:
151 if exclude is None:
152 exclude = []
153 if len(exclude) != 0:
154 print("will not search in ", exclude)
155 gt = ROOT.G4TransportationManager.GetTransportationManager()
156 gn = gt.GetNavigatorForTracking()
157 world = gn.GetWorldVolume().GetLogicalVolume()
158 magnetMass = 0
159 alreadyPrinted: dict[str, str] = {}
160 _dummy, nM = nextLevel(world, magnetMass, onlyWithField, exclude, alreadyPrinted)
161 print("total magnet mass", nM / 1000.0, "t")
162 return
163
164

◆ printWF()

def geomGeant4.printWF (   vl,
  alreadyPrinted,
bool   onlyWithField = True 
)

Definition at line 104 of file geomGeant4.py.

104def printWF(vl, alreadyPrinted, onlyWithField: bool = True):
105 magnetMass = 0
106 vname = vl.GetName().data()
107 if vname in alreadyPrinted:
108 return magnetMass
109 vln = vname + " " + str(vl.GetCopyNo())
110 mvl = vl.GetMotherLogical().GetName().data()
111 alreadyPrinted[vname] = mvl
112 if mvl != "cave":
113 vln = mvl + "/" + vln
114 lvl = vl.GetLogicalVolume()
115 cvol = lvl.GetSolid().GetCubicVolume() / G4Unit.m3
116 M = lvl.GetMass() / G4Unit.kg
117 fm = lvl.GetFieldManager()
118 if not fm and onlyWithField:
119 return magnetMass
120 if M < 5000.0:
121 print("%-35s volume = %5.2Fm3 mass = %5.2F kg" % (vln, cvol, M))
122 else:
123 print("%-35s volume = %5.2Fm3 mass = %5.2F t" % (vln, cvol, M / 1000.0))
124 if fm:
125 fi = fm.GetDetectorField()
126 print(" Magnetic field:", fi.GetConstantFieldValue() / G4Unit.tesla)
127 # if vl.GetName().c_str()[0:3]=='Mag': magnetMass = M # only count volumes starting with Mag
128 name = vl.GetName().c_str()
129 if "_" in name and "Mag" in name.split("_")[1]:
130 magnetMass = M # only count volumes starting with Mag
131 return magnetMass
132
133

◆ setMagnetField()

None geomGeant4.setMagnetField (   flag = None)

Definition at line 53 of file geomGeant4.py.

53def setMagnetField(flag=None) -> None:
54 print("setMagnetField() called. Out of date, does not set field for tau neutrino detector!")
55 fGeo = ROOT.gGeoManager
56 vols = fGeo.GetListOfVolumes()
57 # copy field by hand to geant4
58 listOfFields = {}
59 for v in vols:
60 field = v.GetField()
61 if not field:
62 continue
63 bx = field.GetFieldValue()[0] / u.tesla * G4Unit.tesla
64 by = field.GetFieldValue()[1] / u.tesla * G4Unit.tesla
65 bz = field.GetFieldValue()[2] / u.tesla * G4Unit.tesla
66 magFieldIron = ROOT.G4UniformMagField(ROOT.G4ThreeVector(bx, by, bz))
67 FieldIronMgr = ROOT.G4FieldManager(magFieldIron)
68 FieldIronMgr.CreateChordFinder(magFieldIron)
69 listOfFields[v.GetName()] = FieldIronMgr
70 gt = ROOT.G4TransportationManager.GetTransportationManager()
71 gn = gt.GetNavigatorForTracking()
72 world = gn.GetWorldVolume().GetLogicalVolume()
73 setField = {}
74 for da in range(world.GetNoDaughters()):
75 vl0 = world.GetDaughter(da)
76 vln = vl0.GetName().c_str()
77 lvl0 = vl0.GetLogicalVolume()
78 if vln in listOfFields:
79 setField[lvl0] = vln
80 for dda in range(lvl0.GetNoDaughters()):
81 vl = lvl0.GetDaughter(dda)
82 vln = vl.GetName().c_str()
83 lvl = vl.GetLogicalVolume()
84 if vln in listOfFields:
85 setField[lvl] = vln
86 modified = False
87 for lvl in setField:
88 # check if field already exists
89 fm = lvl.GetFieldManager()
90 if not fm.DoesFieldExist():
91 lvl.SetFieldManager(listOfFields[setField[lvl]], True)
92 modified = True
93 if flag == "dump":
94 constField = listOfFields[setField[lvl]].GetDetectorField().GetConstantFieldValue()
95 print("set field for ", setField[lvl], constField)
96 else:
97 if flag == "dump":
98 print("field already set:", setField[lvl])
99 if modified:
100 g4Run = ROOT.G4RunManager.GetRunManager()
101 g4Run.GeometryHasBeenModified(True)
102
103

◆ startUI()

None geomGeant4.startUI ( )

Definition at line 250 of file geomGeant4.py.

250def startUI() -> None:
251 session = ROOT.G4UIterminal()
252 session.SessionStart()
253
254