Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:40:00

0001 #ifndef ALIGNMENT_OFFLINEVALIDATION_PREPAREDMRTRENDS_H_
0002 #define ALIGNMENT_OFFLINEVALIDATION_PREPAREDMRTRENDS_H_
0003 
0004 #include <iostream>
0005 #include <string>
0006 #include <vector>
0007 #include <algorithm>
0008 #include <map>
0009 #include <iomanip>
0010 #include <fstream>
0011 #include <experimental/filesystem>
0012 #include "TPad.h"
0013 #include "TCanvas.h"
0014 #include "TGraph.h"
0015 #include "TGraphErrors.h"
0016 #include "TMultiGraph.h"
0017 #include "TH1.h"
0018 #include "THStack.h"
0019 #include "TROOT.h"
0020 #include "TFile.h"
0021 #include "TLegend.h"
0022 #include "TLegendEntry.h"
0023 #include "TMath.h"
0024 #include "TRegexp.h"
0025 #include "TPaveLabel.h"
0026 #include "TPaveText.h"
0027 #include "TStyle.h"
0028 #include "TLine.h"
0029 #include "boost/property_tree/ptree.hpp"
0030 #include "boost/property_tree/json_parser.hpp"
0031 
0032 /*!
0033  * \def Dummy value in case a DMR would fail for instance
0034  */
0035 #define DUMMY -999.
0036 /*!
0037  * \def Scale factor value to have mean and sigmas expressed in micrometers.
0038  */
0039 #define DMRFactor 10000.
0040 
0041 /*! \struct Point
0042  *  \brief Structure Point
0043  *         Contains parameters of Gaussian fits to DMRs
0044  *  
0045  * @param run:             run number (IOV boundary)
0046  * @param scale:           scale for the measured quantity: cm->μm for DMRs, 1 for normalized residuals
0047  * @param mu:              mu/mean from Gaussian fit to DMR/DrmsNR
0048  * @param sigma:           sigma/standard deviation from Gaussian fit to DMR/DrmsNR
0049  * @param muplus:          mu/mean for the inward pointing modules
0050  * @param muminus:         mu/mean for outward pointing modules
0051  * @param sigmaplus:       sigma/standard for inward pointing modules 
0052  * @param sigmaminus: //!< sigma/standard for outward pointing modules
0053  */
0054 struct Point {
0055   float run, scale, mu, sigma, muplus, muminus, sigmaplus, sigmaminus;
0056 
0057   /*! \fn Point
0058      *  \brief Constructor of structure Point, initialising all members one by one
0059      */
0060   Point(float Run = DUMMY,
0061         float ScaleFactor = DMRFactor,
0062         float y1 = DUMMY,
0063         float y2 = DUMMY,
0064         float y3 = DUMMY,
0065         float y4 = DUMMY,
0066         float y5 = DUMMY,
0067         float y6 = DUMMY)
0068       : run(Run), scale(ScaleFactor), mu(y1), sigma(y2), muplus(y3), muminus(y5), sigmaplus(y4), sigmaminus(y6) {}
0069 
0070   /*! \fn Point
0071      *  \brief Constructor of structure Point, initialising all members from DMRs directly (with split)
0072      */
0073   Point(float Run, float ScaleFactor, TH1 *histo, TH1 *histoplus, TH1 *histominus)
0074       : Point(Run,
0075               ScaleFactor,
0076               histo->GetMean(),
0077               histo->GetMeanError(),
0078               histoplus->GetMean(),
0079               histoplus->GetMeanError(),
0080               histominus->GetMean(),
0081               histominus->GetMeanError()) {}
0082 
0083   /*! \fn Point
0084      *  \brief Constructor of structure Point, initialising all members from DMRs directly (without split)
0085      */
0086   Point(float Run, float ScaleFactor, TH1 *histo) : Point(Run, ScaleFactor, histo->GetMean(), histo->GetMeanError()) {}
0087 
0088   Point &operator=(const Point &p) {
0089     run = p.run;
0090     mu = p.mu;
0091     muplus = p.muplus;
0092     muminus = p.muminus;
0093     sigma = p.sigma;
0094     sigmaplus = p.sigmaplus;
0095     sigmaminus = p.sigmaminus;
0096     return *this;
0097   }
0098 
0099   inline float GetRun() const { return run; }
0100   inline float GetMu() const { return scale * mu; }
0101   inline float GetMuPlus() const { return scale * muplus; }
0102   inline float GetMuMinus() const { return scale * muminus; }
0103   inline float GetSigma() const { return scale * sigma; }
0104   inline float GetSigmaPlus() const { return scale * sigmaplus; }
0105   inline float GetSigmaMinus() const { return scale * sigmaminus; }
0106   inline float GetDeltaMu() const {
0107     if (muplus == DUMMY && muminus == DUMMY)
0108       return DUMMY;
0109     else
0110       return scale * (muplus - muminus);
0111   }
0112   inline float GetSigmaDeltaMu() const {
0113     if (sigmaplus == DUMMY && sigmaminus == DUMMY)
0114       return DUMMY;
0115     else
0116       return scale * hypot(sigmaplus, sigmaminus);
0117   }
0118 };
0119 
0120 /*! \class Geometry
0121  *  \brief Class Geometry
0122  *         Contains vector for fit parameters (mean, sigma, etc.) obtained from multiple IOVs
0123  *         See Structure Point for description of the parameters.
0124  */
0125 
0126 class Geometry {
0127 public:
0128   std::vector<Point> points;
0129 
0130 private:
0131   //template<typename T> std::vector<T> GetQuantity (T (Point::*getter)() const) const {
0132   std::vector<float> GetQuantity(float (Point::*getter)() const) const {
0133     std::vector<float> v;
0134     for (Point point : points) {
0135       float value = (point.*getter)();
0136       v.push_back(value);
0137     }
0138     return v;
0139   }
0140 
0141 public:
0142   TString title;
0143   Geometry() : title("") {}
0144   Geometry(TString Title) : title(Title) {}
0145   Geometry &operator=(const Geometry &geom) {
0146     title = geom.title;
0147     points = geom.points;
0148     return *this;
0149   }
0150   inline void SetTitle(TString Title) { title = Title; }
0151   inline TString GetTitle() { return title; }
0152   inline std::vector<float> Run() const { return GetQuantity(&Point::GetRun); }
0153   inline std::vector<float> Mu() const { return GetQuantity(&Point::GetMu); }
0154   inline std::vector<float> MuPlus() const { return GetQuantity(&Point::GetMuPlus); }
0155   inline std::vector<float> MuMinus() const { return GetQuantity(&Point::GetMuMinus); }
0156   inline std::vector<float> Sigma() const { return GetQuantity(&Point::GetSigma); }
0157   inline std::vector<float> SigmaPlus() const { return GetQuantity(&Point::GetSigmaPlus); }
0158   inline std::vector<float> SigmaMinus() const { return GetQuantity(&Point::GetSigmaMinus); }
0159   inline std::vector<float> DeltaMu() const { return GetQuantity(&Point::GetDeltaMu); }
0160   inline std::vector<float> SigmaDeltaMu() const { return GetQuantity(&Point::GetSigmaDeltaMu); }
0161 };
0162 
0163 class PrepareDMRTrends {
0164 public:
0165   PrepareDMRTrends(const char *outputFileName, boost::property_tree::ptree &json);
0166   ~PrepareDMRTrends() {}
0167 
0168   TString getName(TString structure, int layer, TString geometry);
0169   void compileDMRTrends(std::vector<int> IOVlist,
0170                         TString Variable,
0171                         std::vector<std::string> inputFiles,
0172                         std::vector<TString> structures,
0173                         const std::map<TString, int> nlayers,
0174                         bool FORCE = false);
0175 
0176 private:
0177   const char *outputFileName_;
0178   std::vector<std::string> geometries;
0179 };
0180 
0181 #endif  // ALIGNMENT_OFFLINEVALIDATION_PREPAREDMRTRENDS_H_