CMS 3D CMS Logo

AutoFillTreeProducer.py
Go to the documentation of this file.
1 from builtins import range
2 from PhysicsTools.Heppy.analyzers.core.TreeAnalyzerNumpy import TreeAnalyzerNumpy
3 from PhysicsTools.Heppy.analyzers.core.AutoHandle import AutoHandle
4 #from ROOT import TriggerBitChecker
7 
8 import six
9 
10 class AutoFillTreeProducer( TreeAnalyzerNumpy ):
11 
12  #-----------------------------------
13  # BASIC TREE PRODUCER
14  #-----------------------------------
15  def __init__(self, cfg_ana, cfg_comp, looperName):
16  super(AutoFillTreeProducer,self).__init__(cfg_ana, cfg_comp, looperName)
17 
18  ## Read whether we want vectors or flat trees
19  self.scalar = not self.cfg_ana.vectorTree
20 
21  ## Read whether we want 4-vectors
22  if not getattr(self.cfg_ana, 'saveTLorentzVectors', False):
23  fourVectorType.removeVariable("p4")
24 
25 
26  self.collections = {}
27  self.globalObjects = {}
28  self.globalVariables = []
29  if hasattr(cfg_ana,"collections"):
30  self.collections.update(cfg_ana.collections)
31  if hasattr(cfg_ana,"globalObjects"):
32  self.globalObjects.update(cfg_ana.globalObjects)
33  if hasattr(cfg_ana,"globalVariables"):
34  self.globalVariables=cfg_ana.globalVariables[:]
35 
36  def beginLoop(self, setup) :
37  super(AutoFillTreeProducer, self).beginLoop(setup)
38 
39  def declareHandles(self):
40  super(AutoFillTreeProducer, self).declareHandles()
41 # self.handles['TriggerResults'] = AutoHandle( ('TriggerResults','','HLT'), 'edm::TriggerResults' )
42  self.mchandles['GenInfo'] = AutoHandle( ('generator','',''), 'GenEventInfoProduct' )
43  for k,v in six.iteritems(self.collections):
44  if isinstance(v, tuple) and isinstance(v[0], AutoHandle):
45  self.handles[k] = v[0]
46 
47  def declareCoreVariables(self, tr, isMC):
48  """Here we declare the variables that we always want and that are hard-coded"""
49  tr.var('run', int, storageType="i")
50  tr.var('lumi', int, storageType="i")
51  tr.var('evt', int, storageType="l")
52  tr.var('isData', int)
53 
54  # self.triggerBitCheckers = []
55  # if hasattr(self.cfg_ana, 'triggerBits'):
56  # for T, TL in six.iteritems(self.cfg_ana.triggerBits):
57  # trigVec = ROOT.vector(ROOT.string)()
58  # for TP in TL:
59  # trigVec.push_back(TP)
60  # tr.var( 'HLT_'+T, int )
61 # self.triggerBitCheckers.append( (T, TriggerBitChecker(trigVec)) )
62 
63  if not isMC:
64  tr.var('intLumi', int, storageType="i")
65 
66  if isMC:
67  ## cross section
68  tr.var('xsec', float)
69  ## PU weights
70  tr.var("puWeight")
71  ## number of true interactions
72  tr.var("nTrueInt")
73  ## generator weight
74  tr.var("genWeight")
75  ## PDF weights
76  self.pdfWeights = []
77  if hasattr(self.cfg_ana, "PDFWeights") and len(self.cfg_ana.PDFWeights) > 0:
78  self.pdfWeights = self.cfg_ana.PDFWeights
79  for (pdf,nvals) in self.pdfWeights:
80  if self.scalar:
81  for i in range(nvals): tr.var('pdfWeight_%s_%d' % (pdf,i))
82  else:
83  tr.vector('pdfWeight_%s' % pdf, nvals)
84 
85  def declareVariables(self,setup):
86  isMC = self.cfg_comp.isMC
87  tree = self.tree
88  self.declareCoreVariables(tree, isMC)
89 
90  if not hasattr(self.cfg_ana,"ignoreAnalyzerBookings") or not self.cfg_ana.ignoreAnalyzerBookings :
91  #import variables declared by the analyzers
92  if hasattr(setup,"globalVariables"):
93  self.globalVariables+=setup.globalVariables
94  if hasattr(setup,"globalObjects"):
95  self.globalObjects.update(setup.globalObjects)
96  if hasattr(setup,"collections"):
97  self.collections.update(setup.collections)
98 
99  for v in self.globalVariables:
100  v.makeBranch(tree, isMC)
101  for o in six.itervalues(self.globalObjects):
102  o.makeBranches(tree, isMC)
103  for c in six.itervalues(self.collections):
104  if isinstance(c, tuple): c = c[-1]
105  if self.scalar:
106  c.makeBranchesScalar(tree, isMC)
107  else:
108  c.makeBranchesVector(tree, isMC)
109 
110  def fillCoreVariables(self, tr, event, isMC):
111  """Here we fill the variables that we always want and that are hard-coded"""
112  tr.fill('run', event.input.eventAuxiliary().id().run())
113  tr.fill('lumi',event.input.eventAuxiliary().id().luminosityBlock())
114  tr.fill('evt', event.input.eventAuxiliary().id().event())
115  tr.fill('isData', 0 if isMC else 1)
116 
117 # triggerResults = self.handles['TriggerResults'].product()
118 # for T,TC in self.triggerBitCheckers:
119 # tr.fill("HLT_"+T, TC.check(event.object(), triggerResults))
120 
121  if not isMC:
122  tr.fill('intLumi', getattr(self.cfg_comp,'intLumi',1.0))
123 
124  if isMC:
125  ## xsection, if available
126  tr.fill('xsec', getattr(self.cfg_comp,'xSection',1.0))
127  ## PU weights, check if a PU analyzer actually filled it
128  if hasattr(event,"nPU"):
129  tr.fill("nTrueInt", event.nPU)
130  tr.fill("puWeight", event.puWeight)
131  else :
132  tr.fill("nTrueInt", -1)
133  tr.fill("puWeight", 1.0)
134 
135  tr.fill("genWeight", self.mchandles['GenInfo'].product().weight())
136  ## PDF weights
137  if hasattr(event,"pdfWeights") :
138  for (pdf,nvals) in self.pdfWeights:
139  if len(event.pdfWeights[pdf]) != nvals:
140  raise RuntimeError("PDF lenght mismatch for %s, declared %d but the event has %d" % (pdf,nvals,event.pdfWeights[pdf]))
141  if self.scalar:
142  for i,w in enumerate(event.pdfWeights[pdf]):
143  tr.fill('pdfWeight_%s_%d' % (pdf,i), w)
144  else:
145  tr.vfill('pdfWeight_%s' % pdf, event.pdfWeights[pdf])
146 
147  def process(self, event):
148  if hasattr(self.cfg_ana,"filter") :
149  if not self.cfg_ana.filter(event) :
150  return True #do not stop processing, just filter myself
151  self.readCollections( event.input)
152  self.fillTree(event)
153 
154  def fillTree(self, event, resetFirst=True):
155  isMC = self.cfg_comp.isMC
156  if resetFirst: self.tree.reset()
157 
158  self.fillCoreVariables(self.tree, event, isMC)
159 
160  for v in self.globalVariables:
161  if not isMC and v.mcOnly: continue
162  v.fillBranch(self.tree, event, isMC)
163 
164  for on, o in six.iteritems(self.globalObjects):
165  if not isMC and o.mcOnly: continue
166  o.fillBranches(self.tree, getattr(event, on), isMC)
167 
168  for cn, c in six.iteritems(self.collections):
169  if isinstance(c, tuple) and isinstance(c[0], AutoHandle):
170  if not isMC and c[-1].mcOnly: continue
171  objects = self.handles[cn].product()
172  setattr(event, cn, [objects[i] for i in range(objects.size())])
173  c = c[-1]
174  if not isMC and c.mcOnly: continue
175  if self.scalar:
176  c.fillBranchesScalar(self.tree, getattr(event, cn), isMC)
177  else:
178  c.fillBranchesVector(self.tree, getattr(event, cn), isMC)
179 
180  self.tree.tree.Fill()
181 
182  def getPythonWrapper(self):
183  """
184  This function produces a string that contains a Python wrapper for the event.
185  The wrapper is automatically generated based on the collections and allows the full
186  event contents to be accessed from subsequent Analyzers using e.g.
187 
188  leps = event.selLeptons #is of type selLeptons
189  pt0 = leps[0].pt
190 
191  One just needs to add the EventAnalyzer to the sequence.
192  """
193 
194  isMC = self.cfg_comp.isMC
195 
196  classes = ""
197  anclass = ""
198  anclass += "from PhysicsTools.HeppyCore.framework.analyzer import Analyzer\n"
199  anclass += "class EventAnalyzer(Analyzer):\n"
200  anclass += " def __init__(self, cfg_ana, cfg_comp, looperName):\n"
201  anclass += " super(EventAnalyzer, self).__init__(cfg_ana, cfg_comp, looperName)\n"
202 
203  anclass += " def process(self, event):\n"
204 
205  for cname, coll in self.collections.items():
206  classes += coll.get_py_wrapper_class(isMC)
207  anclass += " event.{0} = {0}.make_array(event)\n".format(coll.name)
208 
209  return classes + "\n" + anclass
210 
scalar
Read whether we want vectors or flat trees.
Definition: weight.py:1
def __init__(self, cfg_ana, cfg_comp, looperName)
Definition: event.py:1