Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <cstdlib>
0002 #include <string>
0003 #include <tuple>
0004 #include <iostream>
0005 #include <numeric>
0006 #include <functional>
0007 #include <unistd.h>
0008 
0009 #include "TFile.h"
0010 #include "TGraph.h"
0011 #include "TH1.h"
0012 
0013 #include "exceptions.h"
0014 #include "toolbox.h"
0015 #include "Options.h"
0016 
0017 #include "FWCore/ParameterSet/interface/FileInPath.h"
0018 #include "boost/filesystem.hpp"
0019 #include "boost/algorithm/string.hpp"
0020 #include "boost/property_tree/ptree.hpp"
0021 #include "boost/property_tree/json_parser.hpp"
0022 #include "boost/optional.hpp"
0023 #include "fmt/printf.h"
0024 
0025 #include "TString.h"
0026 #include "TColor.h"
0027 
0028 #include "Alignment/OfflineValidation/interface/PreparePVTrends.h"
0029 #include "Alignment/OfflineValidation/interface/Trend.h"
0030 
0031 using namespace std;
0032 using namespace AllInOneConfig;
0033 namespace fs = boost::filesystem;
0034 namespace bc = boost::container;
0035 
0036 static const char *bold = "\e[1m", *normal = "\e[0m";
0037 static const float defaultConvertScale = 1000.;
0038 static const int startRun2016 = 272930;
0039 static const int endRun2018 = 325175;
0040 namespace pt = boost::property_tree;
0041 
0042 int trends(int argc, char *argv[]) {
0043   // parse the command line
0044 
0045   Options options;
0046   options.helper(argc, argv);
0047   options.parser(argc, argv);
0048 
0049   //Read in AllInOne json config
0050   pt::ptree main_tree;
0051   pt::read_json(options.config, main_tree);
0052 
0053   pt::ptree alignments = main_tree.get_child("alignments");
0054   pt::ptree validation = main_tree.get_child("validation");
0055   pt::ptree style = main_tree.get_child("style");
0056 
0057   //Read all configure variables and set default for missing keys
0058   string outputdir = main_tree.get<string>("output");
0059   bool doRMS = validation.count("doRMS") ? validation.get<bool>("doRMS") : true;
0060   bool doUnitTest = validation.count("doUnitTest") ? validation.get<bool>("doUnitTest") : false;
0061   int nWorkers = validation.count("nWorkers") ? validation.get<int>("nWorkers") : 20;
0062   TString lumiInputFile = style.get_child("trends").count("lumiInputFile")
0063                               ? style.get_child("trends").get<string>("lumiInputFile")
0064                               : "Alignment/OfflineValidation/data/lumiPerRun_Run2.txt";
0065 
0066   fs::path lumiFile = lumiInputFile.Data();
0067   edm::FileInPath fip = edm::FileInPath(lumiFile.string());
0068   fs::path pathToLumiFile = "";
0069   if (!fs::exists(lumiFile)) {
0070     pathToLumiFile = fip.fullPath();
0071   } else {
0072     pathToLumiFile = lumiFile;
0073   }
0074   if (!fs::exists(pathToLumiFile)) {
0075     cout << "ERROR: lumi-per-run file (" << lumiFile.string().data() << ") not found!" << endl
0076          << "Please check!" << endl;
0077     exit(EXIT_FAILURE);
0078   } else {
0079     cout << "Found lumi-per-run file: " << pathToLumiFile.string().data() << endl;
0080   }
0081 
0082   string lumiAxisType = "recorded";
0083   if (lumiInputFile.Contains("delivered"))
0084     lumiAxisType = "delivered";
0085 
0086   std::cout << Form("NOTE: using %s luminosity!", lumiAxisType.data()) << std::endl;
0087 
0088   for (auto const &childTreeAlignments : alignments) {
0089     fs::create_directory(childTreeAlignments.first.c_str());
0090     for (auto const &childTreeIOV : validation.get_child("IOV")) {
0091       int iov = childTreeIOV.second.get_value<int>();
0092       TString file = childTreeAlignments.second.get<string>("file");
0093       fs::path input_dir = Form("%s/PVValidation_%s_%d.root",
0094                                 file.ReplaceAll("{}", to_string(iov)).Data(),
0095                                 childTreeAlignments.first.c_str(),
0096                                 iov);
0097       fs::path link_dir = Form(
0098           "./%s/PVValidation_%s_%d.root", childTreeAlignments.first.c_str(), childTreeAlignments.first.c_str(), iov);
0099       if (!fs::exists(link_dir) && fs::exists(input_dir))
0100         fs::create_symlink(input_dir, link_dir);
0101     }
0102   }
0103 
0104   string labels_to_add = "";
0105   if (validation.count("labels")) {
0106     for (auto const &label : validation.get_child("labels")) {
0107       labels_to_add += "_";
0108       labels_to_add += label.second.get_value<string>();
0109     }
0110   }
0111 
0112   fs::path pname = Form("%s/PVtrends%s.root", outputdir.data(), labels_to_add.data());
0113 
0114   PreparePVTrends prepareTrends(pname.c_str(), nWorkers, alignments);
0115   prepareTrends.multiRunPVValidation(doRMS, pathToLumiFile.string().data(), doUnitTest);
0116 
0117   assert(fs::exists(pname));
0118 
0119   float convertUnit = style.get_child("trends").count("convertUnit")
0120                           ? style.get_child("trends").get<float>("convertUnit")
0121                           : defaultConvertScale;
0122   int firstRun = validation.count("firstRun") ? validation.get<int>("firstRun") : startRun2016;
0123   int lastRun = validation.count("lastRun") ? validation.get<int>("lastRun") : endRun2018;
0124 
0125   const Run2Lumi GetLumi(pathToLumiFile.string().data(), firstRun, lastRun, convertUnit);
0126 
0127   vector<string> variables{"dxy_phi_vs_run", "dxy_eta_vs_run", "dz_phi_vs_run", "dz_eta_vs_run"};
0128 
0129   vector<string> titles{"of impact parameter in transverse plane as a function of azimuth angle",
0130                         "of impact parameter in transverse plane as a function of pseudorapidity",
0131                         "of impact parameter along z-axis as a function of azimuth angle",
0132                         "of impact parameter along z-axis as a function of pseudorapidity"};
0133 
0134   vector<string> ytitles{
0135       "of d_{xy}(#phi)   [#mum]", "of d_{xy}(#eta)   [#mum]", "of d_{z}(#phi)   [#mum]", "of d_{z}(#eta)   [#mum]"};
0136 
0137   auto f = TFile::Open(pname.c_str());
0138   for (size_t i = 0; i < variables.size(); i++) {
0139     Trend mean(Form("mean_%s", variables[i].data()),
0140                outputdir.data(),
0141                Form("mean %s", titles[i].data()),
0142                Form("mean %s", ytitles[i].data()),
0143                -4.,
0144                8.,
0145                style,
0146                GetLumi,
0147                lumiAxisType.data());
0148     Trend RMS(Form("RMS_%s", variables[i].data()),
0149               outputdir.data(),
0150               Form("RMS %s", titles[i].data()),
0151               Form("RMS %s", ytitles[i].data()),
0152               0.,
0153               25.,
0154               style,
0155               GetLumi,
0156               lumiAxisType.data());
0157     mean.lgd.SetHeader("p_{T} (track) > 3 GeV");
0158     RMS.lgd.SetHeader("p_{T} (track) > 3 GeV");
0159 
0160     for (auto const &alignment : alignments) {
0161       bool fullRange = true;
0162       if (style.get_child("trends").count("earlyStops")) {
0163         for (auto const &earlyStop : style.get_child("trends.earlyStops")) {
0164           if (earlyStop.second.get_value<string>() == alignment.first)
0165             fullRange = false;
0166         }
0167       }
0168 
0169       TString gname = alignment.second.get<string>("title");
0170       gname.ReplaceAll(" ", "_");
0171 
0172       auto gMean = Get<TGraph>(fmt::sprintf("mean_%s_%s", gname.Data(), variables[i].data()).c_str());
0173       auto hRMS = Get<TH1>(fmt::sprintf("RMS_%s_%s", gname.Data(), variables[i].data()).c_str());
0174       assert(gMean != nullptr);
0175       assert(hRMS != nullptr);
0176 
0177       TString gtitle = alignment.second.get<string>("title");
0178       gMean->SetTitle(gtitle);  // for the legend
0179       //gMean->SetTitle(""); // for the legend
0180       gMean->SetMarkerSize(0.6);
0181       int color = alignment.second.get<int>("color");
0182       int style = floor(alignment.second.get<double>("style") / 100.);
0183       gMean->SetMarkerColor(color);
0184       gMean->SetLineColor(color);  // no need to be set but looks better IMHO
0185       gMean->SetMarkerStyle(style);
0186       gMean->SetMarkerSize(1.6);
0187 
0188       hRMS->SetTitle(gtitle);  // for the legend
0189       //hRMS ->SetTitle(""); // for the legend
0190       hRMS->SetMarkerSize(0.6);
0191       hRMS->SetMarkerColor(color);
0192       hRMS->SetLineColor(color);  // needs to be set, otherwise color is NOT picked up
0193       hRMS->SetMarkerStyle(style);
0194 
0195       mean(gMean, "PZ", "p", fullRange);
0196       RMS(hRMS, "", "p", fullRange);
0197     }
0198   }
0199 
0200   f->Close();
0201   cout << bold << "Done" << normal << endl;
0202 
0203   return EXIT_SUCCESS;
0204 }
0205 
0206 #ifndef DOXYGEN_SHOULD_SKIP_THIS
0207 int main(int argc, char *argv[]) { return exceptions<trends>(argc, argv); }
0208 #endif