CMS 3D CMS Logo

TH2PolyOfflineMaps.py
Go to the documentation of this file.
1 #!/usr/bin/env python
2 
3 import sys
4 import os
5 from ROOT import *
6 from copy import deepcopy
7 from array import array
8 
9 gROOT.SetBatch() # don't pop up canvases
10 
11 #Find 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 
20 # Default values
21 inputFileName = "DQM_V0013_R000292154__StreamExpressCosmics__Commissioning2017-Express-v1__DQMIO.root"
22 limitsFileName = "limits.dat"
23 outputDirectoryName = "OUT/"
24 minMaxFileName = "minmax.out"
25 #detIDsFileName = "DATA/detids.dat"
26 detIDsFileName = getFileInPath('DQM/SiStripMonitorClient/data/detids.dat')
27 #default one
28 baseRootDirs = ["DQMData/Run 292154/PixelPhase1/Run summary/Phase1_MechanicalView"
29  ,"DQMData/Run 292154/PixelPhase1/Run summary/Tracks"
30  ]
31 
32 
33 maxPxBarrel = 4
34 maxPxForward = 3
35 barrelLadderShift = [0, 14, 44, 90]
36 
37 forwardDiskXShift = [25, 75, 125]
38 forwardDiskYShift = 45; # to make +DISK on top in the 'strip-like' layout
39 
40 plotWidth, plotHeight = 3000, 2000
41 extremeBinsNum = 20
42 
43 limits = ["num_digis 0.01 90 1 0",
44  "num_clusters 0.01 25 1 0",
45  "Trechitsize_y 0.01 10 0 0",
46  "Trechitsize_x 0.01 10 0 0",
47  "Tresidual_y 0.0000001 0.004 0 1",
48  "Tresidual_x 0.0000001 0.004 0 1",
49  "Tcharge 2000 80000 0 0",
50  "Thitefficiency 0.95 1 0 0",
51  #"Tmissing 0.01 500 0 0",
52  "Tnum_clusters_ontrack 0.01 15 1 0",
53  "Tsize 0.01 15 0 0",
54  #"Tvalid 0.01 90 0 0",
55  "adc 0.01 256 0 0",
56  "charge 2000 80000 0 0",
57  "size 0.01 15 0 0",]
58 
60 
61  ###
62  # LOTS OF CODE BORROWED FROM: PYTHONBINREADER, PIXELTRACKERMAP
63  ###
64 
65  ############################################################################
66 
67  def __TraverseDirTree(self, dir):
68 
69  try:
70  currPath = (dir.GetPath().split(":/"))[1]
71  except:
72  print("Exception raised: Path not found in the input file")
73  return
74 
75  for obj in dir.GetListOfKeys():
76  if not obj.IsFolder():
77  if obj.ReadObjectAny(TClass.GetClass("TH2")):
78  th2 = deepcopy(obj.ReadObj())
79  name = th2.GetName()
80  if 6 < th2.GetNbinsX() < 10 and name.find("per") != -1 and name.find("Lumisection") == -1: #take only module lvl plots
81  print(''.join([dir.GetPath(), '/', name]))
82 
83  # fix when there are plots starting with the same strings in different directories
84  prefix = ""
85  for i in self.dirs:
86  if currPath.startswith(i):
87  prefix = self.dirsAliases[i]
88  break
89  # print(currPath, prefix)
90  th2.SetName(prefix + th2.GetName())
91  self.listOfNumHistograms.append(th2)
92  else:
93  self.__TraverseDirTree(obj.ReadObj())
94 
95  def __GetPartStr(self, isXlowerThanZero, isYlowerThanZero):
96  if isXlowerThanZero and isYlowerThanZero:
97  return "mO"
98  if isXlowerThanZero and isYlowerThanZero == False:
99  return "mI"
100  if isXlowerThanZero == False and isYlowerThanZero:
101  return "pO"
102  if isXlowerThanZero == False and isYlowerThanZero == False:
103  return "pI"
104 
105  def __GetBarrelSector(self, layer, signedLadder, signedModule): #adapted from PixelBarrelName
106  theLadder = abs(signedLadder)
107  theModule = abs(signedModule)
108 
109  sector = 0
110 
111  if layer == 1:
112 
113  if theLadder == 1:
114  if theModule >= 2:
115  return 1
116  else:
117  return 2
118  if theLadder == 2:
119  if theModule >= 3:
120  return 2
121  else:
122  return 3
123  if theLadder == 3:
124  if theModule >= 4:
125  return 3
126  else:
127  return 4
128  if theLadder == 4:
129  if theModule >= 2:
130  return 5
131  else:
132  return 6
133  if theLadder == 5:
134  if theModule >= 3:
135  return 6
136  else:
137  return 7
138  if theLadder == 6:
139  if theModule >= 4:
140  return 7
141  else:
142  return 8
143  # here is used simplified form of assignment, see source file for reference
144  elif layer == 2:
145  i = theLadder // 5
146  sector = i * 3
147  shortLadder = theLadder - 5 * i
148  for i in range(0, shortLadder, 2):
149  sector = sector + 1
150  return sector
151  elif layer == 3:
152  sector = 1
153  for i in range(2, theLadder, 3):
154  if (i + 1) % 3 == 0:
155  sector = sector + 1
156  return sector
157  elif layer == 4:
158  sector = (theLadder + 3) // 4
159  return sector
160 
161  def __BuildOnlineBarrelName(self, signedModule, signedLadder, layer): #in Phase1 it is assumed that there are only full modules
162  thePart = self.__GetPartStr(signedModule < 0, signedLadder < 0)
163  theSector = str(self.__GetBarrelSector(layer, signedLadder, signedModule))
164  return "BPix_B" + thePart + "_SEC" + theSector + "_LYR" + str(layer) + "_LDR" + str(abs(signedLadder)) + "F_MOD" + str(abs(signedModule))
165 
166  def __BuildOnlineDiskName(self, signedDisk, signedBlade, panel, ring):
167  thePart = self.__GetPartStr(signedDisk < 0, signedBlade < 0)
168  return "FPix_B" + thePart + "_D" + str(abs(signedDisk)) + "_BLD" + str(abs(signedBlade)) + "_PNL" + str(panel) + "_RNG" + str(ring)
169 
170  def __GroupHistograms(self):
171  currentGroupName = ""
172  groupOfHists = []
174 
175  ##### GROUP ALL LAYERS/RINGS HAVING SIMILAR INFORMATION
176  for obj in self.listOfNumHistograms:
177  objName = obj.GetName()
178  objNameSplit = objName.split("_")
179  objNameCollected = ''.join(objNameSplit[0:-1])
180  if objNameCollected != currentGroupName:
181  if len(groupOfHists):
182  self.groupedHistograms.append(groupOfHists)
183  groupOfHists = []
184 
185  currentGroupName = objNameCollected
186  groupOfHists.append(obj)
187  self.groupedHistograms.append(groupOfHists) #the last group
188 
189  def __AddNamedBins(self, geoFile, tX, tY, sX, sY, applyModuleRotation = False):
190 
191  for line in geoFile:
192  lineSpl = line.strip().split("\"")
193 
194  detId = lineSpl[0].split(" ")[0]
195 
196  vertices = lineSpl[1]
197  xy = vertices.split(" ")
198  x, y = array('d'), array('d')
199  verNum = 1
200  for coord in xy:
201  coordSpl = coord.split(",")
202  if applyModuleRotation:
203  x.append(-(float(coordSpl[0]) * sX + tX))
204  y.append((float(coordSpl[1]) * sY + tY))
205  else:
206  x.append(float(coordSpl[0]) * sX + tX)
207  y.append(float(coordSpl[1]) * sY + tY)
208  verNum = verNum + 1
209  #close polygon
210  x.append(x[0])
211  y.append(y[0])
212 
213  # print(detId, vertices)
214  # print(x)
215  # print(y)
216  if applyModuleRotation:
217  bin = TGraph(verNum, y, x)
218  else:
219  bin = TGraph(verNum, x, y)
220  # bin = TGraph(verNum, y, x) # rotation by 90 deg (so that it had the same layout as for the strips)
221  bin.SetName(detId)
222 
223  self.__BaseTrackerMap.AddBin(bin)
224 
226 
227  self.__BaseTrackerMap = TH2Poly("Summary", "", -10, 160, -70, 70)
228  # self.__BaseTrackerMap = TH2Poly("Summary", "Tracker Map", 0, 0, 0, 0)
229  self.__BaseTrackerMap.SetFloat(1)
230  self.__BaseTrackerMap.GetXaxis().SetTitle("")
231  self.__BaseTrackerMap.GetYaxis().SetTitle("")
232  self.__BaseTrackerMap.SetOption("COLZ L")
233  self.__BaseTrackerMap.SetStats(0)
234 
235  # BARREL FIRST
236  for i in range(maxPxBarrel):
237  with open(self.geometryFilenames[i], "r") as geoFile:
238  currBarrelTranslateX = 0
239  currBarrelTranslateY = barrelLadderShift[i]
240 
241  self.__AddNamedBins(geoFile, currBarrelTranslateX, currBarrelTranslateY, 1, 1, True)
242 
243  # break # debug only 1st layer
244 
245  # MINUS FORWARD
246  for i in range(-maxPxForward, 0):
247  with open(self.geometryFilenames[maxPxBarrel + maxPxForward + i], "r") as geoFile:
248  currForwardTranslateX = forwardDiskXShift[-i - 1]
249  currForwardTranslateY = -forwardDiskYShift
250 
251  self.__AddNamedBins(geoFile, currForwardTranslateX, currForwardTranslateY, 1, 1)
252 
253  # PLUS FORWARD
254  for i in range(maxPxForward):
255  with open(self.geometryFilenames[maxPxBarrel + maxPxForward + i], "r") as geoFile:
256  currForwardTranslateX = forwardDiskXShift[i]
257  currForwardTranslateY = forwardDiskYShift
258 
259  self.__AddNamedBins(geoFile, currForwardTranslateX, currForwardTranslateY, 1, 1)
260 
261  # self.__BaseTrackerMap.Fill("305139728", 2)
262 
263  print("Base Tracker Map: constructed")
264 
265  ############################################################################
266  def __init__(self, inputDQMName, outputDirName, minMaxFileName, limits, modDicName, runNumber, dirs, dirsAliases):
267 # def __init__(self, inputDQMName, outputDirName, minMaxFileName, limitsFileName, modDicName, runNumber, dirs, dirsAliases):
268  self.inputFileName = inputDQMName
269  self.outputDirName = outputDirName
270  self.minMaxFileName = minMaxFileName
271 # self.limitsFileName = limitsFileName
272  self.detIDsFileName = modDicName
273  self.limits = limits
274 
275  self.runNumber = runNumber
276  self.dirs = dirs
277  self.dirsAliases = dirsAliases
278 
279  self.inputFile = TFile(self.inputFileName)
281  self.availableNames = []
282 
283  self.maxLadderToLayer = {6:1, 14:2, 22:3, 32:4}
284  self.maxBladeToRing = {11:1, 17:2}
285 
287  for i in range(maxPxBarrel):
288  self.geometryFilenames.append(getFileInPath("DQM/SiStripMonitorClient/data/Geometry/vertices_barrel_" + str(i + 1)))
289 # self.geometryFilenames.append("DATA/Geometry/vertices_barrel_" + str(i + 1))
290  for i in range(-maxPxForward, maxPxForward + 1):
291  if i == 0:
292  continue #there is no 0 disk
293  self.geometryFilenames.append(getFileInPath("DQM/SiStripMonitorClient/data/Geometry/vertices_forward_" + str(i)))
294 # self.geometryFilenames.append("DATA/Geometry/vertices_forward_" + str(i))
295 
296  self.internalData = {}
297 
298  if self.inputFile.IsOpen():
299  print("%s opened successfully!" % (self.inputFileName))
300  #Get all neeeded histograms
301  for dir in self.dirs:
302  self.__TraverseDirTree(self.inputFile.Get(dir))
303  # print("Histograms to read %d" % (len(self.listOfNumHistograms)))
304 
305  self.detDict = {}
306 
307  with open(self.detIDsFileName, "r") as detIDs: # create dictionary online -> rawid
308  for entry in detIDs:
309  items = entry.replace("\n", " ").split(" ")
310  self.detDict.update({items[1] : int(items[0])})
311  # init internal data structure
312  self.internalData.update({int(items[0]) : {}})
313 
314  self.rawToOnlineDict = dict((v,k) for k,v in self.detDict.iteritems())
315 
316  self.__GroupHistograms()
317 
319 
320  else:
321  print("Unable to open file %s" % (self.inputFileName))
322 
323  ### CREATE LIMITS DICTIONARY
324 
325  self.limitsDic = {}
326  for y in limits:
327 
328  lineSpl = y.strip().split(" ")
329 
330  if len(lineSpl) < 5:
331  continue
332 
333  currName = lineSpl[0]
334  zMin = float(lineSpl[1])
335  zMax = float(lineSpl[2])
336  isLog = False if lineSpl[3] == "0" else True
337  isAbs = False if lineSpl[4] == "0" else True
338 
339  self.limitsDic.update({currName : {"zMin" : zMin, "zMax" : zMax, "isLog" : isLog, "isAbs" : isAbs}})
340  # print limitsDic
341 
342  def ReadHistograms(self):
343  if self.inputFile.IsOpen():
344  for group in self.groupedHistograms:
345  # name = ''.join(group[0].GetName().split("_")[0:-1])
346  if len(group) == 0:
347  return
348  print(group[0].GetName())
349  name = ''.join(group[0].GetName().split("_per_")[0])
350  self.availableNames.append(name)
351  # print(name)
352  for obj in group:
353  nbinsX = obj.GetNbinsX()
354  nbinsY = obj.GetNbinsY()
355 
356  if nbinsX == 9: # BARREL
357  maxX = nbinsX // 2
358  maxY = nbinsY // 2
359 
360  for x in range(-maxX, maxX + 1):
361  if x == 0:
362  continue
363  for y in range(-maxY, maxY + 1, 1):
364  if y == 0:
365  continue
366  onlineName = self.__BuildOnlineBarrelName(x, y, self.maxLadderToLayer[maxY])
367  self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, y + maxY + 1)})
368 
369  elif nbinsX == 7: # FORWARD
370  maxX = nbinsX // 2
371  maxY = nbinsY // 4
372 
373  for x in range(-maxX, maxX + 1):
374  if x == 0:
375  continue
376  for y in range(-maxY, maxY + 1):
377  if int(y) == 0:
378  continue
379  for panel in range(1, 3):
380  onlineName = self.__BuildOnlineDiskName(x, y, panel, self.maxBladeToRing[maxY])
381  self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, (y + maxY) * 2 + (3-panel))})
382  else:
383  print("Unrecognized plot")
384  else:
385  print("Histograms saved to internal data structure")
386 
387  def DumpData(self):
388  for key in self.internalData:
389  print("#"*20)
390  print(key)
391  module = self.internalData[key]
392  for d in module:
393  print(d, module[d])
394 
395  print(len(self.internalData))
396 
397  for i in self.availableNames:
398  print(i)
399  print(len(self.availableNames))
400 
401  def PrintTrackerMaps(self):
402  monitoredValues = []
403  gStyle.SetPalette(1)
404  for key in self.internalData:
405  monitoredValues = self.internalData[key].keys()
406  # print(monitoredValues)
407  break
408 
409  if os.path.exists(self.outputDirName) == False: # check whether directory exists
410  os.system("mkdir " + self.outputDirName)
411 
412  with open(self.outputDirName + self.minMaxFileName, "w") as minMaxFile:
413 
414  for mv in monitoredValues:
415  currentHist = deepcopy(self.__BaseTrackerMap)
416  # currentHist.SetTitle("Run " + self.runNumber + ": Tracker Map for " + mv) // to make it compatible between ROOT v.
417  histoTitle = "Run " + self.runNumber + ": Tracker Map for " + mv
418 
419  applyLogScale = False
420  applyAbsValue = False
421  if mv in self.limitsDic:
422  limitsElem = self.limitsDic[mv]
423 
424  print(mv + " found in limits dictionary - applying custom limits...")
425 
426  currentHist.SetMinimum(limitsElem["zMin"])
427  currentHist.SetMaximum(limitsElem["zMax"])
428  applyLogScale = limitsElem["isLog"]
429  applyAbsValue = limitsElem["isAbs"]
430 
431  listOfVals = []
432  onlineName = ""
433  for detId in self.internalData:
434  val = (self.internalData[detId])[mv]
435  onlineName = self.rawToOnlineDict[detId]
436  listOfVals.append([val, detId, onlineName])
437  if applyAbsValue:
438  currentHist.Fill(str(detId), abs(val))
439  else:
440  currentHist.Fill(str(detId), val)
441 
442  listOfVals = sorted(listOfVals, key = lambda item: item[0])
443 
444  minMaxFile.write("\n" + mv + "\n\n")
445 
446  minMaxFile.write("MIN:\n")
447  for i in range(extremeBinsNum):
448  minMaxFile.write("\t" + str(listOfVals[i][1]) + " " + str(listOfVals[i][2]) + " " + str(listOfVals[i][0]) + "\n")
449 
450  minMaxFile.write("MAX:\n")
451  for i in range(extremeBinsNum):
452  minMaxFile.write("\t" + str(listOfVals[-i - 1][1]) + " " + str(listOfVals[-i - 1][2]) + " " + str(listOfVals[-i - 1][0]) + "\n")
453 
454  c1 = TCanvas(mv, mv, plotWidth , plotHeight)
455 
456  if applyLogScale:
457  c1.SetLogz()
458 
459 
460 
461  currentHist.Draw("AC COLZ L")
462 
463  gPad.Update()
464  palette = currentHist.FindObject("palette");
465  palette.SetX1NDC(0.89);
466  palette.SetX2NDC(0.91);
467  palette.SetLabelSize(0.05);
468  gPad.Update()
469 
470 
471 
472  ### IMPORTANT - REALTIVE POSITIONING IS MESSY IN CURRENT VERION OF PYROOT
473  ### IT CAN CHANGE FROM VERSION TO VERSION, SO YOU HAVE TO ADJUST IT FOR YOUR NEEDS
474  ### !!!!!!!!!!!!!
475 
476  # draw axes (z, phi -> BARREL; x, y -> FORWARD)
477  ###################################################
478 
479  ### z arrow
480  arrow = TArrow(0.05, 27.0, 0.05, -30.0, 0.02, "|>")
481  arrow.SetLineWidth(4)
482  arrow.Draw()
483  ### phi arrow
484  phiArrow = TArrow(0.0, 27.0, 30.0, 27.0, 0.02, "|>")
485  phiArrow.SetLineWidth(4)
486  phiArrow.Draw()
487  ### x arror
488  xArrow = TArrow(25.0, 44.5, 50.0, 44.5, 0.02, "|>")
489  xArrow.SetLineWidth(4)
490  xArrow.Draw()
491  ### y arror
492  yArrow = TArrow(25.0, 44.5, 25.0, 69.5, 0.02, "|>")
493  yArrow.SetLineWidth(4)
494  yArrow.Draw()
495 
496  ###################################################
497 
498  # add some captions
499  txt = TLatex()
500  txt.SetNDC()
501  txt.SetTextFont(1)
502  txt.SetTextColor(1)
503  txt.SetTextAlign(22)
504  txt.SetTextAngle(0)
505 
506  # draw new-style title
507  txt.SetTextSize(0.05)
508  txt.DrawLatex(0.5, 0.95, histoTitle)
509 
510  txt.SetTextSize(0.03)
511 
512  txt.DrawLatex(0.5, 0.125, "-DISK")
513  txt.DrawLatex(0.5, 0.075, "NUMBER ->")
514  txt.DrawLatex(0.5, 0.875, "+DISK")
515 
516  txt.DrawLatex(0.17, 0.35, "+z")
517  txt.DrawLatexNDC(0.36, 0.685, "+phi") # WAY TO FORCE IT TO DRAW LATEX CORRECTLY NOT FOUND ('#' DOESN'T WORK)
518  txt.DrawLatex(0.38, 0.73, "+x")
519  txt.DrawLatex(0.26, 0.875, "+y")
520 
521  txt.SetTextAngle(90)
522  txt.DrawLatex(0.17, 0.5, "BARREL")
523 
524  #save to the png
525  c1.Print(self.outputDirName + mv + ".png")
526 
527  def __del__(self):
528  if self.inputFile.IsOpen():
529  self.inputFile.Close()
530 
531 #--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--+--
532 for i in range(1, len(sys.argv), 1):
533  if i == 1:
534  inputFileName = sys.argv[i]
535  elif i == 2:
536  plotWidth = int(sys.argv[i])
537  elif i == 3:
538  plotHeight = int(sys.argv[i])
539 # elif i == 4:
540 # limitsFileName = sys.argv[i]
541 # elif i == 5:
542  elif i == 4:
543  detIDsFileName = sys.argv[i]
544 
545 deductedRunNumber = inputFileName.split("_R000")[1][0:6]
546 print(deductedRunNumber)
547 
548 baseRootDirs = ["DQMData/Run " + deductedRunNumber + "/PixelPhase1/Run summary/Phase1_MechanicalView" #maybe read it from the input file???
549  ,"DQMData/Run " + deductedRunNumber + "/PixelPhase1/Run summary/Tracks"
550  ]
551 
552 baseRootDirsAliases = {baseRootDirs[0]:""
553  , baseRootDirs[1]:"T"
554  }
555 
556 readerObj = TH2PolyOfflineMaps(inputFileName, outputDirectoryName, minMaxFileName, limits, detIDsFileName, deductedRunNumber, baseRootDirs, baseRootDirsAliases)
557 #readerObj = TH2PolyOfflineMaps(inputFileName, outputDirectoryName, minMaxFileName, limitsFileName, detIDsFileName, deductedRunNumber, baseRootDirs, baseRootDirsAliases)
558 readerObj.ReadHistograms()
559 # readerObj.DumpData()
560 readerObj.PrintTrackerMaps()
def __GetBarrelSector(self, layer, signedLadder, signedModule)
def __TraverseDirTree(self, dir)
LOTS OF CODE BORROWED FROM: PYTHONBINREADER, PIXELTRACKERMAP.
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:65
def __BuildOnlineDiskName(self, signedDisk, signedBlade, panel, ring)
def __init__(self, inputDQMName, outputDirName, minMaxFileName, limits, modDicName, runNumber, dirs, dirsAliases)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
def __GetPartStr(self, isXlowerThanZero, isYlowerThanZero)
def __BuildOnlineBarrelName(self, signedModule, signedLadder, layer)
#define update(a, b)
def __AddNamedBins(self, geoFile, tX, tY, sX, sY, applyModuleRotation=False)
double split
Definition: MVATrainer.cc:139