Fork me on GitHub

source: git/readers/DelphesProMC.py@ 154d5da

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 154d5da was 154d5da, checked in by pavel <pavel@…>, 11 years ago

first working version of DelphesProMC.py

  • Property mode set to 100644
File size: 5.6 KB
Line 
1import sys
2
3import urllib2
4import cStringIO
5import zipfile
6
7import ROOT
8
9import ProMCHeader_pb2
10import ProMC_pb2
11
12################################################################################
13# HttpFile class is written by Eli Carter (retracile)
14# http://stackoverflow.com/a/7852229
15
16class HttpFile(object):
17 def __init__(self, url):
18 self.url = url
19 self.offset = 0
20 self.length = -1
21
22 def size(self):
23 if self.length < 0:
24 f = urllib2.urlopen(self.url)
25 self.length = int(f.headers["Content-length"])
26 return self.length
27
28 def read(self, count = -1):
29 req = urllib2.Request(self.url)
30 if count < 0:
31 end = self.size() - 1
32 else:
33 end = self.offset + count - 1
34 req.headers['Range'] = "bytes=%s-%s" % (self.offset, end)
35 f = urllib2.urlopen(req)
36 data = f.read()
37 # FIXME: should check that we got the range expected, etc.
38 chunk = len(data)
39 if count >= 0:
40 assert chunk == count
41 self.offset += chunk
42 return data
43
44 def seek(self, offset, whence = 0):
45 if whence == 0:
46 self.offset = offset
47 elif whence == 1:
48 self.offset += offset
49 elif whence == 2:
50 self.offset = self.size() + offset
51 else:
52 raise Exception("Invalid whence")
53
54 def tell(self):
55 return self.offset
56
57################################################################################
58
59def ConvertInput(name, momentumUnit, lengthUnit, branch, factory,
60 allParticleOutputArray, stableParticleOutputArray, partonOutputArray):
61
62 pdg = ROOT.TDatabasePDG.Instance()
63
64 data = ProMC_pb2.ProMCEvent()
65 data.ParseFromString(zip.read(name))
66
67 event = data.event
68 particles = data.particles
69
70 element = branch.NewEntry()
71
72 element.Number = event.Number[0] if event.Number else 0
73 element.ProcessID = event.Process_ID[0] if event.Process_ID else 0
74 element.MPI = event.MPI[0] if event.MPI else 0
75# element.Weight = event.Weight[0] if event.Weight else 0
76 element.Scale = event.Scale[0] if event.Scale else 0
77 element.AlphaQED = event.Alpha_QED[0] if event.Alpha_QED else 0
78 element.AlphaQCD = event.Alpha_QCD[0] if event.Alpha_QCD else 0
79 element.ID1 = event.ID1[0] if event.ID1 else 0
80 element.ID2 = event.ID2[0] if event.ID2 else 0
81 element.X1 = event.X1[0] if event.X1 else 0
82 element.X2 = event.X2[0] if event.X2 else 0
83 element.ScalePDF = event.Scale_PDF[0] if event.Scale_PDF else 0
84 element.PDF1 = event.PDF1[0] if event.PDF1 else 0
85 element.PDF2 = event.PDF2[0] if event.PDF2 else 0
86
87 for i in range(len(particles.pdg_id)):
88 pid = particles.pdg_id[i]
89 status = particles.status[i]
90 px = particles.Px[i]
91 py = particles.Py[i]
92 pz = particles.Pz[i]
93 mass = particles.mass[i]
94
95 x = particles.X[i]
96 y = particles.Y[i]
97 z = particles.Z[i]
98 t = particles.T[i]
99
100 candidate = factory.NewCandidate()
101
102 candidate.PID = pid
103 pdgCode = ROOT.TMath.Abs(candidate.PID)
104
105 candidate.Status = status
106
107 candidate.M1 = particles.mother1[i]
108 candidate.M2 = particles.mother2[i]
109
110 candidate.D1 = particles.daughter1[i]
111 candidate.D2 = particles.daughter2[i]
112
113 candidate.Mass = mass
114
115 candidate.Momentum.SetXYZM(px, py, pz, mass)
116
117 candidate.Position.SetXYZT(x, y, z, t)
118
119 pdgParticle = pdg.GetParticle(pid)
120 if pdgParticle == None:
121 candidate.Charge = -999
122 continue
123 else:
124 candidate.Charge = int(pdgParticle.Charge()/3.0)
125
126 allParticleOutputArray.Add(candidate)
127
128 if status == 1:
129 stableParticleOutputArray.Add(candidate)
130 elif pdgCode <= 5 or pdgCode == 21 or pdgCode == 15:
131 partonOutputArray.Add(candidate)
132
133################################################################################
134
135if len(sys.argv) < 2:
136 print " Usage: DelphesProMC.py config_file output_file input_file(s)"
137 sys.exit(1)
138
139ROOT.gSystem.Load('libDelphes')
140
141outputFile = ROOT.TFile(sys.argv[2], "CREATE")
142
143if not outputFile.IsOpen():
144 print "** ERROR: can't open", sys.argv[2]
145 sys.exit(1)
146
147treeWriter = ROOT.ExRootTreeWriter(outputFile, "Delphes")
148
149branchEvent = treeWriter.NewBranch("Event", ROOT.HepMCEvent.Class())
150
151confReader = ROOT.ExRootConfReader()
152confReader.ReadFile(sys.argv[1])
153
154modularDelphes = ROOT.Delphes("Delphes")
155modularDelphes.SetConfReader(confReader)
156modularDelphes.SetTreeWriter(treeWriter)
157
158factory = modularDelphes.GetFactory()
159allParticleOutputArray = modularDelphes.ExportArray("allParticles")
160stableParticleOutputArray = modularDelphes.ExportArray("stableParticles")
161partonOutputArray = modularDelphes.ExportArray("partons")
162
163modularDelphes.InitTask()
164
165for fileName in sys.argv[3:]:
166 print "** Reading", fileName
167
168 if fileName.startswith("http://"):
169 file = HttpFile(fileName)
170 else:
171 file = open(fileName)
172
173 zip = zipfile.ZipFile(file)
174
175 numberOfEvents = len(zip.namelist())
176 if numberOfEvents <= 0: continue
177
178 progressBar = ROOT.ExRootProgressBar(numberOfEvents - 1)
179
180 # retrive information from the header file
181 header = ProMCHeader_pb2.ProMCHeader()
182 header.ParseFromString(zip.read("header"))
183 momentumUnit = float(header.MomentumUnit)
184 lengthUnit = float(header.LengthUnit)
185
186 modularDelphes.Clear()
187 treeWriter.Clear()
188 eventCounter = 0
189 for name in zip.namelist():
190 eventCounter += 1
191 if not name.isdigit(): continue
192
193 ConvertInput(name, momentumUnit, lengthUnit, branchEvent, factory,
194 allParticleOutputArray, stableParticleOutputArray, partonOutputArray)
195
196 modularDelphes.ProcessTask()
197
198 treeWriter.Fill()
199
200 modularDelphes.Clear()
201 treeWriter.Clear()
202 progressBar.Update(eventCounter)
203
204 progressBar.Update(eventCounter, eventCounter, ROOT.kTRUE)
205 progressBar.Finish()
206 zip.close()
207
208modularDelphes.FinishTask()
209treeWriter.Write()
210
211del treeWriter
Note: See TracBrowser for help on using the repository browser.