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