CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PlotAlignmentValidation.h
Go to the documentation of this file.
1 #ifndef PLOTALIGNNMENTVALIDATION_H_
2 #define PLOTALIGNNMENTVALIDATION_H_
3 
4 #include <TStyle.h>
5 #include <TSystem.h>
6 #include <vector>
7 #include <memory>
8 #include <string>
9 #include <fstream>
10 #include <iostream>
11 #include <sstream>
12 #include <stdio.h>
13 #include <stdlib.h>
14 #include "TTree.h"
15 #include "TString.h"
16 #include "TDirectory.h"
17 #include "TCanvas.h"
18 #include "TFile.h"
19 #include "TDirectoryFile.h"
20 #include "TLegend.h"
21 #include "THStack.h"
22 #include <exception>
23 
24 // This line works only if we have a CMSSW environment...
26 
28 public:
29  TkOfflineVariables(std::string fileName, std::string baseDir, std::string legName="", int color=1, int style=1);
30  int getLineColor(){ return lineColor; };
31  int getLineStyle(){ return lineStyle; };
33  TTree* getTree(){ return tree; };
34  TFile* getFile(){ return file; };
35 private:
36  TFile* file;
37  TTree* tree;
38  int lineColor;
39  int lineStyle;
41 };
42 
43 
45 {
46  lineColor = lColor;
47  lineStyle = lStyle;
48  if (legName=="") {
49  int start = 0;
50  if (fileName.find('/') ) start =fileName.find_last_of('/')+1;
51  int stop = fileName.find_last_of('.');
52  legendName = fileName.substr(start,stop-start);
53  } else {
54  legendName = legName;
55  }
56 
57  //fill the tree pointer
58  file = TFile::Open( fileName.c_str() );
59  TDirectoryFile *d = 0;
60  if (file->Get( baseDir.c_str() ) ) {
61  d = (TDirectoryFile*)file->Get( baseDir.c_str() );
62  if ((*d).Get("TkOffVal")) {
63  tree = (TTree*)(*d).Get("TkOffVal");
64  } else {
65  std::cout<<"no tree named TkOffVal"<<std::endl;
66  }
67  } else {
68  std::cout<<"no directory named "<<baseDir.c_str()<<std::endl;
69  }
70 }
71 
73 public:
74  //PlotAlignmentValidation(TString *tmp);
75  PlotAlignmentValidation(const char *inputFile,std::string fileName="", int lineColor=1, int lineStyle=1);
77  void loadFileList(const char *inputFile, std::string fileName="", int lineColor=2, int lineStyle=1);
78  void useFitForDMRplots(bool usefit = false);
79  void plotOutlierModules(const char *outputFileName="OutlierModules.ps",std::string plotVariable = "chi2PerDofX" ,float chi2_cut = 10,unsigned int minHits = 50);//method dumps selected modules into ps file
80  void plotSubDetResiduals(bool plotNormHisto=false, unsigned int subDetId=7);//subDetector number :1.TPB, 2.TBE+, 3.TBE-, 4.TIB, 5.TID+, 6.TID-, 7.TOB, 8.TEC+ or 9.TEC-
81  void plotDMR(const std::string& plotVar="medianX",Int_t minHits = 50, const std::string& options = "plain"); // plotVar=mean,meanX,meanY,median,rms etc., comma-separated list can be given; minHits=the minimum hits needed for module to appear in plot; 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")
82  void plotSurfaceShapes(const std::string& options = "layers",const std::string& variable="");
83  // plotSurfaceShapes: options="split","layers"/"layer","subdet"
84  void plotHitMaps();
86  void setTreeBaseDir( std::string dir = "TrackerOfflineValidationStandalone");
87 
88  THStack* addHists(const char *selection, const TString &residType = "xPrime", bool printModuleIds = false);//add hists fulfilling 'selection' on TTree; residType: xPrime,yPrime,xPrimeNorm,yPrimeNorm,x,y,xNorm; if (printModuleIds): cout DetIds
89 
90 private :
91  TList getTreeList();
93 
94  bool useFit_;
95 
96  std::pair<float,float> fitGauss(TH1 *hist,int color);
97  //void plotBoxOverview(TCanvas &c1, TList &treeList,std::string plot_Var1a,std::string plot_Var1b, std::string plot_Var2, Int_t filenumber,Int_t minHits);
98  //void plot1DDetailsSubDet(TCanvas &c1, TList &treeList, std::string plot_Var1a,std::string plot_Var1b, std::string plot_Var2, Int_t minHits);
99  //void plot1DDetailsBarrelLayer(TCanvas &c1, TList &treeList, std::string plot_Var1a,std::string plot_Var1b, Int_t minHits);
100  //void plot1DDetailsDiskWheel(TCanvas &c1, TList &treelist, std::string plot_Var1a,std::string plot_Var1b, Int_t minHits);
101  void plotSS(const std::string& options = "layers",const std::string& variable="");
102  void setHistStyle( TH1& hist,const char* titleX, const char* titleY, int color);
103  void setTitleStyle( TNamed& h,const char* titleX, const char* titleY, int subDetId);
104  void setNiceStyle();
105  void setCanvasStyle( TCanvas& canv );
106  void setLegendStyle( TLegend& leg );
107 
108  TString outputFile;
110  TList *sourcelist;
111  std::vector<TkOfflineVariables*> sourceList;
115 
116  // These are helpers for DMR plotting
117 
118  struct DMRPlotInfo {
120  int nbins;
121  double min, max;
122  int minHits;
125  THStack* hstack;
126  TLegend* legend;
128  float maxY;
129  TH1F* h;
130  TH1F* h1;
131  TH1F* h2;
133  };
134 
135  std::string getSelectionForDMRPlot(int minHits, int subDetId, int direction = 0, int layer = 0);
136  std::string getVariableForDMRPlot(const std::string& histoname, const std::string& variable,
137  int nbins, double min, double max);
138  void setDMRHistStyleAndLegend(TH1F* h, DMRPlotInfo& plotinfo, int direction = 0, int layer = 0);
139  void plotDMRHistogram(DMRPlotInfo& plotinfo, int direction = 0, int layer = 0);
140 
141 };
142 
143 #endif // PLOTALIGNNMENTVALIDATION_H_
void plotSurfaceShapes(const std::string &options="layers", const std::string &variable="")
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
void setHistStyle(TH1 &hist, const char *titleX, const char *titleY, int color)
std::pair< float, float > fitGauss(TH1 *hist, int color)
selection
main part
Definition: corrVsCorr.py:98
void plotSS(const std::string &options="layers", const std::string &variable="")
std::string getSelectionForDMRPlot(int minHits, int subDetId, int direction=0, int layer=0)
#define min(a, b)
Definition: mlp_lapack.h:161
void loadFileList(const char *inputFile, std::string fileName="", int lineColor=2, int lineStyle=1)
void plotOutlierModules(const char *outputFileName="OutlierModules.ps", std::string plotVariable="chi2PerDofX", float chi2_cut=10, unsigned int minHits=50)
void plotDMRHistogram(DMRPlotInfo &plotinfo, int direction=0, int layer=0)
void setTreeBaseDir(std::string dir="TrackerOfflineValidationStandalone")
void setLegendStyle(TLegend &leg)
void plotSubDetResiduals(bool plotNormHisto=false, unsigned int subDetId=7)
std::string getVariableForDMRPlot(const std::string &histoname, const std::string &variable, int nbins, double min, double max)
THStack * addHists(const char *selection, const TString &residType="xPrime", bool printModuleIds=false)
TkOfflineVariables(std::string fileName, std::string baseDir, std::string legName="", int color=1, int style=1)
const T & max(const T &a, const T &b)
void plotDMR(const std::string &plotVar="medianX", Int_t minHits=50, const std::string &options="plain")
void setTitleStyle(TNamed &h, const char *titleX, const char *titleY, int subDetId)
void useFitForDMRplots(bool usefit=false)
constexpr unsigned int subDetId[12]
void setDMRHistStyleAndLegend(TH1F *h, DMRPlotInfo &plotinfo, int direction=0, int layer=0)
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::vector< TkOfflineVariables * > sourceList
void setCanvasStyle(TCanvas &canv)
void setOutputDir(std::string dir)
PlotAlignmentValidation(const char *inputFile, std::string fileName="", int lineColor=1, int lineStyle=1)
tuple cout
Definition: gather_cfg.py:121
dbl *** dir
Definition: mlp_gen.cc:35