File indexing completed on 2024-04-06 12:08:43
0001
0002
0003 from __future__ import print_function
0004 import sys
0005 import os
0006 from ROOT import gROOT, gStyle, gPad, TCanvas, TClass, TGraph, TFile, TArrow, TLatex, TH2Poly, kBlack
0007 from copy import deepcopy
0008 from array import array
0009
0010 gROOT.SetBatch()
0011
0012
0013
0014 def getFileInPath(rfile):
0015 import os
0016 for dir in os.environ['CMSSW_SEARCH_PATH'].split(":"):
0017 if os.path.exists(os.path.join(dir,rfile)): return os.path.join(dir,rfile)
0018 return None
0019
0020
0021
0022 inputFileName = "DQM_V0013_R000292154__StreamExpressCosmics__Commissioning2017-Express-v1__DQMIO.root"
0023 limitsFileName = "limits.dat"
0024 outputDirectoryName = "OUT/"
0025 minMaxFileName = "minmax.out"
0026
0027 detIDsFileName = getFileInPath('DQM/SiStripMonitorClient/data/detids.dat')
0028
0029 baseRootDirs = ["DQMData/Run 292154/PixelPhase1/Run summary/Phase1_MechanicalView"
0030 ,"DQMData/Run 292154/PixelPhase1/Run summary/Tracks"
0031 ]
0032
0033
0034 maxPxBarrel = 4
0035 maxPxForward = 3
0036 barrelLadderShift = [0, 14, 44, 90]
0037
0038 forwardDiskXShift = [25, 75, 125]
0039 forwardDiskYShift = 45;
0040
0041 plotWidth, plotHeight = 3000, 2000
0042 extremeBinsNum = 20
0043
0044 limits = ["num_digis 0.01 90 1 0",
0045 "num_clusters 0.01 25 1 0",
0046 "Trechitsize_y 0.01 10 0 0",
0047 "Trechitsize_x 0.01 10 0 0",
0048 "Tresidual_y 0.0000001 0.004 0 1",
0049 "Tresidual_x 0.0000001 0.004 0 1",
0050 "Tcharge 2000 80000 0 0",
0051 "Thitefficiency 0.95 1 0 0",
0052
0053 "Tnum_clusters_ontrack 0.01 15 1 0",
0054 "Tsize 0.01 15 0 0",
0055
0056 "adc 0.01 256 0 0",
0057 "charge 2000 80000 0 0",
0058 "size 0.01 15 0 0",]
0059
0060 class TH2PolyOfflineMaps:
0061
0062
0063
0064
0065
0066
0067
0068 def __TraverseDirTree(self, dir):
0069
0070 try:
0071 currPath = (dir.GetPath().split(":/"))[1]
0072 except:
0073 print("Exception raised: Path not found in the input file")
0074 return
0075
0076 for obj in dir.GetListOfKeys():
0077 if not obj.IsFolder():
0078 if obj.ReadObjectAny(TClass.GetClass("TH2")):
0079 th2 = deepcopy(obj.ReadObj())
0080 name = th2.GetName()
0081 if 6 < th2.GetNbinsX() < 10 and name.find("per") != -1 and name.find("Lumisection") == -1:
0082 print(''.join([dir.GetPath(), '/', name]))
0083
0084
0085 prefix = ""
0086 for i in self.dirs:
0087 if currPath.startswith(i):
0088 prefix = self.dirsAliases[i]
0089 break
0090
0091 th2.SetName(prefix + th2.GetName())
0092 self.listOfNumHistograms.append(th2)
0093 else:
0094 self.__TraverseDirTree(obj.ReadObj())
0095
0096 def __GetPartStr(self, isXlowerThanZero, isYlowerThanZero):
0097 if isXlowerThanZero and isYlowerThanZero:
0098 return "mO"
0099 if isXlowerThanZero and isYlowerThanZero == False:
0100 return "mI"
0101 if isXlowerThanZero == False and isYlowerThanZero:
0102 return "pO"
0103 if isXlowerThanZero == False and isYlowerThanZero == False:
0104 return "pI"
0105
0106 def __GetBarrelSector(self, layer, signedLadder, signedModule):
0107 theLadder = abs(signedLadder)
0108 theModule = abs(signedModule)
0109
0110 sector = 0
0111
0112 if layer == 1:
0113
0114 if theLadder == 1:
0115 if theModule >= 2:
0116 return 1
0117 else:
0118 return 2
0119 if theLadder == 2:
0120 if theModule >= 3:
0121 return 2
0122 else:
0123 return 3
0124 if theLadder == 3:
0125 if theModule >= 4:
0126 return 3
0127 else:
0128 return 4
0129 if theLadder == 4:
0130 if theModule >= 2:
0131 return 5
0132 else:
0133 return 6
0134 if theLadder == 5:
0135 if theModule >= 3:
0136 return 6
0137 else:
0138 return 7
0139 if theLadder == 6:
0140 if theModule >= 4:
0141 return 7
0142 else:
0143 return 8
0144
0145 elif layer == 2:
0146 i = theLadder // 5
0147 sector = i * 3
0148 shortLadder = theLadder - 5 * i
0149 for i in range(0, shortLadder, 2):
0150 sector = sector + 1
0151 return sector
0152 elif layer == 3:
0153 sector = 1
0154 for i in range(2, theLadder, 3):
0155 if (i + 1) % 3 == 0:
0156 sector = sector + 1
0157 return sector
0158 elif layer == 4:
0159 sector = (theLadder + 3) // 4
0160 return sector
0161
0162 def __BuildOnlineBarrelName(self, signedModule, signedLadder, layer):
0163 thePart = self.__GetPartStr(signedModule < 0, signedLadder < 0)
0164 theSector = str(self.__GetBarrelSector(layer, signedLadder, signedModule))
0165 return "BPix_B" + thePart + "_SEC" + theSector + "_LYR" + str(layer) + "_LDR" + str(abs(signedLadder)) + "F_MOD" + str(abs(signedModule))
0166
0167 def __BuildOnlineDiskName(self, signedDisk, signedBlade, panel, ring):
0168 thePart = self.__GetPartStr(signedDisk < 0, signedBlade < 0)
0169 return "FPix_B" + thePart + "_D" + str(abs(signedDisk)) + "_BLD" + str(abs(signedBlade)) + "_PNL" + str(panel) + "_RNG" + str(ring)
0170
0171 def __GroupHistograms(self):
0172 currentGroupName = ""
0173 groupOfHists = []
0174 self.groupedHistograms = []
0175
0176
0177 for obj in self.listOfNumHistograms:
0178 objName = obj.GetName()
0179 objNameSplit = objName.split("_")
0180 objNameCollected = ''.join(objNameSplit[0:-1])
0181 if objNameCollected != currentGroupName:
0182 if len(groupOfHists):
0183 self.groupedHistograms.append(groupOfHists)
0184 groupOfHists = []
0185
0186 currentGroupName = objNameCollected
0187 groupOfHists.append(obj)
0188 self.groupedHistograms.append(groupOfHists)
0189
0190 def __AddNamedBins(self, geoFile, tX, tY, sX, sY, applyModuleRotation = False):
0191
0192 for line in geoFile:
0193 lineSpl = line.strip().split("\"")
0194
0195 detId = str(lineSpl[0].split(" ")[1])+"_("+str(lineSpl[0].split(" ")[0])+")"
0196 vertices = lineSpl[1]
0197 xy = vertices.split(" ")
0198 x, y = array('d'), array('d')
0199 verNum = 1
0200 for coord in xy:
0201 coordSpl = coord.split(",")
0202 if applyModuleRotation:
0203 x.append(-(float(coordSpl[0]) * sX + tX))
0204 y.append((float(coordSpl[1]) * sY + tY))
0205 else:
0206 x.append(float(coordSpl[0]) * sX + tX)
0207 y.append(float(coordSpl[1]) * sY + tY)
0208 verNum = verNum + 1
0209
0210 x.append(x[0])
0211 y.append(y[0])
0212
0213 if applyModuleRotation:
0214 bin = TGraph(verNum, y, x)
0215 else:
0216 bin = TGraph(verNum, x, y)
0217
0218 bin.SetName(detId)
0219
0220 self.__BaseTrackerMap.AddBin(bin)
0221
0222 def __CreateTrackerBaseMap(self):
0223
0224 self.__BaseTrackerMap = TH2Poly("Summary", "", -10, 160, -70, 70)
0225
0226 self.__BaseTrackerMap.SetFloat(1)
0227 self.__BaseTrackerMap.GetXaxis().SetTitle("")
0228 self.__BaseTrackerMap.GetYaxis().SetTitle("")
0229 self.__BaseTrackerMap.SetOption("COLZ L")
0230 self.__BaseTrackerMap.SetStats(0)
0231
0232
0233 for i in range(maxPxBarrel):
0234 with open(self.geometryFilenames[i], "r") as geoFile:
0235 currBarrelTranslateX = 0
0236 currBarrelTranslateY = barrelLadderShift[i]
0237
0238 self.__AddNamedBins(geoFile, currBarrelTranslateX, currBarrelTranslateY, 1, 1, True)
0239
0240
0241
0242
0243 for i in range(-maxPxForward, 0):
0244 with open(self.geometryFilenames[maxPxBarrel + maxPxForward + i], "r") as geoFile:
0245 currForwardTranslateX = forwardDiskXShift[-i - 1]
0246 currForwardTranslateY = -forwardDiskYShift
0247
0248 self.__AddNamedBins(geoFile, currForwardTranslateX, currForwardTranslateY, 1, 1)
0249
0250
0251 for i in range(maxPxForward):
0252 with open(self.geometryFilenames[maxPxBarrel + maxPxForward + i], "r") as geoFile:
0253 currForwardTranslateX = forwardDiskXShift[i]
0254 currForwardTranslateY = forwardDiskYShift
0255
0256 self.__AddNamedBins(geoFile, currForwardTranslateX, currForwardTranslateY, 1, 1)
0257
0258
0259
0260 print("Base Tracker Map: constructed")
0261
0262
0263 def __init__(self, inputDQMName, outputDirName, minMaxFileName, limits, modDicName, runNumber, dirs, dirsAliases):
0264
0265 self.inputFileName = inputDQMName
0266 self.outputDirName = outputDirName
0267 self.minMaxFileName = minMaxFileName
0268
0269 self.detIDsFileName = modDicName
0270 self.limits = limits
0271
0272 self.runNumber = runNumber
0273 self.dirs = dirs
0274 self.dirsAliases = dirsAliases
0275
0276 self.inputFile = TFile(self.inputFileName)
0277 self.listOfNumHistograms = []
0278 self.availableNames = []
0279
0280 self.maxLadderToLayer = {6:1, 14:2, 22:3, 32:4}
0281 self.maxBladeToRing = {11:1, 17:2}
0282
0283 self.geometryFilenames = []
0284 for i in range(maxPxBarrel):
0285 self.geometryFilenames.append(getFileInPath("DQM/SiStripMonitorClient/data/Geometry/vertices_barrel_" + str(i + 1)))
0286
0287 for i in range(-maxPxForward, maxPxForward + 1):
0288 if i == 0:
0289 continue
0290 self.geometryFilenames.append(getFileInPath("DQM/SiStripMonitorClient/data/Geometry/vertices_forward_" + str(i)))
0291
0292
0293 self.internalData = {}
0294
0295 if self.inputFile.IsOpen():
0296 print("%s opened successfully!" % (self.inputFileName))
0297
0298 for dir in self.dirs:
0299 self.__TraverseDirTree(self.inputFile.Get(dir))
0300
0301
0302 self.detDict = {}
0303
0304 with open(self.detIDsFileName, "r") as detIDs:
0305 for entry in detIDs:
0306 items = entry.replace("\n", " ").split(" ")
0307 self.detDict.update({items[1] : int(items[0])})
0308
0309 self.internalData.update({int(items[0]) : {}})
0310
0311 self.rawToOnlineDict = dict((v,k) for k,v in self.detDict.items())
0312
0313 self.__GroupHistograms()
0314
0315 self.__CreateTrackerBaseMap()
0316
0317 else:
0318 print("Unable to open file %s" % (self.inputFileName))
0319
0320
0321
0322 self.limitsDic = {}
0323 for y in limits:
0324
0325 lineSpl = y.strip().split(" ")
0326
0327 if len(lineSpl) < 5:
0328 continue
0329
0330 currName = lineSpl[0]
0331 zMin = float(lineSpl[1])
0332 zMax = float(lineSpl[2])
0333 isLog = False if lineSpl[3] == "0" else True
0334 isAbs = False if lineSpl[4] == "0" else True
0335
0336 self.limitsDic.update({currName : {"zMin" : zMin, "zMax" : zMax, "isLog" : isLog, "isAbs" : isAbs}})
0337
0338
0339 def ReadHistograms(self):
0340 if self.inputFile.IsOpen():
0341 for group in self.groupedHistograms:
0342
0343 if len(group) == 0:
0344 return
0345 print(group[0].GetName())
0346 name = ''.join(group[0].GetName().split("_per_")[0])
0347 self.availableNames.append(name)
0348
0349 for obj in group:
0350 nbinsX = obj.GetNbinsX()
0351 nbinsY = obj.GetNbinsY()
0352
0353 if nbinsX == 9:
0354 maxX = nbinsX // 2
0355 maxY = nbinsY // 2
0356
0357 for x in range(-maxX, maxX + 1):
0358 if x == 0:
0359 continue
0360 for y in range(-maxY, maxY + 1, 1):
0361 if y == 0:
0362 continue
0363 onlineName = self.__BuildOnlineBarrelName(x, y, self.maxLadderToLayer[maxY])
0364 self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, y + maxY + 1)})
0365
0366 elif nbinsX == 7:
0367 maxX = nbinsX // 2
0368 maxY = nbinsY // 4
0369
0370 for x in range(-maxX, maxX + 1):
0371 if x == 0:
0372 continue
0373 for y in range(-maxY, maxY + 1):
0374 if int(y) == 0:
0375 continue
0376 for panel in range(1, 3):
0377 onlineName = self.__BuildOnlineDiskName(x, y, panel, self.maxBladeToRing[maxY])
0378 self.internalData[self.detDict[onlineName]].update({name : obj.GetBinContent(x + maxX + 1, (y + maxY) * 2 + (3-panel))})
0379 else:
0380 print("Unrecognized plot")
0381 else:
0382 print("Histograms saved to internal data structure")
0383
0384 def DumpData(self):
0385 for key in self.internalData:
0386 print("#"*20)
0387 print(key)
0388 module = self.internalData[key]
0389 for d in module:
0390 print((d, module[d]))
0391
0392 print(len(self.internalData))
0393
0394 for i in self.availableNames:
0395 print(i)
0396 print(len(self.availableNames))
0397
0398 def PrintTrackerMaps(self):
0399 monitoredValues = []
0400 gStyle.SetPalette(1)
0401 for key in self.internalData:
0402 monitoredValues = self.internalData[key].keys()
0403
0404 break
0405
0406 if os.path.exists(self.outputDirName) == False:
0407 os.system("mkdir " + self.outputDirName)
0408
0409 with open(self.outputDirName + self.minMaxFileName, "w") as minMaxFile:
0410
0411 for mv in monitoredValues:
0412 currentHist = deepcopy(self.__BaseTrackerMap)
0413
0414 histoTitle = "Run " + self.runNumber + ": Tracker Map for " + mv
0415
0416 applyLogScale = False
0417 applyAbsValue = False
0418 if mv in self.limitsDic:
0419 limitsElem = self.limitsDic[mv]
0420
0421 print(mv + " found in limits dictionary - applying custom limits...")
0422
0423 currentHist.SetMinimum(limitsElem["zMin"])
0424 currentHist.SetMaximum(limitsElem["zMax"])
0425 applyLogScale = limitsElem["isLog"]
0426 applyAbsValue = limitsElem["isAbs"]
0427
0428 listOfVals = []
0429 onlineName = ""
0430 nameId = ""
0431
0432 for detId in self.internalData:
0433 val = (self.internalData[detId])[mv]
0434 onlineName = self.rawToOnlineDict[detId]
0435 listOfVals.append([val, detId, onlineName])
0436
0437 nameId = str(onlineName)+"_("+str(detId)+")"
0438
0439 if applyAbsValue:
0440 currentHist.Fill(str(nameId), abs(val))
0441 else:
0442 currentHist.Fill(str(nameId), val)
0443
0444 listOfVals = sorted(listOfVals, key = lambda item: item[0])
0445
0446 minMaxFile.write("\n" + mv + "\n\n")
0447
0448 minMaxFile.write("MIN:\n")
0449 for i in range(extremeBinsNum):
0450 minMaxFile.write("\t" + str(listOfVals[i][1]) + " " + str(listOfVals[i][2]) + " " + str(listOfVals[i][0]) + "\n")
0451
0452 minMaxFile.write("MAX:\n")
0453 for i in range(extremeBinsNum):
0454 minMaxFile.write("\t" + str(listOfVals[-i - 1][1]) + " " + str(listOfVals[-i - 1][2]) + " " + str(listOfVals[-i - 1][0]) + "\n")
0455
0456
0457 c1 = TCanvas("MyT", "MyT", plotWidth , plotHeight)
0458
0459 if applyLogScale:
0460 c1.SetLogz()
0461
0462 currentHist.Draw("AC COLZ L")
0463
0464 gPad.Update()
0465 palette = currentHist.FindObject("palette");
0466 palette.SetX1NDC(0.89);
0467 palette.SetX2NDC(0.91);
0468 palette.SetLabelSize(0.05);
0469 gPad.Update()
0470
0471
0472
0473
0474
0475
0476
0477
0478
0479 arrow = TArrow(0.05, 27.0, 0.05, -30.0, 0.02, "|>")
0480 arrow.SetLineWidth(4)
0481 arrow.Draw()
0482
0483 phiArrow = TArrow(0.0, 27.0, 30.0, 27.0, 0.02, "|>")
0484 phiArrow.SetLineWidth(4)
0485 phiArrow.Draw()
0486
0487 xArrow = TArrow(25.0, 44.5, 50.0, 44.5, 0.02, "|>")
0488 xArrow.SetLineWidth(4)
0489 xArrow.Draw()
0490
0491 yArrow = TArrow(25.0, 44.5, 25.0, 69.5, 0.02, "|>")
0492 yArrow.SetLineWidth(4)
0493 yArrow.Draw()
0494
0495
0496
0497
0498 txt = TLatex()
0499 txt.SetNDC()
0500 txt.SetTextFont(1)
0501 txt.SetTextColor(1)
0502 txt.SetTextAlign(22)
0503 txt.SetTextAngle(0)
0504
0505
0506 txt.SetTextSize(0.05)
0507 txt.DrawLatex(0.5, 0.95, histoTitle)
0508
0509 txt.SetTextSize(0.03)
0510
0511 txt.DrawLatex(0.5, 0.125, "-DISK")
0512 txt.DrawLatex(0.5, 0.075, "NUMBER ->")
0513 txt.DrawLatex(0.5, 0.875, "+DISK")
0514
0515 txt.DrawLatex(0.17, 0.35, "+z")
0516 txt.DrawLatexNDC(0.36, 0.685, "+phi")
0517 txt.DrawLatex(0.38, 0.73, "+x")
0518 txt.DrawLatex(0.26, 0.875, "+y")
0519
0520 txt.SetTextAngle(90)
0521 txt.DrawLatex(0.17, 0.5, "BARREL")
0522
0523
0524 c1.Print(self.outputDirName + mv + ".png")
0525
0526
0527 c1.Clear()
0528 c1.SetLogz(False)
0529 currentHist.GetZaxis().UnZoom()
0530 currentHist.SetLineColor(kBlack)
0531 currentHist.Draw("AL COLZ")
0532 currentHist.GetXaxis().SetRangeUser(-10,155)
0533
0534 palette.SetX1NDC(0.92);
0535 palette.SetX2NDC(0.94);
0536 palette.SetY1NDC(0.02);
0537 palette.SetY2NDC(0.91);
0538 gPad.SetRightMargin(0.08);
0539 gPad.SetLeftMargin(0.01);
0540 gPad.SetTopMargin(0.09);
0541 gPad.SetBottomMargin(0.02);
0542 gPad.Update()
0543
0544 zarrow = TArrow(0, 27, 0, -30, 0.02, "|>")
0545 zarrow.SetLineWidth(3)
0546 zarrow.Draw()
0547 phiArrow.SetLineWidth(3)
0548 phiArrow.Draw()
0549 xArrow.SetLineWidth(3)
0550 xArrow.Draw()
0551 yArrow.SetLineWidth(3)
0552 yArrow.Draw()
0553
0554 txt.Clear()
0555 txt.SetTextAngle(0)
0556 txt.SetTextSize(0.05)
0557 PixelTitle = "Run " + self.runNumber + ": Pixel " + mv
0558 txt.DrawLatex(0.5, 0.95, PixelTitle)
0559
0560 txt.SetTextSize(0.04)
0561 txt.SetNDC(False)
0562 txt.DrawLatex(75, -65, "-DISK")
0563 txt.DrawLatex(75, 65, "+DISK")
0564 txt.DrawLatex(50, -60, "NUMBER ->")
0565
0566 txt.DrawLatex(-5, -20, "+z")
0567 txt.DrawLatex(35, 30, "+phi")
0568 txt.DrawLatex(55, 45, "+x")
0569 txt.DrawLatex(30, 65, "+y")
0570
0571 txt.SetTextAngle(90)
0572 txt.DrawLatex(-5, 0, "BARREL")
0573
0574 c1.SaveAs(self.outputDirName + mv + ".root")
0575 c1.Close()
0576
0577 def __del__(self):
0578 if self.inputFile :
0579 if self.inputFile.IsOpen():
0580 self.inputFile.Close()
0581
0582
0583 for i in range(1, len(sys.argv), 1):
0584 if i == 1:
0585 inputFileName = sys.argv[i]
0586 elif i == 2:
0587 plotWidth = int(sys.argv[i])
0588 elif i == 3:
0589 plotHeight = int(sys.argv[i])
0590
0591
0592
0593 elif i == 4:
0594 detIDsFileName = sys.argv[i]
0595
0596 deductedRunNumber = inputFileName.split("_R000")[1][0:6]
0597 print(deductedRunNumber)
0598
0599 baseRootDirs = ["DQMData/Run " + deductedRunNumber + "/PixelPhase1/Run summary/Phase1_MechanicalView"
0600 ,"DQMData/Run " + deductedRunNumber + "/PixelPhase1/Run summary/Tracks"
0601 ]
0602
0603 baseRootDirsAliases = {baseRootDirs[0]:""
0604 , baseRootDirs[1]:"T"
0605 }
0606
0607 readerObj = TH2PolyOfflineMaps(inputFileName, outputDirectoryName, minMaxFileName, limits, detIDsFileName, deductedRunNumber, baseRootDirs, baseRootDirsAliases)
0608
0609 readerObj.ReadHistograms()
0610
0611 readerObj.PrintTrackerMaps()