CMS 3D CMS Logo

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