Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:38:34

0001 import ROOT
0002 ROOT.gROOT.SetBatch(True)
0003 from setTDRStyle import setTDRStyle
0004 
0005 
0006 from granularity import *
0007 
0008 
0009 class ValidationPlotter:
0010     def __init__(self):
0011         setTDRStyle()
0012         self.inFiles = {}
0013         self.labels = {}
0014         self.colors = {}
0015         self.markers = {}
0016         self.outPath = None
0017         self.granularity = standardGranularity
0018         self.order = []
0019         
0020     def addInputFile(self, label, inFile, color=None, marker=20):
0021         self.order.append(label)
0022         self.inFiles[label] = inFile
0023         self.labels[label] = label
0024         self.markers[label] = marker
0025         if color != None:
0026             self.colors[label] = color
0027         else:
0028             # choose first not occupied color (other than white)
0029             for autoColor in range(1,100):
0030                 if autoColor not in self.colors.values() and not autoColor == 10:
0031                     self.colors[label] = autoColor
0032                     break
0033     
0034     def setGranularity(self, granularity):
0035         self.granularity = granularity
0036     
0037     def setOutputPath(self, outPath):
0038         self.outPath = outPath
0039         
0040     def convertName(self, name):
0041         out = name.replace("Bpix", "BPIX")
0042         out = out.replace("Fpix", "FPIX")
0043         out = out.replace("Plus", "+")
0044         out = out.replace("Minus", "-")
0045         out = out.replace("Fpix", "FPIX")
0046         out = out.replace("Tib", "TIB")
0047         out = out.replace("Tob", "TOB")
0048         out = out.replace("Tid", "TID")
0049         out = out.replace("Tec", "TEC")
0050         out = out.replace("Layer", " L")
0051         out = out.replace("Ring", " R")
0052         out = out.replace("Stereo", "S")
0053         out = out.replace("Rphi", "R") # different from Ring, this one does not add a space in front
0054         out = out.replace("In", "i")
0055         out = out.replace("Out", "o")
0056         return out
0057     
0058     def plotHist(self, folder, name, title, hists, twoDimensional=False):
0059         self.canvas = ROOT.TCanvas("canvas", "canvas", ROOT.gStyle.GetCanvasDefW(),ROOT.gStyle.GetCanvasDefH())
0060         ROOT.gPad.SetRightMargin(0.10)
0061         if twoDimensional:
0062             ROOT.gPad.SetRightMargin(0.2)
0063             
0064         
0065         legend = ROOT.TLegend(0.2,0.7,0.9,0.90)
0066         legend.SetFillColor(0)
0067         legend.SetFillStyle(0)
0068         legend.SetTextSize(0.03)
0069         legend.SetMargin(0.15)
0070         legend.SetNColumns(2)
0071         legend.SetBorderSize(0)
0072         
0073         normalize = False
0074         if len (hists) > 1:
0075             normalize = True
0076         
0077 
0078         firstHist = True
0079         scaleHist = None
0080         maximum = 0
0081         for hist, label in hists:
0082             
0083             n = int(hist.Integral())
0084             mu = hist.GetMean()
0085             sigma = hist.GetRMS()
0086             
0087             if normalize:
0088                 hist.Scale(1./hist.Integral())
0089             
0090             if hist.GetMaximum() > maximum:
0091                 maximum = hist.GetMaximum()
0092             
0093             if firstHist:
0094                 scaleHist = hist
0095                 addDraw = ""
0096                 firstHist = False
0097             else:
0098                 addDraw = "same"
0099             
0100             if not twoDimensional:
0101                 if self.markers[label] != 0:
0102                     drawMode = "P0%s"%(addDraw)
0103                 else:
0104                     drawMode = "hist%s"%(addDraw)
0105             else:
0106                 drawMode = "COLZ"
0107                 
0108             scaleHist.SetMaximum(maximum*1.5)
0109             scaleHist.SetMinimum(0)
0110             
0111             hist.Draw(drawMode)
0112             
0113             legText = "#splitline{{{label}}}{{N={n:.2E},#mu={mu:.2f},RMS={sigma:.2f}}}".format(label=label,n=n, mu=mu,sigma=sigma)
0114             
0115             if self.markers[label] != 0:
0116                 legend.AddEntry(hist, legText, "p")
0117             else:
0118                 legend.AddEntry(hist, legText, "l")
0119         
0120         if not twoDimensional:
0121             legend.Draw()
0122         
0123         cmsText = ROOT.TLatex(0.16,0.96,title)
0124         cmsText.SetTextFont(42)
0125         cmsText.SetNDC()
0126         cmsText.Draw("same")
0127         
0128         import os
0129         if not os.path.isdir(self.outPath):
0130             os.makedirs(self.outPath)
0131         if not os.path.isdir(self.outPath+"/"+folder):
0132             os.makedirs(self.outPath+"/"+folder)
0133         self.canvas.SaveAs("{}/{}/{}.pdf".format(self.outPath, folder, name))
0134         self.canvas = None
0135             
0136         
0137         
0138     def makeResidualPlot(self, sectorNumber, coordinate):
0139         # residual
0140         hists = []
0141         for label in self.order:
0142             fi = ROOT.TFile(self.inFiles[label], "READ")
0143             hist = fi.Get("ApeEstimator1/Sector_{sectorNumber}/Results/h_Res{coordinate}".format(sectorNumber=sectorNumber, coordinate=coordinate))
0144             hist.SetLineColor(self.colors[label])
0145             hist.SetMarkerColor(self.colors[label])
0146             hist.SetMarkerStyle(self.markers[label])
0147             hist.SetDirectory(0)
0148             if self.markers[label] == 0:
0149                 hist.SetMarkerSize(0)
0150                 hist.SetLineWidth(2)
0151                 
0152             hists.append((hist, label))
0153             nameHist = fi.Get("ApeEstimator1/Sector_{sectorNumber}/z_name".format(sectorNumber=sectorNumber))
0154             nameHist.SetDirectory(0)
0155             title = self.convertName(nameHist.GetTitle())
0156             fi.Close()
0157             
0158         name = "Sector_{sectorNumber}_Res{coordinate}".format(sectorNumber=sectorNumber, coordinate=coordinate)
0159         
0160         self.plotHist("residuals", name, title, hists)
0161         
0162         
0163         # normalized residual
0164         hists = []
0165         for label in self.order:
0166             fi = ROOT.TFile(self.inFiles[label], "READ")
0167             hist = fi.Get("ApeEstimator1/Sector_{sectorNumber}/Results/h_NorRes{coordinate}".format(sectorNumber=sectorNumber, coordinate=coordinate))
0168             hist.SetLineColor(self.colors[label])
0169             hist.SetMarkerColor(self.colors[label])
0170             hist.SetMarkerStyle(self.markers[label])
0171             hist.SetDirectory(0)
0172             if self.markers[label] == 0:
0173                 hist.SetMarkerSize(0)
0174                 hist.SetLineWidth(2)
0175                 
0176             hists.append((hist, label))
0177             nameHist = fi.Get("ApeEstimator1/Sector_{sectorNumber}/z_name".format(sectorNumber=sectorNumber))
0178             nameHist.SetDirectory(0)
0179             title = self.convertName(nameHist.GetTitle())
0180             fi.Close()
0181         name = "Sector_{sectorNumber}_NorRes{coordinate}".format(sectorNumber=sectorNumber, coordinate=coordinate)
0182         
0183         self.plotHist("residuals", name, title, hists)
0184         
0185     def makeTrackPlot(self,histName, twoDimensional=False):
0186         hists = []
0187         for label in self.order:
0188             fi = ROOT.TFile(self.inFiles[label], "READ")
0189             hist = fi.Get("ApeEstimator2/TrackVariables/{histName}".format(histName=histName))
0190             hist.SetLineColor(self.colors[label])
0191             hist.SetMarkerColor(self.colors[label])
0192             hist.SetMarkerStyle(self.markers[label])
0193             hist.SetDirectory(0)
0194             if self.markers[label] == 0:
0195                 hist.SetMarkerSize(0)
0196                 hist.SetLineWidth(2)
0197             
0198             if twoDimensional:
0199                 self.plotHist("tracks", histName+"_"+label, label, [(hist, label),], twoDimensional=True)
0200             else:
0201                 hists.append((hist, label))
0202         if len(hists) > 0:
0203             self.plotHist("tracks", histName, histName, hists, twoDimensional=twoDimensional)
0204             
0205     
0206     
0207     def draw(self):
0208         for coordinate in self.granularity.sectors.keys():
0209             rangeList = self.granularity.sectors[coordinate]
0210             for first, last in rangeList:
0211                 for i in range(first, last+1):
0212                     self.makeResidualPlot(i, coordinate)
0213                     
0214         self.makeTrackPlot("h_hitsSize")
0215         self.makeTrackPlot("h_hitsValid")
0216         self.makeTrackPlot("h_hitsInvalid")
0217         self.makeTrackPlot("h_hits2D")
0218         self.makeTrackPlot("h_layersMissed")
0219         self.makeTrackPlot("h_hitsPixel")
0220         self.makeTrackPlot("h_hitsStrip")
0221         self.makeTrackPlot("h_charge")
0222         self.makeTrackPlot("h_chi2")
0223         self.makeTrackPlot("h_ndof")
0224         self.makeTrackPlot("h_norChi2")
0225         self.makeTrackPlot("h_prob")
0226         self.makeTrackPlot("h_eta")
0227         self.makeTrackPlot("h_etaErr")
0228         self.makeTrackPlot("h_theta")
0229         self.makeTrackPlot("h_phi")
0230         self.makeTrackPlot("h_phiErr")
0231         self.makeTrackPlot("h_d0Beamspot")
0232         self.makeTrackPlot("h_d0BeamspotErr")
0233         self.makeTrackPlot("h_dz")
0234         self.makeTrackPlot("h_dzErr")
0235         self.makeTrackPlot("h_pt")
0236         self.makeTrackPlot("h_ptErr")
0237         self.makeTrackPlot("h_meanAngle")
0238         self.makeTrackPlot("h_hitsGood")
0239         # these don't have several histograms in one plot
0240         self.makeTrackPlot("h2_meanAngleVsHits", twoDimensional=True)
0241         self.makeTrackPlot("h2_hitsGoodVsHitsValid", twoDimensional=True)
0242         self.makeTrackPlot("h2_hitsPixelVsEta", twoDimensional=True)
0243         self.makeTrackPlot("h2_hitsPixelVsTheta", twoDimensional=True)
0244         self.makeTrackPlot("h2_hitsStripVsEta", twoDimensional=True)
0245         self.makeTrackPlot("h2_hitsStripVsTheta", twoDimensional=True)
0246         self.makeTrackPlot("h2_ptVsTheta", twoDimensional=True)
0247         
0248         
0249 
0250 def main():
0251     pass
0252 
0253 
0254 if __name__ == "__main__":
0255     main()