CMS 3D CMS Logo

PhaseITreeProducer.py
Go to the documentation of this file.
1 #!/usr/bin/env python3
2 
3 from __future__ import print_function
4 import sys
5 from ROOT import *
6 from array import array
7 from copy import deepcopy
8 
9 gROOT.SetBatch() # don't pop up canvases
10 
11 #get data files
12 
13 def getFileInPath(rfile):
14  import os
15  for dir in os.environ['CMSSW_SEARCH_PATH'].split(":"):
16  if os.path.exists(os.path.join(dir,rfile)): return os.path.join(dir,rfile)
17  return None
18 
19 # Default values
20 inputFileName = "DQM.root"
21 outputFileName = "DQMTree.root"
22 #detIDsFileName = "DATA/detids.dat"
23 detIDsFileName = getFileInPath('DQM/SiStripMonitorClient/data/detids.dat')
24 
26 
27 
28 
29  def __TraverseDirTree(self, dir):
30 
31  currPath = (dir.GetPath().split(":/"))[1]
32 
33  for obj in dir.GetListOfKeys():
34  if not obj.IsFolder():
35  if obj.ReadObjectAny(TClass.GetClass("TH2")):
36  th2 = deepcopy(obj.ReadObj())
37  name = th2.GetName()
38  if 6 < th2.GetNbinsX() < 10 and name.find("per") != -1 and name.find("Lumisection") == -1: #take only module lvl plots
39  print(''.join([dir.GetPath(), '/', name]))
40 
41  # fix when there are plots starting with the same strings in different directories
42  prefix = ""
43  for i in self.dirs:
44  if currPath.startswith(i):
45  prefix = self.dirsAliases[i]
46  break
47  # print(currPath, prefix)
48  th2.SetName(prefix + th2.GetName())
49  self.listOfNumHistograms.append(th2)
50  else:
51  self.__TraverseDirTree(obj.ReadObj())
52 
53  def __GetPartStr(self, isXlowerThanZero, isYlowerThanZero):
54  if isXlowerThanZero and isYlowerThanZero:
55  return "mO"
56  if isXlowerThanZero and isYlowerThanZero == False:
57  return "mI"
58  if isXlowerThanZero == False and isYlowerThanZero:
59  return "pO"
60  if isXlowerThanZero == False and isYlowerThanZero == False:
61  return "pI"
62 
63  def __GetBarrelSector(self, layer, signedLadder, signedModule): #adapted from PixelBarrelName
64  theLadder = abs(signedLadder)
65  theModule = abs(signedModule)
66 
67  sector = 0
68 
69  if layer == 1:
70 
71  if theLadder == 1:
72  if theModule >= 2:
73  return 1
74  else:
75  return 2
76  if theLadder == 2:
77  if theModule >= 3:
78  return 2
79  else:
80  return 3
81  if theLadder == 3:
82  if theModule >= 4:
83  return 3
84  else:
85  return 4
86  if theLadder == 4:
87  if theModule >= 2:
88  return 5
89  else:
90  return 6
91  if theLadder == 5:
92  if theModule >= 3:
93  return 6
94  else:
95  return 7
96  if theLadder == 6:
97  if theModule >= 4:
98  return 7
99  else:
100  return 8
101  # here is used simplified form of assignment, see source file for reference
102  elif layer == 2:
103  i = theLadder // 5
104  sector = i * 3
105  shortLadder = theLadder - 5 * i
106  for i in range(0, shortLadder, 2):
107  sector = sector + 1
108  return sector
109  elif layer == 3:
110  sector = 1
111  for i in range(2, theLadder, 3):
112  if (i + 1) % 3 == 0:
113  sector = sector + 1
114  return sector
115  elif layer == 4:
116  sector = (theLadder + 3) // 4
117  return sector
118 
119  def __BuildOnlineBarrelName(self, signedModule, signedLadder, layer): #in Phase1 it is assumed that there are only full modules
120  thePart = self.__GetPartStr(signedModule < 0, signedLadder < 0)
121  theSector = str(self.__GetBarrelSector(layer, signedLadder, signedModule))
122  return "BPix_B" + thePart + "_SEC" + theSector + "_LYR" + str(layer) + "_LDR" + str(abs(signedLadder)) + "F_MOD" + str(abs(signedModule))
123 
124  def __BuildOnlineDiskName(self, signedDisk, signedBlade, panel, ring):
125  thePart = self.__GetPartStr(signedDisk < 0, signedBlade < 0)
126  return "FPix_B" + thePart + "_D" + str(abs(signedDisk)) + "_BLD" + str(abs(signedBlade)) + "_PNL" + str(panel) + "_RNG" + str(ring)
127 
128  def __GroupHistograms(self):
129  currentGroupName = ""
130  groupOfHists = []
132 
133 
134  for obj in self.listOfNumHistograms:
135  objName = obj.GetName()
136  objNameSplit = objName.split("_")
137  objNameCollected = ''.join(objNameSplit[0:-1])
138  if objNameCollected != currentGroupName:
139  if len(groupOfHists):
140  self.groupedHistograms.append(groupOfHists)
141  groupOfHists = []
142 
143  currentGroupName = objNameCollected
144  groupOfHists.append(obj)
145  self.groupedHistograms.append(groupOfHists) #the last group
146 
147  def __CreateDummyStructAsStr(self, dicData):
148  str = "struct MyStruct{"
149  str = str + "Int_t key;"
150  leafStr = "key/I"
151  for k in dicData:
152  str = str + "Float_t " + k + ";"
153  leafStr = leafStr + ":" + k + "/F"
154  str = str + "};"
155  return str, leafStr
156 
157 
158 
159  def __init__(self, inputDQMName, outputDQMName, modDicName):
160  self.inputFileName = inputDQMName
161  self.outputFileName = outputDQMName
162  self.detIDsFileName = modDicName
163 
164  index = self.inputFileName.find('R000')
165  runNumber = self.inputFileName[index+4:index+10]
166 
167  self.dirs = ['DQMData/Run %s/PixelPhase1/Run summary/Phase1_MechanicalView' % (runNumber),
168  'DQMData/Run %s/PixelPhase1/Run summary/Tracks' % (runNumber)]
169  self.dirsAliases = {self.dirs[0]:"", self.dirs[1]: "T"}
170 
171  self.inputFile = TFile(self.inputFileName)
173  self.availableNames = []
174 
175  self.maxLadderToLayer = {6:1, 14:2, 22:3, 32:4}
176  self.maxBladeToRing = {11:1, 17:2}
177 
178  self.internalData = {}
179 
180  if self.inputFile.IsOpen():
181  print("%s opened successfully!" % (self.inputFileName))
182  #Get all neeeded histograms
183  for dir in self.dirs:
184  self.__TraverseDirTree(self.inputFile.Get(dir))
185  print("Histograms to read %d" % (len(self.listOfNumHistograms)))
186 
187  self.detDict = {}
188 
189  with open(self.detIDsFileName, "r") as detIDs: # create dictionary online -> rawid
190  for entry in detIDs:
191  items = entry.replace("\n", " ").split(" ")
192  self.detDict.update({items[1] : int(items[0])})
193  # init internal data structure
194  self.internalData.update({int(items[0]) : {}})
195 
196  self.__GroupHistograms()
197 
198  else:
199  print("Unable to open file %s" % (self.inputFileName))
200 
201  def ReadHistograms(self):
202  if self.inputFile.IsOpen():
203  for group in self.groupedHistograms:
204  # name = ''.join(group[0].GetName().split("_")[0:-1])
205  name = ''.join(group[0].GetName().split("_per_")[0])
206  self.availableNames.append(name)
207  # print(name)
208  for obj in group:
209  nbinsX = obj.GetNbinsX()
210  nbinsY = obj.GetNbinsY()
211 
212  if nbinsX == 9: # BARREL
213  maxX = nbinsX // 2
214  maxY = nbinsY // 2
215 
216  for x in range(-maxX, maxX + 1):
217  if x == 0:
218  continue
219  for y in range(-maxY, maxY + 1, 1):
220  if y == 0:
221  continue
222  onlineName = self.__BuildOnlineBarrelName(x, y, self.maxLadderToLayer[maxY])
223  self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, y + maxY + 1)})
224 
225  elif nbinsX == 7: # FORWARD
226  maxX = nbinsX // 2
227  maxY = nbinsY // 4
228 
229  for x in range(-maxX, maxX + 1):
230  if x == 0:
231  continue
232  for y in range(-maxY, maxY + 1):
233  if int(y) == 0:
234  continue
235  for panel in range(1, 3):
236  onlineName = self.__BuildOnlineDiskName(x, y, panel, self.maxBladeToRing[maxY])
237  self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, (y + maxY) * 2 + panel)})
238  else:
239  print("Unrecognized plot")
240  else:
241  print("Histograms saved to internal data structure")
242  def CreateTree(self): # too much complication lvl, use CreateTree2
243  if len(self.internalData):
244  # <- TOTAL WORKAROUND: START -> #
245  # SEE: http://wlav.web.cern.ch/wlav/pyroot/tpytree.html
246 
247  leafStr = ""
248  for key in self.internalData:
249  s, leafStr = self.__CreateDummyStructAsStr(self.internalData[key])
250  print(s)
251  print(leafStr)
252  gROOT.ProcessLine(s)
253  break #all modules are assumed to have the same set of measured parameters
254 
255  from ROOT import MyStruct
256  ms = MyStruct()
257  # <- TOTAL WORKAROUND: END -> #
258 
259  self.outputFile = TFile(self.outputFileName, "recreate")
260  tree = TTree("tree", "readData")
261  tree.Branch("b", ms, leafStr)
262 
263  tree.SetBranchAddress("b", tree)
264 
265  for key in self.internalData:
266  setattr(ms, "key", key)
267  for d in self.internalData[key]:
268  setattr(ms, d, (self.internalData[key])[d])
269  tree.Fill()
270 
271  # break
272  tree.Write()
273  self.outputFile.Close()
274 
275  print("Data saved as TTree object")
276 
277  def CreateTree2(self): # use for TkCommissioner
278  if len(self.internalData):
279 
280  self.outputFile = TFile(self.outputFileName, "recreate")
281  tree = TTree("tree", "readData")
282 
283  key = array("i", [0])
284  tree.Branch("detid", key, "detid/i")
285 
286  dat = {}
287  for k in self.internalData:
288  for d in self.internalData[k]:
289  dat.update({d : array("f", [0])})
290  tree.Branch(d, dat[d], d + "/F")
291  break
292 
293  for k in self.internalData:
294  key[0] = k
295  for d in self.internalData[k]:
296  (dat[d])[0] = (self.internalData[k])[d]
297  tree.Fill()
298 
299  tree.Write()
300  self.outputFile.Close()
301 
302  print("Data saved as TTree object")
303 
304  def DumpData(self):
305  for key in self.internalData:
306  print("#"*20)
307  print(key)
308  module = self.internalData[key]
309  for d in module:
310  print((d, module[d]))
311 
312  for i in self.availableNames:
313  print(i)
314  print(len(self.availableNames))
315 
316  def __del__(self):
317  if self.inputFile:
318  if self.inputFile.IsOpen():
319  self.inputFile.Close()
320 
321 #--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--
322 for i in range(1, len(sys.argv), 1):
323  if i == 1:
324  inputFileName = sys.argv[i]
325  elif i == 2:
326 # detIDsFileName = sys.argv[i]
327 # elif i == 3:
328  outputFileName = sys.argv[i]
329 
330 #readerObj = ModuleLvlValuesReader(inputFileName, outputFileName)
331 readerObj = ModuleLvlValuesReader(inputFileName, outputFileName, detIDsFileName)
332 readerObj.ReadHistograms()
333 readerObj.CreateTree2()
334 # readerObj.DumpData()
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
def __init__(self, inputDQMName, outputDQMName, modDicName)
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def __GetBarrelSector(self, layer, signedLadder, signedModule)
static std::string join(char **cmd)
Definition: RemoteFile.cc:19
def __GetPartStr(self, isXlowerThanZero, isYlowerThanZero)
#define update(a, b)
def __BuildOnlineDiskName(self, signedDisk, signedBlade, panel, ring)
#define str(s)
def __BuildOnlineBarrelName(self, signedModule, signedLadder, layer)