CMS 3D CMS Logo

haddnano.py
Go to the documentation of this file.
1 #!/bin/env python
2 import ROOT
3 import numpy
4 import sys
5 
6 if len(sys.argv) < 3:
7  print("Syntax: haddnano.py out.root input1.root input2.root ...")
8 ofname = sys.argv[1]
9 files = sys.argv[2:]
10 
11 
12 def zeroFill(tree, brName, brObj, allowNonBool=False):
13  # typename: (numpy type code, root type code)
14  branch_type_dict = {'Bool_t': ('?', 'O'), 'Float_t': ('f4', 'F'), 'UInt_t': (
15  'u4', 'i'), 'Long64_t': ('i8', 'L'), 'Double_t': ('f8', 'D')}
16  brType = brObj.GetLeaf(brName).GetTypeName()
17  if (not allowNonBool) and (brType != "Bool_t"):
18  print(("Did not expect to back fill non-boolean branches", tree, brName, brObj.GetLeaf(br).GetTypeName()))
19  else:
20  if brType not in branch_type_dict:
21  raise RuntimeError('Impossible to backfill branch of type %s' % brType)
22  buff = numpy.zeros(1, dtype=numpy.dtype(branch_type_dict[brType][0]))
23  b = tree.Branch(brName, buff, brName + "/" +
24  branch_type_dict[brType][1])
25  # be sure we do not trigger flushing
26  b.SetBasketSize(tree.GetEntries() * 2)
27  for x in range(0, tree.GetEntries()):
28  b.Fill()
29  b.ResetAddress()
30 
31 
32 fileHandles = []
33 goFast = True
34 for fn in files:
35  print("Adding file" + str(fn))
36  fileHandles.append(ROOT.TFile.Open(fn))
37  if fileHandles[-1].GetCompressionSettings() != fileHandles[0].GetCompressionSettings():
38  goFast = False
39  print("Disabling fast merging as inputs have different compressions")
40 of = ROOT.TFile(ofname, "recreate")
41 if goFast:
42  of.SetCompressionSettings(fileHandles[0].GetCompressionSettings())
43 of.cd()
44 
45 for e in fileHandles[0].GetListOfKeys():
46  name = e.GetName()
47  print("Merging" + str(name))
48  obj = e.ReadObj()
49  cl = ROOT.TClass.GetClass(e.GetClassName())
50  inputs = ROOT.TList()
51  isTree = obj.IsA().InheritsFrom(ROOT.TTree.Class())
52  if isTree:
53  obj = obj.CloneTree(-1, "fast" if goFast else "")
54  branchNames = set([x.GetName() for x in obj.GetListOfBranches()])
55  for fh in fileHandles[1:]:
56  otherObj = fh.GetListOfKeys().FindObject(name).ReadObj()
57  inputs.Add(otherObj)
58  if isTree and obj.GetName() == 'Events':
59  otherObj.SetAutoFlush(0)
60  otherBranches = set([x.GetName()
61  for x in otherObj.GetListOfBranches()])
62  missingBranches = list(branchNames - otherBranches)
63  additionalBranches = list(otherBranches - branchNames)
64  print("missing: " + str(missingBranches) + "\n Additional:" + str(additionalBranches))
65  for br in missingBranches:
66  # fill "Other"
67  zeroFill(otherObj, br, obj.GetListOfBranches().FindObject(br))
68  for br in additionalBranches:
69  # fill main
70  branchNames.add(br)
71  zeroFill(obj, br, otherObj.GetListOfBranches().FindObject(br))
72  # merge immediately for trees
73  if isTree and obj.GetName() == 'Runs':
74  otherObj.SetAutoFlush(0)
75  otherBranches = set([x.GetName()
76  for x in otherObj.GetListOfBranches()])
77  missingBranches = list(branchNames - otherBranches)
78  additionalBranches = list(otherBranches - branchNames)
79  print("missing: " + str(missingBranches) + "\n Additional:" + str(additionalBranches))
80  for br in missingBranches:
81  # fill "Other"
82  zeroFill(otherObj, br, obj.GetListOfBranches(
83  ).FindObject(br), allowNonBool=True)
84  for br in additionalBranches:
85  # fill main
86  branchNames.add(br)
87  zeroFill(obj, br, otherObj.GetListOfBranches(
88  ).FindObject(br), allowNonBool=True)
89  # merge immediately for trees
90  if isTree:
91  obj.Merge(inputs, "fast" if goFast else "")
92  inputs.Clear()
93 
94  if isTree:
95  obj.Write()
96  elif obj.IsA().InheritsFrom(ROOT.TH1.Class()):
97  obj.Merge(inputs)
98  obj.Write()
99  elif obj.IsA().InheritsFrom(ROOT.TObjString.Class()):
100  for st in inputs:
101  if st.GetString() != obj.GetString():
102  print("Strings are not matching")
103  obj.Write()
104  else:
105  print("Cannot handle " + str(obj.IsA().GetName()))
def zeroFill(tree, brName, brObj, allowNonBool=False)
Definition: haddnano.py:12
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
#define str(s)