Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:42

0001 #!/usr/bin/env python3
0002 
0003 from __future__ import print_function
0004 import ROOT
0005 from ROOT import TBufferFile, TH1F, TProfile, TProfile2D, TH2F, TFile, TH1D, TH2D
0006 import re
0007 import os
0008 
0009 
0010 def draw_line(lineList,x1,x2,y1,y2,width=1,style=1,color=1):
0011     from ROOT import TLine
0012     l=TLine(x1,y1,x2,y2)
0013     l.SetBit(ROOT.kCanDelete)
0014     l.SetLineWidth(width)
0015     l.SetLineStyle(style)
0016     l.SetLineColor(color)
0017     l.Draw()
0018     lineList.append(l)
0019 
0020 def draw_box(boxList,xl,xr,yl,yr,opacity=1,color=1,style=1001,lstyle=1,lw=3):
0021     from ROOT import TBox
0022     b=TBox(xl,yl,xr,yr)
0023     b.SetBit(ROOT.kCanDelete)
0024     b.SetFillStyle(style)
0025     b.SetFillColorAlpha(color, opacity)
0026     b.SetLineColor(color)
0027     b.SetLineWidth(lw)
0028     b.SetLineStyle(lstyle)
0029     b.Draw()
0030     boxList.append(b)
0031 
0032 def renderPluginBPIX(lineList,layer) :
0033     from ROOT import TCanvas,TLine
0034     nlad=[6,14,22,32]
0035     coordSign=[(-1,-1),(-1,1),(1,-1),(1,1)]
0036     for xsign,ysign in coordSign:
0037         xlow = xsign*0.5
0038         xhigh= xsign*(0.5+4)
0039         ylow = ysign*0.5
0040         yhigh= ysign*(0.5 + nlad[layer-1])
0041         # Outside Box
0042         draw_line(lineList,xlow,  xhigh,  ylow,  ylow) # bottom
0043         draw_line(lineList,xlow,  xhigh, yhigh, yhigh) # top
0044         draw_line(lineList,xlow,   xlow,  ylow, yhigh) # left
0045         draw_line(lineList,xhigh, xhigh,  ylow, yhigh) # right
0046         # Inner Horizontal lines
0047         for lad in range(nlad[layer-1]):
0048             lad+=1
0049             if lad != nlad[layer-1]:
0050                 y = ysign * (lad+0.5)
0051                 draw_line(lineList,xlow, xhigh,  y,  y)
0052             y = ysign * (lad);
0053             draw_line(lineList,xlow, xhigh,  y,  y, 1, 3);
0054         # Inner Vertical lines
0055         for mod in range(3) : 
0056             mod+=1
0057             x = xsign * (mod + 0.5);
0058             draw_line(lineList,x, x,  ylow,  yhigh);
0059 
0060         # Draw ROC0
0061         for mod in range(4):
0062             mod+=1
0063             for lad in range(nlad[layer-1]):
0064                 lad+=1
0065                 if ysign==1:
0066                     flipped = not(lad%2==0)
0067                 else :
0068                     flipped = not(lad%2==1)
0069                 if flipped : roc0_orientation = -1
0070                 else : roc0_orientation = 1
0071                 if xsign==-1 : roc0_orientation *= -1
0072                 if ysign==-1 : roc0_orientation *= -1
0073                 x1 = xsign * (mod+0.5)
0074                 x2 = xsign * (mod+0.5 - 1./8);
0075                 y1 = ysign * (lad)
0076                 y2 = ysign * (lad + roc0_orientation*1./2)
0077                 if layer == 1 and xsign == -1 :
0078                     x1 = xsign * (mod-0.5)
0079                     x2 = xsign * (mod-0.5 + 1./8)
0080                     y1 = ysign * (lad)
0081                     y2 = ysign * (lad - roc0_orientation*1./2)
0082 
0083                     draw_line(lineList,x1, x2, y1, y1, 1)
0084                     draw_line(lineList,x2, x2, y1, y2, 1)
0085                   
0086                 else:
0087                     draw_line(lineList,x1, x2, y1, y1, 1)
0088                     draw_line(lineList,x2, x2, y1, y2, 1)
0089 
0090 def maskBPixROC(boxList,xsign,ysign,layer,lad,mod,roc):
0091     if roc<8 : 
0092         rocShiftX=roc*1./8
0093         rocShiftY=0
0094     else : 
0095         rocShiftX=(15-roc)*1./8
0096         rocShiftY=1./2
0097     if ysign==1:
0098         flipped = not(lad%2==0)
0099     else :
0100         flipped = not(lad%2==1)
0101     if flipped : roc0_orientation = -1
0102     else : roc0_orientation = 1
0103     if xsign==-1 : roc0_orientation *= -1
0104     if ysign==-1 : roc0_orientation *= -1
0105     x1 = xsign * (mod+0.5-rocShiftX)
0106     x2 = xsign * (mod+0.5 - 1./8-rocShiftX);
0107     y1 = ysign * (lad-roc0_orientation*rocShiftY)
0108     y2 = ysign * (lad + roc0_orientation*1./2-roc0_orientation*rocShiftY)
0109     if layer == 1 and xsign == -1 :
0110         x1 = xsign * (mod-0.5)-rocShiftX
0111         x2 = xsign * (mod-0.5 + 1./8)-rocShiftX
0112         y1 = ysign * (lad +rocShiftY)
0113         y2 = ysign * (lad - roc0_orientation*1./2+rocShiftY)
0114     draw_box(boxList,min(x1,x2),max(x1,x2),min(y1, y2),max(y1,y2),0.75)
0115 
0116 
0117                   
0118 def renderPluginFPIX(lineList,ring) :
0119     from ROOT import TCanvas,TLine
0120     coordSign=[(-1,-1),(-1,1),(1,-1),(1,1)]
0121     for dsk in range(3) :
0122         dsk+=1
0123         for xsign,ysign in coordSign:            
0124             for bld in range(5+ring*6):
0125                 bld+=1
0126                 # Panel 2 has dashed mid-plane
0127                 x1      = xsign * (0.5 + dsk - 1)
0128                 x2      = xsign * (0.5 + dsk)
0129                 sign = ysign
0130                 y1      = ysign * (bld + sign*0.5)
0131                 y2      = ysign * (bld)
0132                 yp2_mid = ysign * (bld - sign*0.25)
0133                 y3      = ysign * (bld - sign*0.5)
0134                 draw_line(lineList,x1, x2, y1, y1)
0135                 draw_line(lineList,x1, x2, y2, y2)
0136                 draw_line(lineList,x1, x2, yp2_mid, yp2_mid,1,2)
0137                 draw_line(lineList,x1, x2, y3, y3)
0138                 # Vertical lines
0139                 x = xsign * (0.5 + dsk - 1)
0140                 draw_line(lineList,x,  x,  y1,  y2)
0141                 draw_line(lineList,x,  x,  y2,  y3)
0142                 if ring==2 :
0143                     x = xsign * (0.5 + dsk)
0144                     draw_line(lineList,x,  x,  y1,  y2)
0145                     draw_line(lineList,x,  x,  y2,  y3)
0146                 #Make a BOX around ROC 0
0147                 x1 = xsign * (0.5 + dsk - 1/8.)
0148                 x2 = xsign * (0.5 + dsk)
0149                 y1_p1 = ysign * (bld + sign*0.25)
0150                 y2_p1 = ysign * (bld + sign*0.25 + xsign*ysign*0.25)
0151                 draw_line(lineList,x1, x2, y1_p1, y1_p1, 1)
0152                 draw_line(lineList,x1, x1, y1_p1, y2_p1, 1)
0153                 y1_p2 = ysign * (bld - sign*0.25)
0154                 y2_p2 = ysign * (bld - sign*0.25 - xsign*ysign*0.25)
0155                 draw_line(lineList,x1, x2, y1_p2, y1_p2)
0156                 draw_line(lineList,x1, x1, y1_p2, y2_p2)
0157 
0158 def maskFPixROC(boxList,xsign,ysign,dsk,bld,pnl,roc) :
0159     from ROOT import TCanvas,TLine       
0160     if roc<8 : 
0161         rocShiftX=roc*1./8
0162         rocShiftY=0
0163     else : 
0164         rocShiftX=(15-roc)*1./8
0165         rocShiftY=1./4
0166     sign=ysign
0167     x1 = xsign * (0.5 + dsk - 1/8.-rocShiftX)
0168     x2 = xsign * (0.5 + dsk-rocShiftX)
0169     if pnl==1:
0170         y1 = ysign * (bld + sign*0.25)-xsign*rocShiftY
0171         y2 = ysign * (bld + sign*0.25 + xsign*ysign*0.25)-xsign*rocShiftY
0172     else:
0173         y1 = ysign * (bld - sign*0.25)+xsign*rocShiftY
0174         y2 = ysign * (bld - sign*0.25 - xsign*ysign*0.25)+xsign*rocShiftY
0175     draw_box(boxList,min(x1,x2),max(x1,x2),min(y1,y2),max(y1,y2),0.75)
0176 
0177 
0178 def dqm_get_dataset(server, match, run, type="offline_data"):
0179     datareq = urllib2.Request(('%s/data/json/samples?match=%s') % (server, match))
0180     datareq.add_header('User-agent', ident)
0181     # Get data                                                                                                                              
0182     data = eval(re.sub(r"\bnan\b", "0", urllib2.build_opener(X509CertOpen()).open(datareq).read()),
0183                        { "__builtins__": None }, {})
0184     ret = ""
0185     for l in data['samples']:
0186         if l['type'] == type:
0187             for x in l['items']:
0188                 if int(x['run']) == int(run):
0189                     ret=x['dataset']
0190                     break
0191     print(ret)
0192     return ret
0193 
0194 
0195 
0196 
0197 def main():
0198     import sys
0199     import os
0200     import ROOT
0201 
0202     if len(sys.argv) != 3:
0203         print("input files needed!")
0204         return
0205     else:
0206         filename=sys.argv[1]
0207         pclfile=sys.argv[2]
0208     print(filename+" -- "+filename[19:25])
0209     print(pclfile+" -- "+pclfile[19:25])
0210     runNum=filename[19:25]
0211 
0212     dir="DQMData/Run " + runNum + "/PixelPhase1/Run summary/Pahse1_MechanicalView/"
0213     dirERROR="DQMData/Run " + runNum + "/PixelPhase1/Run summary/SiPixelQualityPCL/BadROC_PCL/"
0214 
0215 
0216     dirBPix=dir + "PXBarrel/"
0217     dirFPix=dir + "PXForward/"
0218 
0219     hoccB="digi_occupancy_per_SignedModuleCoord_per_SignedLadderCoord_PXLayer_"
0220     hoccF="digi_occupancy_per_SignedDiskCoord_per_SignedBladePanelCoord_PXRing_"
0221     hdeadB="Dead Channels per ROC_per_SignedModuleCoord_per_SignedLadderCoord_PXLayer_"
0222     hdeadF="Dead Channels per ROC_per_SignedDiskCoord_per_SignedBladePanelCoord_PXRing_"
0223         
0224     ROOT.gROOT.SetBatch(1)    
0225     ROOT.gStyle.SetOptStat(0)
0226     ROOT.gStyle.SetPalette(52) #104 kTemperatureMap // 55 kRainBow // 97 kRust // 57 kBird
0227     color=[]
0228     for i in range(0,255):
0229         color.append(ROOT.TColor.GetColorPalette(i))
0230     ROOT.gStyle.SetPalette(1) #104 kTemperatureMap // 55 kRainBow
0231     ROOT.gStyle.SetNumberContours(128)
0232     rootf=ROOT.TFile(filename)
0233     rootp=ROOT.TFile(pclfile)
0234     
0235 
0236 
0237     c=ROOT.TCanvas("c","c",1250,1000)
0238     #BPIX
0239     print("----> Build maps for BPix")
0240     histOccList=[]
0241     histDeadList=[]
0242     for lyr in range(1,5):
0243         histOccList.append(rootf.FindObjectAny(hoccB+str(lyr)))
0244         histDeadList.append(rootp.Get(dirERROR+hdeadB+str(lyr)))
0245     for hist1, hist2 in zip(histOccList, histDeadList):
0246         if hist1 != None or hist2 !=None:
0247             hist1.Draw("colz")
0248             match=re.search('(?<=PXLayer_)[0-9]',hist1.GetName())
0249             if match != None and "per_SignedModuleCoord_per_SignedLadderCoord" in hist1.GetName():
0250                 lyr=int(match.group(0))
0251                 hist1.SetTitle("Digi Occupancy Layer {0}".format(lyr))
0252                 boxList=[]
0253                 lineList=[]
0254                 renderPluginBPIX(lineList,lyr)
0255                 lineWd=3
0256                 if lyr==4 :
0257                     lineWd=2
0258                 if lyr==1: 
0259                     tbmRoc=4
0260                 else:
0261                     tbmRoc=8
0262                 binTBM=[]
0263                 singleROC=0
0264                 maxx=hist2.GetMaximum()
0265                 for biny in range(1,hist2.GetNbinsY()+1):
0266                     if len(binTBM)!=0:
0267                         x1=hist2.GetXaxis().GetBinLowEdge(binTBM[0])
0268                         x2=hist2.GetXaxis().GetBinUpEdge(binTBM[len(binTBM)-1])
0269                         y1=hist2.GetYaxis().GetBinLowEdge(biny-1)
0270                         y2=hist2.GetYaxis().GetBinUpEdge(biny-1)
0271                         draw_box(boxList,x1,x2,y1,y2,0.2,color[100+int((224-100)*(1-hist2.GetBinContent(binTBM[0],biny-1)/maxx))],0,1,lineWd)
0272                     binTBM=[]
0273                     singleROC=0
0274                     for binx in range(1,hist2.GetNbinsX()+1):
0275                         if len(binTBM)==0:
0276                             check=True
0277                         elif hist2.GetBinContent(binx,biny)==hist2.GetBinContent(binTBM[len(binTBM)-1],biny):
0278                             check=True
0279                         else:
0280                             check=False
0281                         if hist2.GetBinContent(binx,biny)!=0 and check:
0282                             if len(binTBM)==0:
0283                                 binTBM.append(binx)
0284                             else:
0285                                 if len(binTBM)==(tbmRoc): 
0286                                     x1=hist2.GetXaxis().GetBinLowEdge(binTBM[0])
0287                                     x2=hist2.GetXaxis().GetBinUpEdge(binTBM[len(binTBM)-1])
0288                                     y1=hist2.GetYaxis().GetBinLowEdge(biny)
0289                                     y2=hist2.GetYaxis().GetBinUpEdge(biny)
0290                                     draw_box(boxList,x1,x2,y1,y2,0.2,color[100+int((224-100)*(1-hist2.GetBinContent(binTBM[0],biny)/maxx))],0,1,lineWd)
0291                                     binTBM=[]
0292                                     singleROC=0
0293                                     binTBM.append(binx)
0294                                 else:
0295                                     binTBM.append(binx)
0296                         else:
0297                             if len(binTBM)!=0:
0298                                 x1=hist2.GetXaxis().GetBinLowEdge(binTBM[0])
0299                                 x2=hist2.GetXaxis().GetBinUpEdge(binTBM[len(binTBM)-1])
0300                                 y1=hist2.GetYaxis().GetBinLowEdge(biny)
0301                                 y2=hist2.GetYaxis().GetBinUpEdge(biny)
0302                                 draw_box(boxList,x1,x2,y1,y2,0.2,color[100+int((224-100)*(1-hist2.GetBinContent(binTBM[0],biny)/maxx))],0,1,lineWd)
0303                                 binTBM=[]
0304                                 if hist2.GetBinContent(binx,biny)!=0:
0305                                     binTBM.append(binx)
0306                 c.SaveAs('MergedPCLDeadROC_BPix_Layer{0}_TBM.pdf'.format(lyr))
0307                 os.system('gs -dBATCH -dNOPAUSE -sDEVICE=png16m -dUseCropBox -sOutputFile=MergedPCLDeadROC_BPix_Layer{0}_TBM.png -r144 -q MergedPCLDeadROC_BPix_Layer{0}_TBM.pdf'.format(lyr))
0308                 os.system('rm -f MergedPCLDeadROC_BPix_Layer{0}_TBM.pdf'.format(lyr))
0309         else :
0310             print("Some Error in get the histograms for FPIX")
0311     #FPIX
0312     print("----> Build maps for FPix")
0313     for rng in range(1,3):
0314         histOccList.append(rootf.FindObjectAny(hoccF+str(rng)))
0315         histDeadList.append(rootp.Get(dirERROR+hdeadF+str(rng)))
0316     for hist1, hist2 in zip(histOccList, histDeadList):
0317         if hist1 != None or hist2 !=None:
0318             hist1.Draw("colz")
0319             match=re.search('(?<=PXRing_)[0-9]',hist1.GetName())
0320             if match != None and "per_SignedDiskCoord_per_SignedBladePanelCoord" in hist1.GetName():
0321                 ring=int(match.group(0))
0322                 hist1.SetTitle("Digi Occupancy Ring {0}".format(ring))
0323                 boxList=[]
0324                 lineList=[]
0325                 renderPluginFPIX(lineList,ring)
0326                 lineWd=3
0327                 if ring==2 :
0328                     lineWd=2
0329                 tbmRoc=8
0330                 binTBM=[]
0331                 maxx=hist2.GetMaximum()
0332                 for biny in range(1,hist2.GetNbinsY()+1):
0333                     if len(binTBM)!=0:
0334                         x1=hist2.GetXaxis().GetBinLowEdge(binTBM[0])
0335                         x2=hist2.GetXaxis().GetBinUpEdge(binTBM[len(binTBM)-1])
0336                         y1=hist2.GetYaxis().GetBinLowEdge(biny-1)
0337                         y2=hist2.GetYaxis().GetBinUpEdge(biny-1)
0338                         draw_box(boxList,x1,x2,y1,y2,0.2,color[100+int((224-100)*(1-hist2.GetBinContent(binTBM[0],biny-1)/maxx))],0,1,lineWd)
0339                     binTBM=[]
0340                     for binx in range(1,hist2.GetNbinsX()+1):
0341                         if len(binTBM)==0:
0342                             check=True
0343                         elif hist2.GetBinContent(binx,biny)==hist2.GetBinContent(binTBM[len(binTBM)-1],biny):
0344                             check=True
0345                         else:
0346                             check=False
0347                         if hist2.GetBinContent(binx,biny)!=0 and check:
0348                             if len(binTBM)==0:
0349                                 binTBM.append(binx)
0350                             else:
0351                                 if len(binTBM)==tbmRoc: 
0352                                     x1=hist2.GetXaxis().GetBinLowEdge(binTBM[0])
0353                                     x2=hist2.GetXaxis().GetBinUpEdge(binTBM[len(binTBM)-1])
0354                                     y1=hist2.GetYaxis().GetBinLowEdge(biny)
0355                                     y2=hist2.GetYaxis().GetBinUpEdge(biny)
0356                                     draw_box(boxList,x1,x2,y1,y2,0.2,color[100+int((224-100)*(1-hist2.GetBinContent(binTBM[0],biny)/maxx))],0,1,lineWd)
0357                                     binTBM=[]
0358                                     binTBM.append(binx)
0359                                 else:
0360                                     binTBM.append(binx)
0361                         else:
0362                             if len(binTBM)!=0:
0363                                 x1=hist2.GetXaxis().GetBinLowEdge(binTBM[0])
0364                                 x2=hist2.GetXaxis().GetBinUpEdge(binTBM[len(binTBM)-1])
0365                                 y1=hist2.GetYaxis().GetBinLowEdge(biny)
0366                                 y2=hist2.GetYaxis().GetBinUpEdge(biny)
0367                                 draw_box(boxList,x1,x2,y1,y2,0.2,color[100+int((224-100)*(1-hist2.GetBinContent(binTBM[0],biny)/maxx))],0,1,lineWd)
0368                                 binTBM=[]
0369                                 if hist2.GetBinContent(binx,biny)!=0:
0370                                     binTBM.append(binx)
0371                 c.SaveAs('MergedPCLDeadROC_FPix_Ring{0}_TBM.pdf'.format(ring))
0372                 os.system('gs -dBATCH -dNOPAUSE -sDEVICE=png16m -dUseCropBox -sOutputFile=MergedPCLDeadROC_FPix_Ring{0}_TBM.png -r144 -q MergedPCLDeadROC_FPix_Ring{0}_TBM.pdf'.format(ring))
0373                 os.system('rm -f MergedPCLDeadROC_FPix_Ring{0}_TBM.pdf'.format(ring))
0374 
0375         else :
0376             print("Some Error in get the histograms for FPIX")
0377 
0378 
0379 
0380 
0381 if __name__ == '__main__':
0382     main()