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