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 % 100;
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);
76  PlotAlignmentValidation(const char *inputFile,std::string fileName="", int lineColor=1, int lineStyle=1, bool bigtext=false);
78  void loadFileList(const char *inputFile, std::string fileName="", int lineColor=2, int lineStyle=1);
79  void useFitForDMRplots(bool usefit = false);
80  void legendOptions(TString options);
81  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
82  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-
83  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")
84  void plotSurfaceShapes(const std::string& options = "layers",const std::string& variable="");
85  void plotChi2(const char *inputFile);
86  // plotSurfaceShapes: options="split","layers"/"layer","subdet"
87  void plotHitMaps();
89  void setTreeBaseDir( std::string dir = "TrackerOfflineValidationStandalone");
90 
91  THStack* addHists(const TString& selection, const TString &residType = "xPrime", TLegend **myLegend = 0, bool printModuleIds = false);//add hists fulfilling 'selection' on TTree; residType: xPrime,yPrime,xPrimeNorm,yPrimeNorm,x,y,xNorm; if (printModuleIds): cout DetIds
92 
93  // These are helpers for DMR plotting
94 
95  struct DMRPlotInfo {
97  int nbins;
98  double min, max;
99  int minHits;
102  THStack* hstack;
103  TLegend* legend;
105  float maxY;
106  TH1F* h;
107  TH1F* h1;
108  TH1F* h2;
110  };
111 
112 private :
113  TList* getTreeList();
115 
116  bool useFit_;
117  bool showMean_;
118  bool showRMS_;
123  bool twolines_;
124  bool bigtext_;
125 
126  TF1 *fitGauss(TH1 *hist,int color);
127  //void plotBoxOverview(TCanvas &c1, TList &treeList,std::string plot_Var1a,std::string plot_Var1b, std::string plot_Var2, Int_t filenumber,Int_t minHits);
128  //void plot1DDetailsSubDet(TCanvas &c1, TList &treeList, std::string plot_Var1a,std::string plot_Var1b, std::string plot_Var2, Int_t minHits);
129  //void plot1DDetailsBarrelLayer(TCanvas &c1, TList &treeList, std::string plot_Var1a,std::string plot_Var1b, Int_t minHits);
130  //void plot1DDetailsDiskWheel(TCanvas &c1, TList &treelist, std::string plot_Var1a,std::string plot_Var1b, Int_t minHits);
131  void plotSS(const std::string& options = "layers",const std::string& variable="");
132  void setHistStyle( TH1& hist,const char* titleX, const char* titleY, int color);
133  void setTitleStyle( TNamed& h,const char* titleX, const char* titleY, int subDetId, bool isSurfaceDeformation=false, TString secondline="");
134  void setNiceStyle();
135  void setCanvasStyle( TCanvas& canv );
136  void setLegendStyle( TLegend& leg );
137  void scaleXaxis(TH1* hist, Int_t scale);
138  TObject* findObjectFromCanvas(TCanvas* canv, const char* className, Int_t n=1);
139 
140  TString outputFile;
142  TList *sourcelist;
143  std::vector<TkOfflineVariables*> sourceList;
147 
148  std::string getSelectionForDMRPlot(int minHits, int subDetId, int direction = 0, int layer = 0);
149  std::string getVariableForDMRPlot(const std::string& histoname, const std::string& variable,
150  int nbins, double min, double max);
151  void setDMRHistStyleAndLegend(TH1F* h, DMRPlotInfo& plotinfo, int direction = 0, int layer = 0);
152  void plotDMRHistogram(DMRPlotInfo& plotinfo, int direction = 0, int layer = 0);
153  void modifySSHistAndLegend(THStack* hs, TLegend* legend);
154 };
155 
156 #endif // PLOTALIGNNMENTVALIDATION_H_
void plotSurfaceShapes(const std::string &options="layers", const std::string &variable="")
void scaleXaxis(TH1 *hist, Int_t scale)
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
void setHistStyle(TH1 &hist, const char *titleX, const char *titleY, int color)
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
TObject * findObjectFromCanvas(TCanvas *canv, const char *className, Int_t n=1)
void legendOptions(TString options)
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)
THStack * addHists(const TString &selection, const TString &residType="xPrime", TLegend **myLegend=0, bool printModuleIds=false)
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 setTitleStyle(TNamed &h, const char *titleX, const char *titleY, int subDetId, bool isSurfaceDeformation=false, TString secondline="")
void setTreeBaseDir(std::string dir="TrackerOfflineValidationStandalone")
void plotChi2(const char *inputFile)
void setLegendStyle(TLegend &leg)
tuple d
Definition: ztail.py:151
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)
TkOfflineVariables(std::string fileName, std::string baseDir, std::string legName="", int color=1, int style=1)
void modifySSHistAndLegend(THStack *hs, TLegend *legend)
void plotDMR(const std::string &plotVar="medianX", Int_t minHits=50, const std::string &options="plain")
unsigned int subDetId[18]
void useFitForDMRplots(bool usefit=false)
T min(T a, T b)
Definition: MathUtil.h:58
void setDMRHistStyleAndLegend(TH1F *h, DMRPlotInfo &plotinfo, int direction=0, int layer=0)
std::vector< TkOfflineVariables * > sourceList
void setCanvasStyle(TCanvas &canv)
void setOutputDir(std::string dir)
tuple cout
Definition: gather_cfg.py:145
dbl *** dir
Definition: mlp_gen.cc:35
TF1 * fitGauss(TH1 *hist, int color)
std::string className(const T &t)
Definition: ClassName.h:30