Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef ALIGNMENT_OFFLINEVALIDATION_PLOTALIGNNMENTVALIDATION_H_
0002 #define ALIGNMENT_OFFLINEVALIDATION_PLOTALIGNNMENTVALIDATION_H_
0003 
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "Alignment/OfflineValidation/interface/TkAlStyle.h"
0006 #include "Alignment/OfflineValidation/interface/TkOffTreeVariables.h"
0007 
0008 #include "Math/ProbFunc.h"
0009 
0010 #include "TAxis.h"
0011 #include "TCanvas.h"
0012 #include "TDirectory.h"
0013 #include "TDirectoryFile.h"
0014 #include "TF1.h"
0015 #include "TFile.h"
0016 #include "TFrame.h"
0017 #include "TGaxis.h"
0018 #include "TH2F.h"
0019 #include "THStack.h"
0020 #include "TKey.h"
0021 #include "TLatex.h"
0022 #include "TLegend.h"
0023 #include "TLegendEntry.h"
0024 #include "TPad.h"
0025 #include "TPaveStats.h"
0026 #include "TPaveText.h"
0027 #include "TProfile.h"
0028 #include "TRandom3.h"
0029 #include "TRegexp.h"
0030 #include "TROOT.h"
0031 #include "TString.h"
0032 #include "TStyle.h"
0033 #include "TSystem.h"
0034 #include "TTree.h"
0035 #include "TTreeReader.h"
0036 
0037 #include <algorithm>
0038 #include <cmath>
0039 #include <cstdio>
0040 #include <cstdlib>
0041 #include <exception>
0042 #include <fstream>
0043 #include <iostream>
0044 #include <memory>
0045 #include <sstream>
0046 #include <string>
0047 #include <vector>
0048 
0049 class TkOfflineVariables {
0050 public:
0051   TkOfflineVariables(std::string fileName, std::string baseDir, std::string legName = "", int color = 1, int style = 1);
0052   ~TkOfflineVariables();
0053   int getLineColor() { return lineColor; }
0054   int getLineStyle() { return lineStyle; }
0055   std::string getName() { return legendName; }
0056   TTree* getTree() { return tree; }
0057   TFile* getFile() { return file; }
0058   int getPhase() { return phase; }
0059 
0060 private:
0061   TFile* file;
0062   TTree* tree;
0063   int lineColor;
0064   int lineStyle;
0065   int phase;
0066   std::string legendName;
0067 };
0068 
0069 inline TkOfflineVariables::TkOfflineVariables(
0070     std::string fileName, std::string baseDir, std::string legName, int lColor, int lStyle) {
0071   lineColor = lColor;
0072   lineStyle = lStyle % 100;
0073   if (legName.empty()) {
0074     int start = 0;
0075     if (fileName.find('/'))
0076       start = fileName.find_last_of('/') + 1;
0077     int stop = fileName.find_last_of('.');
0078     legendName = fileName.substr(start, stop - start);
0079   } else {
0080     legendName = legName;
0081   }
0082 
0083   //fill the tree pointer
0084   file = TFile::Open(fileName.c_str());
0085   TDirectoryFile* d = nullptr;
0086   if (file->Get(baseDir.c_str())) {
0087     d = (TDirectoryFile*)file->Get(baseDir.c_str());
0088     if ((*d).Get("TkOffVal")) {
0089       tree = (TTree*)(*d).Get("TkOffVal");
0090     } else {
0091       std::cout << "no tree named TkOffVal" << std::endl;
0092       assert(false);
0093     }
0094     TDirectoryFile* d2 = (TDirectoryFile*)d->Get("Pixel");
0095     assert(d2);
0096     phase = (int)((bool)d2->Get("P1PXBBarrel_1"));
0097   } else {
0098     std::cout << "no directory named " << baseDir.c_str() << std::endl;
0099     assert(false);
0100   }
0101 }
0102 
0103 inline TkOfflineVariables::~TkOfflineVariables() { delete file; }
0104 
0105 class PlotAlignmentValidation {
0106 public:
0107   //PlotAlignmentValidation(TString *tmp);
0108   PlotAlignmentValidation(bool bigtext = false);
0109   PlotAlignmentValidation(
0110       const char* inputFile, std::string fileName = "", int lineColor = 1, int lineStyle = 1, bool bigtext = false);
0111   ~PlotAlignmentValidation();
0112   void loadFileList(const char* inputFile, std::string fileName = "", int lineColor = 2, int lineStyle = 1);
0113   void useFitForDMRplots(bool usefit = false);
0114   void legendOptions(TString options);
0115   void plotOutlierModules(const char* outputFileName = "OutlierModules.ps",
0116                           std::string plotVariable = "chi2PerDofX",
0117                           float chi2_cut = 10,
0118                           unsigned int minHits = 50);  //method dumps selected modules into ps file
0119   void plotSubDetResiduals(
0120       bool plotNormHisto = false,
0121       unsigned int subDetId =
0122           7); /**<subDetector number :1.TPB, 2.TBE+, 3.TBE-, 4.TIB, 5.TID+, 6.TID-, 7.TOB, 8.TEC+ or 9.TEC- */
0123   void plotDMR(const std::string& plotVar = "medianX",
0124                Int_t minHits = 50,
0125                const std::string& options = "plain",
0126                const std::string& filterName = "");
0127   /**<
0128   * plotVar=mean,meanX,meanY,median,rms etc., comma-separated list can be given; 
0129   * minHits=the minimum hits needed for module to appear in plot; 
0130   * options="plain" for regular DMR, "split" for inwards/outwards split, "layers" for layerwise DMR, "layer=N" for Nth layer, or combination of the previous (e.g. "split layers")
0131   * filterName=rootfile containing tree with module ids to be skipped in plotting (to be used for averaged plots or in debugging)
0132   */
0133   void plotSurfaceShapes(const std::string& options = "layers", const std::string& variable = "");
0134   void plotChi2(const char* inputFile);
0135   /**< plotSurfaceShapes: options="split","layers"/"layer","subdet" */
0136   void plotHitMaps();
0137   void setOutputDir(std::string dir);
0138   void setTreeBaseDir(std::string dir = "TrackerOfflineValidation");
0139 
0140   void residual_by_moduleID(unsigned int moduleid);
0141   int numberOfLayers(int phase, int subdetector);
0142   int maxNumberOfLayers(int subdetector);
0143 
0144   THStack* addHists(
0145       const TString& selection,
0146       const TString& residType = "xPrime",
0147       TLegend** myLegend = nullptr,
0148       bool printModuleIds = false,
0149       bool validforphase0 =
0150           false); /**<add hists fulfilling 'selection' on TTree; residType: xPrime,yPrime,xPrimeNorm,yPrimeNorm,x,y,xNorm; if (printModuleIds): cout DetIds */
0151 
0152   float twotailedStudentTTestEqualMean(float t, float v);
0153 
0154   /**< These are helpers for DMR plotting */
0155 
0156   struct DMRPlotInfo {
0157     std::string variable;
0158     int nbins;
0159     double min, max;
0160     int minHits;
0161     bool plotPlain, plotSplits, plotLayers;
0162     int subDetId, nLayers;
0163     THStack* hstack;
0164     TLegend* legend;
0165     TkOfflineVariables* vars;
0166     float maxY;
0167     TH1F* h;
0168     TH1F* h1;
0169     TH1F* h2;
0170     bool firsthisto;
0171     std::string filterName;
0172   };
0173 
0174 private:
0175   TList* getTreeList();
0176   std::string treeBaseDir;
0177 
0178   bool useFit_;
0179   bool showMean_;
0180   bool showRMS_;
0181   bool showMeanError_;
0182   bool showRMSError_;
0183   bool showModules_;
0184   bool showUnderOverFlow_;
0185   bool twolines_;
0186   bool bigtext_;
0187   const static TString summaryfilename;
0188   std::ofstream summaryfile;
0189   bool openedsummaryfile = false;
0190   TFile* rootsummaryfile;
0191 
0192   std::vector<double> vmean, vdeltamean, vrms, vmeanerror, vPValueEqualSplitMeans, vPValueMeanEqualIdeal,
0193       vPValueRMSEqualIdeal, vAlignmentUncertainty;
0194   double resampleTestOfEqualMeans(TH1F* h1, TH1F* h2, int numSamples);
0195   double resampleTestOfEqualRMS(TH1F* h1, TH1F* h2, int numSamples);
0196 
0197   void storeHistogramInRootfile(TH1* hist);
0198   TF1* fitGauss(TH1* hist, int color);
0199   /**< 
0200   * void plotBoxOverview(TCanvas &c1, TList &treeList,std::string plot_Var1a,std::string plot_Var1b, std::string plot_Var2, Int_t filenumber,Int_t minHits);
0201   * void plot1DDetailsSubDet(TCanvas &c1, TList &treeList, std::string plot_Var1a,std::string plot_Var1b, std::string plot_Var2, Int_t minHits);
0202   * void plot1DDetailsBarrelLayer(TCanvas &c1, TList &treeList, std::string plot_Var1a,std::string plot_Var1b, Int_t minHits);
0203   * void plot1DDetailsDiskWheel(TCanvas &c1, TList &treelist, std::string plot_Var1a,std::string plot_Var1b, Int_t minHits);
0204   */
0205   void plotSS(const std::string& options = "layers", const std::string& variable = "");
0206   void setHistStyle(TH1& hist, const char* titleX, const char* titleY, int color);
0207   void setTitleStyle(TNamed& h,
0208                      const char* titleX,
0209                      const char* titleY,
0210                      int subDetId,
0211                      bool isSurfaceDeformation = false,
0212                      TString secondline = "");
0213   void setNiceStyle();
0214   void setCanvasStyle(TCanvas& canv);
0215   void setLegendStyle(TLegend& leg);
0216   void scaleXaxis(TH1* hist, Int_t scale);
0217   TObject* findObjectFromCanvas(TCanvas* canv, const char* className, Int_t n = 1);
0218 
0219   TString outputFile;
0220   std::string outputDir;
0221   TList* sourcelist;
0222   std::vector<TkOfflineVariables*> sourceList;
0223   bool moreThanOneSource;
0224   std::string fileNames[10];
0225   int fileCounter;
0226 
0227   std::string getSelectionForDMRPlot(int minHits, int subDetId, int direction = 0, int layer = 0);
0228   std::string getVariableForDMRPlot(
0229       const std::string& histoname, const std::string& variable, int nbins, double min, double max);
0230   void setDMRHistStyleAndLegend(TH1F* h, DMRPlotInfo& plotinfo, int direction = 0, int layer = 0);
0231   void plotDMRHistogram(DMRPlotInfo& plotinfo, int direction = 0, int layer = 0, std::string subdet = "");
0232   void modifySSHistAndLegend(THStack* hs, TLegend* legend);
0233   void openSummaryFile();
0234   std::vector<TH1*> findmodule(TFile* f, unsigned int moduleid);
0235 };
0236 
0237 #endif  // ALIGNMENT_OFFLINEVALIDATION_PLOTALIGNNMENTVALIDATION_H_