Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:28

0001 from __future__ import print_function
0002 from __future__ import absolute_import
0003 import ROOT
0004 from .fitResidual import fitResidual
0005 from .drawHistoAllChambers import drawHisto
0006 
0007 def plot(fileName,sl,dir='DQMData/Run 1/DT/Run summary/DTCalibValidation',option="HISTOPE1",draw=True):
0008 
0009     slType = sl
0010     slStr = "SL%d" % slType
0011     verbose = False
0012     nSigmas = 2
0013 
0014     ROOT.TH1.AddDirectory(False)
0015 
0016     file = ROOT.TFile(fileName,'read')
0017 
0018     wheels = (-2,-1,0,1,2)
0019     stations = (1,2,3,4)
0020 
0021     # (Wh-2 MB1 Sec1 ... Wh-2 MB1 Sec12 ... Wh-1 MB1 Sec1 ... Wh-1 MB1 Sec12 ...)
0022     # (Wh-2 MB2 Sec1 ... Wh-2 MB2 Sec12 ... Wh-1 MB2 Sec1 ... Wh-1 MB1 Sec12 ...) ...  
0023     nBins = 250
0024     if slType == 2: nBins = 180
0025     histoMean = ROOT.TH1F("h_ResMeanAll","Mean of residuals",nBins,0,nBins)
0026     histoSigma = ROOT.TH1F("h_ResSigmaAll","Sigma of residuals",nBins,0,nBins)
0027     for st in stations:
0028         nSectors = 12
0029         if st == 4: nSectors = 14
0030         if st == 4 and slType == 2: continue 
0031         if verbose: print("Station",st)
0032         for wh in wheels:
0033             if verbose: print("Wheel",wh) 
0034             for sec in range(1,nSectors+1):
0035                 if verbose: print("Sector",sec)
0036                 # Get histogram
0037                 histoName = "%s/Wheel%d/Station%d/Sector%d/hResDist_STEP3_W%d_St%d_Sec%d_%s" % (dir,wh,st,sec,wh,st,sec,slStr) 
0038                 print("Accessing",histoName)
0039                 histo = file.Get(histoName)
0040                 (histo,fitFunc) = fitResidual(histo,nSigmas,verbose)
0041                 fitMean = fitFunc.GetParameter(1)
0042                 fitMeanErr = fitFunc.GetParError(1)
0043                 fitSigma = fitFunc.GetParameter(2)
0044                 fitSigmaErr = fitFunc.GetParError(2)
0045 
0046                 binHistoNew = (st - 1)*60 + (wh + 2)*nSectors + sec
0047                 if verbose: print("Bin in summary histo",binHistoNew)
0048                 histoMean.SetBinContent(binHistoNew,fitMean)
0049                 histoMean.SetBinError(binHistoNew,fitMeanErr)
0050                 histoSigma.SetBinContent(binHistoNew,fitSigma)
0051                 histoSigma.SetBinError(binHistoNew,fitSigmaErr)  
0052 
0053                 if sec == 1:
0054                     label = "Wheel %d" % wh
0055                     if wh == -2: label += " MB%d" % st  
0056                     histoMean.GetXaxis().SetBinLabel(binHistoNew,label) 
0057                     histoSigma.GetXaxis().SetBinLabel(binHistoNew,label)
0058 
0059     objectsMean = drawHisto(histoMean,title="Mean of residuals (cm)",
0060                                       ymin=-0.1,ymax=0.1,option=option,draw=draw)
0061     objectsSigma = drawHisto(histoSigma,title="Sigma of residuals (cm)",
0062                                       ymin=0.,ymax=0.15,option=option,draw=draw)
0063 
0064     return (objectsMean,objectsSigma)