Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:35

0001 ##########################################################################
0002 # Creates histograms of the modules of one structure. and returns them as
0003 # a list of PlotData objects.
0004 ##
0005 
0006 from builtins import range
0007 import logging
0008 
0009 import ROOT
0010 ROOT.PyConfig.IgnoreCommandLineOptions = True
0011 ROOT.gROOT.SetBatch()
0012 
0013 import Alignment.MillePedeAlignmentAlgorithm.mpsvalidate.subModule as mpsv_subModule
0014 import Alignment.MillePedeAlignmentAlgorithm.mpsvalidate.classes as mpsv_classes
0015 import Alignment.MillePedeAlignmentAlgorithm.mpsvalidate.style as mpsv_style
0016 
0017 
0018 def plot(MillePedeUser, alignables, config):
0019     logger = logging.getLogger("mpsvalidate")
0020 
0021     alignables.create_list(MillePedeUser)
0022 
0023     # number of bins to start
0024     numberOfBins = 10000
0025 
0026     ######################################################################
0027     # initialize data hierarchy
0028     # plots[mode][struct]
0029     #
0030 
0031     plots = []
0032     # loop over mode
0033     for modeNumber, mode in enumerate(["xyz", "rot", "dist"]):
0034         plots.append([])
0035         # loop over structures
0036         for structNumber, struct in enumerate(alignables.structures):
0037             plots[modeNumber].append(mpsv_classes.PlotData(mode))
0038 
0039     # initialize histograms
0040     for modeNumber, mode in enumerate(["xyz", "rot", "dist"]):
0041         for structNumber, struct in enumerate(alignables.structures):
0042             # use a copy for shorter name
0043             plot = plots[modeNumber][structNumber]
0044 
0045             for i in range(3):
0046                 if (mode == "xyz"):
0047                     plot.histo.append(ROOT.TH1F("{0} {1} {2}".format(struct.get_name(), plot.xyz[
0048                                       i], mode), "", numberOfBins, -1000, 1000))
0049                 else:
0050                     plot.histo.append(ROOT.TH1F("{0} {1} {2}".format(struct.get_name(), plot.xyz[
0051                                       i], mode), "", numberOfBins, -0.1, 0.1))
0052 
0053                 if (plot.unit!=""):
0054                     plot.histo[i].SetXTitle("#Delta"+plot.xyz[i]+" ["+plot.unit+"]")
0055                 else:
0056                     plot.histo[i].SetXTitle("#Delta"+plot.xyz[i])
0057                 plot.histo[i].SetYTitle("number of alignables")
0058                 plot.histo[i].GetXaxis().SetTitleOffset(0.8)
0059                 plot.histoAxis.append(plot.histo[i].GetXaxis())
0060 
0061             # add labels
0062             plot.title = ROOT.TPaveLabel(
0063                 0.1, 0.8, 0.9, 0.9, "Module: {0} {1}".format(struct.get_name(), mode))
0064             plot.text = ROOT.TPaveText(0.05, 0.1, 0.95, 0.75)
0065             plot.text.SetTextAlign(12)
0066             plot.text.SetTextSizePixels(20)
0067 
0068             # save copy
0069             plots[modeNumber][structNumber] = plot
0070 
0071     ######################################################################
0072     # fill histogram
0073     #
0074 
0075     for line in MillePedeUser:
0076         # is module ?
0077         if (line.ObjId == 1):
0078             for modeNumber, mode in enumerate(["xyz", "rot", "dist"]):
0079                 for structNumber, struct in enumerate(alignables.structures):
0080                     # use a copy for shorter name
0081                     plot = plots[modeNumber][structNumber]
0082 
0083                     # module in struct ?
0084                     if (struct.contains_detid(line.Id)):
0085                         for i in range(3):
0086                             if (abs(line.Par[plot.data[i]]) != 999999):
0087                                 # transform xyz data from cm to #mu m
0088                                 if (mode == "xyz"):
0089                                     plot.histo[i].Fill(
0090                                         10000 * line.Par[plot.data[i]])
0091                                 else:
0092                                     plot.histo[i].Fill(line.Par[plot.data[i]])
0093 
0094                     # save copy
0095                     plots[modeNumber][structNumber] = plot
0096 
0097     ######################################################################
0098     # find the best range
0099     #
0100 
0101     for modeNumber, mode in enumerate(["xyz", "rot", "dist"]):
0102         for structNumber, struct in enumerate(alignables.structures):
0103             # use a copy for shorter name
0104             plot = plots[modeNumber][structNumber]
0105 
0106             for i in range(3):
0107                 # get first and last bin with content and chose the one which
0108                 # has a greater distance to the center
0109                 if (abs(numberOfBins // 2 - plot.histo[i].FindFirstBinAbove()) > abs(plot.histo[i].FindLastBinAbove() - numberOfBins // 2)):
0110                     plot.maxBinShift[i] = abs(
0111                         numberOfBins // 2 - plot.histo[i].FindFirstBinAbove())
0112                     # set the maxShift value
0113                     plot.maxShift[i] = plot.histo[i].GetBinCenter(
0114                         plot.histo[i].FindFirstBinAbove())
0115                 else:
0116                     plot.maxBinShift[i] = abs(
0117                         plot.histo[i].FindLastBinAbove() - numberOfBins // 2)
0118                     # set the maxShift value
0119                     plot.maxShift[i] = plot.histo[i].GetBinCenter(
0120                         plot.histo[i].FindLastBinAbove())
0121                 # skip empty histogram
0122                 if (abs(plot.maxBinShift[i]) == numberOfBins // 2 + 1):
0123                     plot.maxBinShift[i] = 0
0124 
0125             # three types of ranges
0126 
0127             # 1. multiple of standard dev
0128             if (config.rangemode == "stddev"):
0129                 for i in range(3):
0130                     if (plot.histo[i].GetEntries() != 0 and plot.histo[i].GetStdDev() != 0):
0131                         # if the plotrange is much bigger than the standard
0132                         # deviation use config.widthstdev * StdDev als Range
0133                         if (max(plot.maxShift) / plot.histo[i].GetStdDev() > config.defpeak):
0134                             # corresponding bin config.widthstdev*StdDev
0135                             binShift = int(plot.histo[i].FindBin(
0136                                 config.widthstddev * plot.histo[i].GetStdDev()) - numberOfBins / 2)
0137                         else:
0138                             binShift = max(plot.maxBinShift)
0139 
0140                         # save used binShift
0141                         plot.binShift[i] = binShift
0142 
0143             # 2. show all
0144             if (config.rangemode == "all"):
0145                 for i in range(3):
0146                     plot.binShift[i] = plot.maxBinShift[i]
0147 
0148             # 3. use given ranges
0149             if (config.rangemode == "given"):
0150                 for i in range(3):
0151                     if (mode == "xyz"):
0152                         valuelist = config.rangexyzM
0153                     if (mode == "rot"):
0154                         valuelist = config.rangerotM
0155                     if (mode == "dist"):
0156                         valuelist = config.rangedistM
0157 
0158                     for value in valuelist:
0159                         # maximum smaller than given value
0160                         if (abs(plot.maxShift[i]) < value):
0161                             binShift = value
0162                             break
0163                     # if not possible, force highest
0164                     if (abs(plot.maxShift[i]) > valuelist[-1]):
0165                         binShift = valuelist[-1]
0166                     # calculate binShift
0167                     plot.binShift[i] = int(
0168                         binShift / plot.histo[i].GetBinWidth(1))
0169 
0170             # all plot the same range
0171             if (config.samerange == 1):
0172                 for i in range(3):
0173                     plot.binShift[i] = max(plot.binShift)
0174 
0175             # save used range
0176             for i in range(3):
0177                 plot.usedRange[i] = plot.binShift[i]
0178 
0179             # count entries which are not shown anymore
0180             for i in range(3):
0181                 # bin 1 to begin of histogram
0182                 for j in range(1, numberOfBins // 2 - plot.binShift[i]):
0183                     plot.hiddenEntries[i] += plot.histo[i].GetBinContent(j)
0184                 # from the end of shown bins to the end of histogram
0185                 for j in range(numberOfBins // 2 + plot.binShift[i], plot.histo[i].GetNbinsX()):
0186                     plot.hiddenEntries[i] += plot.histo[i].GetBinContent(j)
0187 
0188             # apply new range
0189             for i in range(3):
0190                 if (plot.histo[i].GetEntries() != 0):
0191                     # merge bins, ca. 100 should be visible in the resulting
0192                     # plot
0193                     mergeNumberBins = plot.binShift[i]
0194                     # skip empty histogram
0195                     if (mergeNumberBins != 0):
0196                         # the 2*maxBinShift bins should shrink to 100 bins
0197                         mergeNumberBins = int(
0198                             2. * mergeNumberBins / config.numberofbins)
0199                         # the total number of bins should be dividable by the
0200                         # bins shrinked together
0201                         if (mergeNumberBins == 0):
0202                             mergeNumberBins = 1
0203                         while (numberOfBins % mergeNumberBins != 0 and mergeNumberBins != 1):
0204                             mergeNumberBins -= 1
0205 
0206                         # Rebin and save new created histogram and axis
0207                         plot.histo[i] = plot.histo[i].Rebin(mergeNumberBins)
0208                         plot.histoAxis[i] = plot.histo[i].GetXaxis()
0209 
0210                         # set view range. it is important to note that the number of bins have changed with the rebinning
0211                         # the total number and the number of shift must be
0212                         # corrected with / mergeNumberBins
0213                         plot.histoAxis[i].SetRange(int(numberOfBins / (2 * mergeNumberBins) - plot.binShift[
0214                                                    i] / mergeNumberBins), int(numberOfBins / (2 * mergeNumberBins) + plot.binShift[i] / mergeNumberBins))
0215 
0216             # error if shift is bigger than limit
0217             limit = config.limit[mode]
0218             for i in range(3):
0219                 # skip empty
0220                 if (plot.histo[i].GetEntries() > 0):
0221                     if (plot.unit!=""):
0222                         plot.text.AddText("max. shift {0}: {1:.2} {2}".format(plot.xyz[i], plot.maxShift[i], plot.unit))
0223                         if (abs(plot.maxShift[i]) > limit):
0224                             plot.text.AddText("! {0} shift bigger than {1} {2}".format(plot.xyz[i], limit, plot.unit))
0225                     else:
0226                         plot.text.AddText("max. shift {0}: {1:.2}".format(plot.xyz[i], plot.maxShift[i]))
0227                         if (abs(plot.maxShift[i]) > limit):
0228                             plot.text.AddText("! {0} shift bigger than {1}".format(plot.xyz[i], limit))
0229 
0230                     if (plot.hiddenEntries[i] != 0):
0231                         plot.text.AddText("! {0}: {1} outlier !".format(
0232                             plot.xyz[i], int(plot.hiddenEntries[i])))
0233 
0234             # save copy
0235             plots[modeNumber][structNumber] = plot
0236 
0237     ######################################################################
0238     # make the plots
0239     #
0240 
0241     # show the skewness in the legend
0242     ROOT.gStyle.SetOptStat("emrs")
0243 
0244     for modeNumber, mode in enumerate(["xyz", "rot", "dist"]):
0245         for structNumber, struct in enumerate(alignables.structures):
0246             # use a copy for shorter name
0247             plot = plots[modeNumber][structNumber]
0248 
0249             canvas = ROOT.TCanvas("canvasModules{0}_{1}".format(
0250                 struct.get_name(), mode), "Parameter", 300, 0, 800, 600)
0251             canvas.Divide(2, 2)
0252 
0253             canvas.cd(1)
0254             plot.title.Draw()
0255             plot.text.Draw()
0256 
0257             # draw identification
0258             ident = mpsv_style.identification(config)
0259             ident.Draw()
0260 
0261             # is there any plot?
0262             plotNumber = 0
0263 
0264             # loop over coordinates
0265             for i in range(3):
0266                 if(plot.histo[i].GetEntries() > 0):
0267                     plotNumber += 1
0268                     canvas.cd(i + 2)
0269                     mpsv_style.setstatsize(canvas, plot.histo[i], config)
0270                     plot.histo[i].DrawCopy()
0271 
0272             if (plotNumber == 0):
0273                 break
0274 
0275             canvas.Update()
0276 
0277             # save as pdf
0278             canvas.Print(
0279                 "{0}/plots/pdf/modules_{1}_{2}.pdf".format(config.outputPath, mode, struct.get_name()))
0280 
0281             # export as png
0282             image = ROOT.TImage.Create()
0283             image.FromPad(canvas)
0284             image.WriteImage(
0285                 "{0}/plots/png/modules_{1}_{2}.png".format(config.outputPath, mode, struct.get_name()))
0286 
0287             # add to output list
0288             output = mpsv_classes.OutputData(plottype="mod", name=struct.get_name(),
0289                                              parameter=mode, filename="modules_{0}_{1}".format(mode, struct.get_name()))
0290             config.outputList.append(output)
0291 
0292     ######################################################################
0293     # make plots with substructure
0294     #
0295 
0296     if config.showsubmodule:
0297         alignables.create_children_list()
0298         for modeNumber, mode in enumerate(["xyz", "rot", "dist"]):
0299             for structNumber, struct in enumerate(alignables.structures):
0300                 # use a copy for shorter name
0301                 plot = plots[modeNumber][structNumber]
0302 
0303                 mpsv_subModule.plot(MillePedeUser, alignables,
0304                                     mode, struct, plot, config)