SND@LHC Software
Loading...
Searching...
No Matches
run_simSND.py
Go to the documentation of this file.
1#!/usr/bin/env python
2import os
3import sys
4import ROOT
5import numpy as np
6
7import shipunit as u
8import shipRoot_conf
9import rootUtils as ut
10from ShipGeoConfig import ConfigRegistry
11from argparse import ArgumentParser
12
13mcEngine = "TGeant4"
14simEngine = "Pythia8" # "Genie" # Ntuple
15inactivateMuonProcesses = False
16
17MCTracksWithHitsOnly = False # copy particles which produced a hit and their history
18MCTracksWithEnergyCutOnly = True # copy particles above a certain kin energy cut
19MCTracksWithHitsOrEnergyCut = False # or of above, factor 2 file size increase compared to MCTracksWithEnergyCutOnly
20
21parser = ArgumentParser()
22
23parser.add_argument("--H6", dest="testbeam", help="use geometry of H8/H6 testbeam setup", action="store_true")
24parser.add_argument("--HX", dest="testbeam2023", help="use geometry of 2023 testbeam setup", action="store_true")
25parser.add_argument("--Genie", dest="genie", help="Genie for reading and processing neutrino interactions (1 standard, 2 FLUKA, 3 Pythia, 4 GENIE geometry driver)", required=False, default = 0, type = int)
26parser.add_argument("--Ntuple", dest="ntuple", help="Use ntuple as input", required=False, action="store_true")
27parser.add_argument("--MuonBack",dest="muonback", help="Generate events from muon background file, --Cosmics=0 for cosmic generator data", required=False, action="store_true")
28parser.add_argument("--Pythia8", dest="pythia8", help="Use Pythia8", required=False, action="store_true")
29parser.add_argument("--PG", dest="pg", help="Use Particle Gun", required=False, action="store_true")
30parser.add_argument("--pID", dest="pID", help="id of particle used by the gun (default=22)", required=False, default=22, type=int)
31parser.add_argument("--Estart", dest="Estart", help="start of energy range of particle gun for muflux detector (default=10 GeV)", required=False, default=10, type=float)
32parser.add_argument("--Eend", dest="Eend", help="end of energy range of particle gun for muflux detector (default=10 GeV)", required=False, default=10, type=float)
33
34parser.add_argument("--multiplePGSources", help="Multiple particle guns in a x-y plane at a fixed z or in a 3D volume", action="store_true")
35parser.add_argument("--EVx", dest="EVx", help="particle gun start xpos", required=False, default=0, type=float)
36parser.add_argument("--EVy", dest="EVy", help="particle gun start ypos", required=False, default=0, type=float)
37parser.add_argument("--EVz", dest="EVz", help="particle gun start zpos", required=False, default=0, type=float)
38parser.add_argument("--Dx", help="size of the full uniform spread of PG xpos", type=float)
39parser.add_argument("--Dy", help="size of the full uniform spread of PG ypos", type=float)
40parser.add_argument("--nZSlices", help="number of z slices for the PG sources", type=int)
41parser.add_argument("--zSliceStep", help="distance between the z slices for the PG sources", type=float)
42
43parser.add_argument("--FollowMuon",dest="followMuon", help="Make muonshield active to follow muons", required=False, action="store_true")
44parser.add_argument("--FastMuon", dest="fastMuon", help="Only transport muons for a fast muon only background estimate", required=False, action="store_true")
45parser.add_argument('--eMin', type=float, help="energy cut", dest='ecut', default=-1.)
46parser.add_argument('--zMax', type=float, help="max distance to apply energy cut", dest='zmax', default=70000.)
47parser.add_argument("--Nuage", dest="nuage", help="Use Nuage, neutrino generator of OPERA", required=False, action="store_true")
48parser.add_argument("--MuDIS", dest="mudis", help="Use muon deep inelastic scattering generator", required=False, action="store_true")
49parser.add_argument("-n", "--nEvents",dest="nEvents", help="Number of events to generate", required=False, default=100, type=int)
50parser.add_argument("-i", "--firstEvent",dest="firstEvent", help="First event of input file to use", required=False, default=0, type=int)
51parser.add_argument("-s", "--seed",dest="theSeed", help="Seed for random number. Only for experts, see TRrandom::SetSeed documentation", required=False, default=0, type=int)
52parser.add_argument("-f", dest="inputFile", help="Input file if not default file", required=False, default=False)
53parser.add_argument("-g", dest="geofile", help="geofile for muon shield geometry, for experts only", required=False, default=None)
54parser.add_argument("-o", "--output",dest="outputDir", help="Output directory", required=False, default=".")
55parser.add_argument("--boostFactor", dest="boostFactor", help="boost mu brems", required=False, type=float,default=0)
56parser.add_argument("--enhancePiKaDecay", dest="enhancePiKaDecay", help="decrease charged pion and kaon lifetime", required=False, type=float,default=0.)
57parser.add_argument("--debug", dest="debug", help="debugging mode, check for overlaps", required=False, action="store_true")
58parser.add_argument("-D", "--display", dest="eventDisplay", help="store trajectories", required=False, action="store_true")
59parser.add_argument("--EmuDet","--nuTargetActive",dest="nuTargetPassive",help="activate emulsiondetector", required=False,action="store_false")
60parser.add_argument("--NagoyaEmu","--useNagoyaEmulsions",dest="useNagoyaEmulsions",help="use bricks of 57 Nagoya emulsion films instead of 60 Slavich", required=False,action="store_true")
61parser.add_argument("-y", dest="year", help="specify the year to generate the respective TI18 detector setup", required=False, type=int, default=2024)
62
63options = parser.parse_args()
64
65# user hook
66userTask = False
67
68class MyTask(ROOT.FairTask):
69 "user task"
70
71 def Exec(self,opt):
72 ioman = ROOT.FairRootManager.Instance()
73 MCTracks = ioman.GetObject("MCTrack")
74 print('Hello',opt,MCTracks.GetEntries())
75 fMC = ROOT.TVirtualMC.GetMC()
76 if MCTracks.GetEntries()>100: fMC.StopRun()
77
78checking4overlaps = False
79if options.debug: checking4overlaps = True
80
81if options.pythia8: simEngine = "Pythia8"
82if options.pg: simEngine = "PG"
83if options.genie: simEngine = "Genie"
84if options.ntuple: simEngine = "Ntuple"
85if options.muonback: simEngine = "MuonBack"
86if options.nuage: simEngine = "Nuage"
87if options.mudis: simEngine = "muonDIS"
88
89if options.inputFile:
90 if options.inputFile == "none": options.inputFile = None
91 inputFile = options.inputFile
92 defaultInputFile = False
93
94if simEngine == "Genie" and defaultInputFile:
95 print('GENIE input file missing, exit')
96 sys.exit()
97if simEngine == "muonDIS" and defaultInputFile:
98 print('input file required if simEngine = muonDIS. Example:')
99 print("/eos/experiment/sndlhc/MonteCarlo/Pythia6/MuonDIS /muonDis_XXXX.root")
100 print(" XXXX = run+cycle*100+k, k: 0 or 1000 for mu+ or mu-, c: number of cycles: 10 events per incoming muon in each cycle, run: 1...10")
101 print(" c = 0 - 2: mu->proton, 5 - 7: mu->neutron")
102 sys.exit()
103if simEngine == "Nuage" and not inputFile:
104 inputFile = 'Numucc.root'
105
106if (simEngine == "Ntuple") and defaultInputFile :
107 print('input file required if simEngine = Ntuple or MuonBack. Examples:')
108 print ("crossing angle up: /eos/experiment/sndlhc/MonteCarlo/FLUKA/muons_up/version1/unit30_Nm.root (unit30_Pm.root)")
109 print ("crossing angle down: /eos/experiment/sndlhc/MonteCarlo/FLUKA/muons_down/muons_VCdown_IR1-LHC.root")
110 sys.exit()
111
112print("SND@LHC setup for",simEngine,"to produce",options.nEvents,"events")
113
114ROOT.gRandom.SetSeed(options.theSeed) # this should be propagated via ROOT to Pythia8 and Geant4VMC
115shipRoot_conf.configure(0) # load basic libraries, prepare atexit for python
116
117if options.testbeam:
118 snd_geo = ConfigRegistry.loadpy("$SNDSW_ROOT/geometry/sndLHC_H6geom_config.py")
119elif options.testbeam2023:
120 snd_geo = ConfigRegistry.loadpy("$SNDSW_ROOT/geometry/sndLHC_HXgeom_config.py",
121 tb_2023_mc = options.testbeam2023)
122else:
123 snd_geo = ConfigRegistry.loadpy("$SNDSW_ROOT/geometry/sndLHC_geom_config.py",
124 nuTargetPassive = options.nuTargetPassive,
125 useNagoyaEmulsions = options.useNagoyaEmulsions,
126 year=options.year)
127
128if simEngine == "PG": tag = simEngine + "_"+str(options.pID)+"-"+mcEngine
129else: tag = simEngine+"-"+mcEngine
130
131if not os.path.exists(options.outputDir):
132 os.makedirs(options.outputDir)
133if options.boostFactor>1:
134 tag+='_boost'+str(options.boostFactor)
135outFile = "%s/sndLHC.%s.root" % (options.outputDir, tag)
136
137# rm older files !!!
138for x in os.listdir(options.outputDir):
139 if not x.find(tag)<0: os.system("rm %s/%s" % (options.outputDir, x) )
140# Parameter file name
141parFile="%s/ship.params.%s.root" % (options.outputDir, tag)
142
143# In general, the following parts need not be touched, except for user task
144# ========================================================================
145
146# -----Timer--------------------------------------------------------
147timer = ROOT.TStopwatch()
148timer.Start()
149# ------------------------------------------------------------------------
150# -----Create simulation run----------------------------------------
151run = ROOT.FairRunSim()
152run.SetName(mcEngine) # Transport engine
153run.SetSink(ROOT.FairRootFileSink(outFile)) # Output file
154run.SetUserConfig("g4Config.C") # user configuration file default g4Config.C
155rtdb = run.GetRuntimeDb()
156# add user task
157if userTask:
158 userTask = MyTask()
159 run.AddTask(userTask)
160
161# -----Create geometry----------------------------------------------
162import shipLHC_conf as sndDet_conf
163modules = sndDet_conf.configure(run,snd_geo)
164
165# -----Create PrimaryGenerator--------------------------------------
166primGen = ROOT.FairPrimaryGenerator()
167
168# -----Particle Gun-----------------------
169if simEngine == "PG":
170 myPgun = ROOT.FairBoxGenerator(options.pID,1)
171 myPgun.SetPRange(options.Estart,options.Eend)
172 myPgun.SetPhiRange(0, 360) # // Azimuth angle range [degree]
173 myPgun.SetThetaRange(0,0) # // Polar angle in lab system range [degree]
174 if options.multiplePGSources:
175 # multiple PG sources in the x-y plane; z is always the same!
176 myPgun.SetBoxXYZ(options.EVx*u.cm,
177 options.EVy*u.cm,
178 (options.EVx+options.Dx)*u.cm,
179 (options.EVy+options.Dy)*u.cm,
180 options.EVz*u.cm)
181 else:
182 # point source
183 myPgun.SetXYZ(options.EVx*u.cm, options.EVy*u.cm, options.EVz*u.cm)
184 primGen.AddGenerator(myPgun)
185 # To generate particle guns along the z axis, create z *target* layers with a set step
186 # For an **unknown** reason simply setting target z thickness doesn't produce the expected result
187 if options.multiplePGSources:
188 targetZpos = np.array(np.arange(options.nZSlices)*options.zSliceStep*u.cm, dtype='d')
189 primGen.SetMultTarget(len(targetZpos), targetZpos, 0*u.cm) # dummy thickness set to 0
190 print('===> Setting particle gun sources starting at (x1,y1,z1)='
191 f'({options.EVx},{options.EVy},{options.EVz})[cm × cm × cm] \n'
192 f'with a uniform x-y spread of (Dx,Dy)=({options.Dx},{options.Dy})[cm × cm]'
193 f' and {options.nZSlices} z slices in steps of {options.zSliceStep}[cm].')
194 ROOT.FairLogger.GetLogger().SetLogScreenLevel("WARNING") # otherwise stupid printout for each event
195# -----muon DIS Background------------------------
196if simEngine == "muonDIS":
197 ut.checkFileExists(inputFile)
198 primGen.SetTarget(0., 0.)
199 DISgen = ROOT.MuDISGenerator()
200 mu_start, mu_end = (-3.7-2.0)*u.m , -0.3*u.m # tunnel wall -30cm in front of SND
201 DISgen.SetPositions(0, mu_start, mu_end)
202 DISgen.Init(inputFile,options.firstEvent)
203 primGen.AddGenerator(DISgen)
204 options.nEvents = min(options.nEvents,DISgen.GetNevents())
205 inactivateMuonProcesses = True # avoid unwanted hadronic events of "incoming" muon flying backward
206 print('MuDIS position info input=',mu_start, mu_end)
207 print('Generate ',options.nEvents,' with DIS input', ' first event',options.firstEvent)
208
209# -----neutrino interactions from nuage------------------------
210if simEngine == "Nuage":
211 primGen.SetTarget(0., 0.)
212 Nuagegen = ROOT.NuageGenerator()
213 Nuagegen.EnableExternalDecayer(1) #with 0 external decayer is disable, 1 is enabled
214 Nuagegen.SetPositions(0., -snd_geo.EmulsionDet.zdim/2, snd_geo.EmulsionDet.zdim/2, -snd_geo.EmulsionDet.xdim/2, snd_geo.EmulsionDet.xdim/2, -snd_geo.EmulsionDet.ydim/2, snd_geo.EmulsionDet.ydim/2)
215 ut.checkFileExists(inputFile)
216 Nuagegen.Init(inputFile,options.firstEvent)
217 primGen.AddGenerator(Nuagegen)
218 options.nEvents = min(options.nEvents,Nuagegen.GetNevents())
219 run.SetPythiaDecayer("DecayConfigNuAge.C")
220 print('Generate ',options.nEvents,' with Nuage input', ' first event',options.firstEvent)
221
222# -----neutrino interactions from GENIE------------------------
223if simEngine=="Genie":
224 ut.checkFileExists(inputFile)
225 primGen.SetTarget(0., 0.) # do not interfere with GenieGenerator
226 Geniegen = ROOT.GenieGenerator()
227 Geniegen.SetGenerationOption(options.genie - 1) # 0 standard, 1 FLUKA,2 Pythia
228 Geniegen.Init(inputFile,options.firstEvent)
229 Geniegen.SetCrossingAngle(150e-6) #used only in option 2
230
231 # Neutrino vertex generation range in z:
232 # Tolerance for neutrino vertex generation range. Mostly to account for tilt in geometry alignment. Take difference in z coordinate of vertical fibres of around 0.5 cm over the fibre length, 39 cm. Assume maximum difference in z is 1 m * 0.5/39.
233 tolerance_vtx_z = 1*u.m * 0.5/39
234 # From first veto bar
235 neutrino_vtx_start_z = snd_geo.MuFilter.Veto1Dy - snd_geo.MuFilter.VetoBarZ/2. - tolerance_vtx_z
236 # To last Scifi plane
237 neutrino_vtx_end_z = snd_geo.Scifi.Ypos4 + snd_geo.Scifi.zdim/2. + tolerance_vtx_z
238
239 Geniegen.SetPositions(-480*u.m, neutrino_vtx_start_z, neutrino_vtx_end_z)
240
241 Geniegen.SetDeltaE_Matching_FLUKAGenie(10.) #energy range for the search of a GENIE interaction with similar energy of FLUKA neutrino
242 primGen.AddGenerator(Geniegen)
243 options.nEvents = min(options.nEvents,Geniegen.GetNevents())
244 run.SetPythiaDecayer('DecayConfigPy8.C')
245 print('Generate ',options.nEvents,' with Genie input for Ship@LHC', ' first event',options.firstEvent)
246
247if simEngine == "Ntuple":
248 ut.checkFileExists(inputFile)
249 Ntuplegen = ROOT.NtupleGenerator_FLUKA()
250 Ntuplegen.SetZ(snd_geo.Floor.z)
251 Ntuplegen.Init(inputFile,options.firstEvent)
252 primGen.AddGenerator(Ntuplegen)
253 options.nEvents = min(options.nEvents,Ntuplegen.GetNevents())
254
255if simEngine == "MuonBack":
256# reading muon tracks from FLUKA
257 fileType = ut.checkFileExists(inputFile)
258 if fileType == 'tree':
259 # 2018 background production
260 primGen.SetTarget(snd_geo.target.z0+70.845*u.m,0.)
261 else:
262 primGen.SetTarget(snd_geo.target.z0+50*u.m,0.)
263 #
264 MuonBackgen = ROOT.MuonBackGenerator()
265 # MuonBackgen.FollowAllParticles() # will follow all particles after hadron absorber, not only muons
266 MuonBackgen.Init(inputFile,options.firstEvent,options.phiRandom)
267 primGen.AddGenerator(MuonBackgen)
268 options.nEvents = min(options.nEvents,MuonBackgen.GetNevents())
269 MCTracksWithHitsOnly = True # otherwise, output file becomes too big
270 print('Process ',options.nEvents,' from input file, with Phi random=',options.phiRandom, ' with MCTracksWithHitsOnly',MCTracksWithHitsOnly)
271
272if options.ecut > 0:
273 modules['Floor'].SetEmin(options.ecut)
274 modules['Floor'].SetZmax(options.zmax)
275
276#
277run.SetGenerator(primGen)
278# ------------------------------------------------------------------------
279if options.followMuon :
280 if 'Veto' in modules:
281 options.fastMuon = True
282 modules['Veto'].SetFollowMuon()
283 if 'Floor' in modules:
284 modules['Floor'].MakeSensitive()
285 print('make floor sensitive')
286if options.fastMuon :
287 if 'Veto' in modules: modules['Veto'].SetFastMuon()
288 elif 'Floor' in modules:
289 modules['Floor'].SetFastMuon()
290 modules['Floor'].SetZmax(options.zmax)
291 print('transport only-muons up to z=',options.zmax)
292# ------------------------------------------------------------------------
293#---Store the visualiztion info of the tracks, this make the output file very large!!
294#--- Use it only to display but not for production!
295if options.eventDisplay: run.SetStoreTraj(ROOT.kTRUE)
296else: run.SetStoreTraj(ROOT.kFALSE)
297
298
299
300# -----Initialize simulation run------------------------------------
301run.Init()
302gMC = ROOT.TVirtualMC.GetMC()
303fStack = gMC.GetStack()
304if MCTracksWithHitsOnly:
305 fStack.SetMinPoints(1)
306 fStack.SetEnergyCut(-100.*u.MeV)
307elif MCTracksWithEnergyCutOnly:
308 fStack.SetMinPoints(-1)
309 fStack.SetEnergyCut(100.*u.MeV)
310elif MCTracksWithHitsOrEnergyCut:
311 fStack.SetMinPoints(1)
312 fStack.SetEnergyCut(100.*u.MeV)
313elif options.deepCopy:
314 fStack.SetMinPoints(0)
315 fStack.SetEnergyCut(0.*u.MeV)
316
317#
318if options.boostFactor > 1:
319 ROOT.gROOT.ProcessLine('#include "Geant4/G4ProcessTable.hh"')
320 ROOT.gROOT.ProcessLine('#include "Geant4/G4MuBremsstrahlung.hh"')
321 ROOT.gROOT.ProcessLine('#include "Geant4/G4GammaConversionToMuons.hh"')
322 ROOT.gROOT.ProcessLine('#include "Geant4/G4MuPairProduction.hh"')
323 ROOT.gROOT.ProcessLine('#include "Geant4/G4AnnihiToMuPair.hh"')
324 ROOT.gROOT.ProcessLine('#include "Geant4/G4MuonToMuonPairProduction.hh"')
325 ROOT.gROOT.ProcessLine('#include "Geant4/G4MuonPlus.hh"')
326 ROOT.gROOT.ProcessLine('#include "Geant4/G4MuonMinus.hh"')
327
328 gProcessTable = ROOT.G4ProcessTable.GetProcessTable()
329 # only muon interaction
330 # procBrems = gProcessTable.FindProcess(ROOT.G4String('muBrems'),ROOT.G4String('mu+'))
331 # muPairProd = gProcessTable.FindProcess(ROOT.G4String('muPairProd'),ROOT.G4String('mu+'))
332 # muPairProd.SetCrossSectionBiasingFactor(options.boostFactor)
333 # procBrems.SetCrossSectionBiasingFactor(options.boostFactor)
334 # muon pair production
335 gammaToMuPair = gProcessTable.FindProcess(ROOT.G4String('GammaToMuPair'),ROOT.G4String('gamma'))
336 gammaToMuPair.SetCrossSecFactor(options.boostFactor)
337 AnnihiToMuPair = gProcessTable.FindProcess(ROOT.G4String('AnnihiToMuPair'),ROOT.G4String('e+'))
338 AnnihiToMuPair.SetCrossSecFactor(options.boostFactor)
339 MuonToMuonPair = gProcessTable.FindProcess(ROOT.G4String('muToMuonPairProd'),ROOT.G4String('mu+'))
340 MuonToMuonPair.SetCrossSectionBiasingFactor(options.boostFactor)
341
342 mygMC = ROOT.TGeant4.GetMC()
343 if options.debug:
344 mygMC.ProcessGeantCommand("/run/particle/dumpOrderingParam")
345 mygMC.ProcessGeantCommand("/particle/select mu+")
346 mygMC.ProcessGeantCommand("/particle/process/dump")
347 mygMC.ProcessGeantCommand("/particle/select gamma")
348 mygMC.ProcessGeantCommand("/particle/process/dump")
349 mygMC.ProcessGeantCommand("/particle/select e+")
350 mygMC.ProcessGeantCommand("/particle/process/dump")
351#
352if options.enhancePiKaDecay:
353 ROOT.gROOT.ProcessLine('#include "Geant4/G4ParticleTable.hh"')
354 ROOT.gROOT.ProcessLine('#include "Geant4/G4DecayTable.hh"')
355 ROOT.gROOT.ProcessLine('#include "Geant4/G4PhaseSpaceDecayChannel.hh"')
356 pt = ROOT.G4ParticleTable.GetParticleTable()
357 for pid in [211,-211,321,-321]:
358 particleG4 = pt.FindParticle(pid)
359 lt = particleG4.GetPDGLifeTime()
360 particleG4.SetPDGLifeTime(lt/options.enhancePiKaDecay)
361 print('### pion kaon lifetime decreased by the factor:',options.enhancePiKaDecay)
362
363if inactivateMuonProcesses :
364 ROOT.gROOT.ProcessLine('#include "Geant4/G4ProcessTable.hh"')
365 mygMC = ROOT.TGeant4.GetMC()
366 mygMC.ProcessGeantCommand("/process/inactivate muPairProd")
367 mygMC.ProcessGeantCommand("/process/inactivate muBrems")
368 #mygMC.ProcessGeantCommand("/process/inactivate muIoni") Temporary fix for DIS Simulations (incoming and outgoing muon hits)
369 mygMC.ProcessGeantCommand("/process/inactivate muonNuclear")
370 mygMC.ProcessGeantCommand("/particle/select mu+")
371 mygMC.ProcessGeantCommand("/particle/process/dump")
372 gProcessTable = ROOT.G4ProcessTable.GetProcessTable()
373 procmu = gProcessTable.FindProcess(ROOT.G4String('muIoni'),ROOT.G4String('mu+'))
374 procmu.SetVerboseLevel(2)
375
376if options.debug: ROOT.fair.Logger.SetConsoleSeverity("debug")
377# -----Start run----------------------------------------------------
378run.Run(options.nEvents)
379# -----Runtime database---------------------------------------------
380kParameterMerged = ROOT.kTRUE
381parOut = ROOT.FairParRootFileIo(kParameterMerged)
382parOut.open(parFile)
383rtdb.setOutput(parOut)
384rtdb.saveOutput()
385rtdb.printParamContexts()
386getattr(rtdb,"print")()
387# ------------------------------------------------------------------------
388geoFile = "%s/geofile_full.%s.root" % (options.outputDir, tag)
389run.CreateGeometryFile(geoFile)
390# save detector parameters dictionary in geofile
391import saveBasicParameters
392saveBasicParameters.execute(geoFile,snd_geo)
393
394# ------------------------------------------------------------------------
395# If using GENIE option 4 (geometry driver) copy GST TTree to the
396# output file. This will make it easy to access the FLUKA variables for
397# each neutrino event.
398if options.genie == 4 :
399
400 f_input = ROOT.TFile(inputFile)
401 gst = f_input.gst
402
403 selection_string = "(Entry$ >= "+str(options.firstEvent)+")"
404 if (options.firstEvent + options.nEvents) < gst.GetEntries() :
405 selection_string += "&&(Entry$ < "+str(options.firstEvent + options.nEvents)+")"
406
407 # Reopen output file
408 f_output = ROOT.TFile(outFile, "UPDATE")
409
410 # Copy only the events used in this run
411 gst_copy = gst.CopyTree(selection_string)
412 gst_copy.Write()
413
414 f_input.Close()
415 f_output.Close()
416
417# -----Finish-------------------------------------------------------
418timer.Stop()
419rtime = timer.RealTime()
420ctime = timer.CpuTime()
421print(' ')
422print("Macro finished succesfully.")
423
424print("Output file is ", outFile)
425print("Geometry file is ",geoFile)
426print("Real time ",rtime, " s, CPU time ",ctime,"s")
427
428# ------------------------------------------------------------------------
429def checkOverlaps(removeScifi=False):
430 sGeo = ROOT.gGeoManager
431 if removeScifi:
432 for n in range(1,6):
433 Hscifi = sGeo.FindVolumeFast('ScifiVolume'+str(n))
434 removalList = []
435 for x in Hscifi.GetNodes():
436 if x.GetName().find('Scifi')==0: removalList.append(x)
437 for x in removalList: Hscifi.RemoveNode(x)
438 sGeo.SetNmeshPoints(10000)
439 sGeo.CheckOverlaps(0.1) # 1 micron takes 5minutes
440 sGeo.PrintOverlaps()
441# check subsystems in more detail
442 for x in sGeo.GetTopNode().GetNodes():
443 x.CheckOverlaps(0.0001)
444 sGeo.PrintOverlaps()
445
447 # after /run/initialize, but prints warning messages, problems with TGeo volume
448 mygMC = ROOT.TGeant4.GetMC()
449 mygMC.ProcessGeantCommand("/geometry/test/recursion_start 0")
450 mygMC.ProcessGeantCommand("/geometry/test/recursion_depth 2")
451 mygMC.ProcessGeantCommand("/geometry/test/run")
452
453# checking for overlaps
454if checking4overlaps:
Exec(self, opt)
Definition run_simSND.py:71
checkOverlaps(removeScifi=False)
checkOverlapsWithGeant4()
execute(f, ox, name='ShipGeo')
configure(darkphoton=None)