CMS 3D CMS Logo

PrepareDMRTrends.h
Go to the documentation of this file.
1 #ifndef ALIGNMENT_OFFLINEVALIDATION_PREPAREDMRTRENDS_H_
2 #define ALIGNMENT_OFFLINEVALIDATION_PREPAREDMRTRENDS_H_
3 
4 #include <iostream>
5 #include <string>
6 #include <vector>
7 #include <algorithm>
8 #include <map>
9 #include <iomanip>
10 #include <fstream>
11 #include <experimental/filesystem>
12 #include "TPad.h"
13 #include "TCanvas.h"
14 #include "TGraph.h"
15 #include "TGraphErrors.h"
16 #include "TMultiGraph.h"
17 #include "TH1.h"
18 #include "THStack.h"
19 #include "TROOT.h"
20 #include "TFile.h"
21 #include "TLegend.h"
22 #include "TLegendEntry.h"
23 #include "TMath.h"
24 #include "TRegexp.h"
25 #include "TPaveLabel.h"
26 #include "TPaveText.h"
27 #include "TStyle.h"
28 #include "TLine.h"
29 #include "boost/property_tree/ptree.hpp"
30 #include "boost/property_tree/json_parser.hpp"
31 
35 #define DUMMY -999.
36 
39 #define DMRFactor 10000.
40 
54 struct Point {
56 
60  Point(float Run = DUMMY,
61  float ScaleFactor = DMRFactor,
62  float y1 = DUMMY,
63  float y2 = DUMMY,
64  float y3 = DUMMY,
65  float y4 = DUMMY,
66  float y5 = DUMMY,
67  float y6 = DUMMY)
68  : run(Run), scale(ScaleFactor), mu(y1), sigma(y2), muplus(y3), muminus(y5), sigmaplus(y4), sigmaminus(y6) {}
69 
73  Point(float Run, float ScaleFactor, TH1 *histo, TH1 *histoplus, TH1 *histominus)
74  : Point(Run,
75  ScaleFactor,
76  histo->GetMean(),
77  histo->GetMeanError(),
78  histoplus->GetMean(),
79  histoplus->GetMeanError(),
80  histominus->GetMean(),
81  histominus->GetMeanError()) {}
82 
86  Point(float Run, float ScaleFactor, TH1 *histo) : Point(Run, ScaleFactor, histo->GetMean(), histo->GetMeanError()) {}
87 
88  Point &operator=(const Point &p) {
89  run = p.run;
90  mu = p.mu;
91  muplus = p.muplus;
92  muminus = p.muminus;
93  sigma = p.sigma;
94  sigmaplus = p.sigmaplus;
95  sigmaminus = p.sigmaminus;
96  return *this;
97  }
98 
99  inline float GetRun() const { return run; }
100  inline float GetMu() const { return scale * mu; }
101  inline float GetMuPlus() const { return scale * muplus; }
102  inline float GetMuMinus() const { return scale * muminus; }
103  inline float GetSigma() const { return scale * sigma; }
104  inline float GetSigmaPlus() const { return scale * sigmaplus; }
105  inline float GetSigmaMinus() const { return scale * sigmaminus; }
106  inline float GetDeltaMu() const {
107  if (muplus == DUMMY && muminus == DUMMY)
108  return DUMMY;
109  else
110  return scale * (muplus - muminus);
111  }
112  inline float GetSigmaDeltaMu() const {
113  if (sigmaplus == DUMMY && sigmaminus == DUMMY)
114  return DUMMY;
115  else
116  return scale * hypot(sigmaplus, sigmaminus);
117  }
118 };
119 
126 class Geometry {
127 public:
128  std::vector<Point> points;
129 
130 private:
131  //template<typename T> std::vector<T> GetQuantity (T (Point::*getter)() const) const {
132  std::vector<float> GetQuantity(float (Point::*getter)() const) const {
133  std::vector<float> v;
134  for (Point point : points) {
135  float value = (point.*getter)();
136  v.push_back(value);
137  }
138  return v;
139  }
140 
141 public:
142  TString title;
143  Geometry() : title("") {}
144  Geometry(TString Title) : title(Title) {}
146  title = geom.title;
147  points = geom.points;
148  return *this;
149  }
150  inline void SetTitle(TString Title) { title = Title; }
151  inline TString GetTitle() { return title; }
152  inline std::vector<float> Run() const { return GetQuantity(&Point::GetRun); }
153  inline std::vector<float> Mu() const { return GetQuantity(&Point::GetMu); }
154  inline std::vector<float> MuPlus() const { return GetQuantity(&Point::GetMuPlus); }
155  inline std::vector<float> MuMinus() const { return GetQuantity(&Point::GetMuMinus); }
156  inline std::vector<float> Sigma() const { return GetQuantity(&Point::GetSigma); }
157  inline std::vector<float> SigmaPlus() const { return GetQuantity(&Point::GetSigmaPlus); }
158  inline std::vector<float> SigmaMinus() const { return GetQuantity(&Point::GetSigmaMinus); }
159  inline std::vector<float> DeltaMu() const { return GetQuantity(&Point::GetDeltaMu); }
160  inline std::vector<float> SigmaDeltaMu() const { return GetQuantity(&Point::GetSigmaDeltaMu); }
161 };
162 
164 public:
165  PrepareDMRTrends(const char *outputFileName, boost::property_tree::ptree &json);
167 
168  TString getName(TString structure, int layer, TString geometry);
169  void compileDMRTrends(std::vector<int> IOVlist,
170  TString Variable,
171  std::vector<std::string> inputFiles,
172  std::vector<TString> structures,
173  const std::map<TString, int> nlayers,
174  bool FORCE = false);
175 
176 private:
177  const char *outputFileName_;
178  std::vector<std::string> geometries;
179 };
180 
181 #endif // ALIGNMENT_OFFLINEVALIDATION_PREPAREDMRTRENDS_H_
Geometry(TString Title)
std::vector< float > MuPlus() const
float GetDeltaMu() const
float mu
Geometry & operator=(const Geometry &geom)
std::vector< float > DeltaMu() const
void compileDMRTrends(std::vector< int > IOVlist, TString Variable, std::vector< std::string > inputFiles, std::vector< TString > structures, const std::map< TString, int > nlayers, bool FORCE=false)
float GetMuMinus() const
nlohmann::json json
float GetRun() const
float sigma
#define DMRFactor
TString getName(TString structure, int layer, TString geometry)
Point(float Run=-999., float ScaleFactor=10000., float y1=-999., float y2=-999., float y3=-999., float y4=-999., float y5=-999., float y6=-999.)
std::vector< Point > points
TString title
std::vector< std::string > geometries
Point(float Run, float ScaleFactor, TH1 *histo, TH1 *histoplus, TH1 *histominus)
void SetTitle(TString Title)
std::vector< float > MuMinus() const
Point & operator=(const Point &p)
std::vector< float > Sigma() const
PrepareDMRTrends(const char *outputFileName, boost::property_tree::ptree &json)
std::vector< float > SigmaPlus() const
float GetMu() const
Class Geometry Contains vector for fit parameters (mean, sigma, etc.) obtained from multiple IOVs See...
std::vector< float > SigmaDeltaMu() const
float run
std::vector< float > Mu() const
float GetMuPlus() const
float muminus
float GetSigmaMinus() const
Definition: value.py:1
float sigmaplus
std::vector< float > GetQuantity(float(Point::*getter)() const) const
#define DUMMY
std::vector< float > Run() const
Point(float Run, float ScaleFactor, TH1 *histo)
Structure Point Contains parameters of Gaussian fits to DMRs.
float GetSigmaDeltaMu() const
float sigmaminus
const char * outputFileName_
float GetSigmaPlus() const
float GetSigma() const
std::vector< float > SigmaMinus() const
TString GetTitle()
float scale
float muplus
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5