CMS 3D CMS Logo

PVtrends.cc
Go to the documentation of this file.
1 #include <cstdlib>
2 #include <string>
3 #include <tuple>
4 #include <iostream>
5 #include <numeric>
6 #include <functional>
7 #include <unistd.h>
8 
9 #include "TFile.h"
10 #include "TGraph.h"
11 #include "TH1.h"
12 
13 #include "exceptions.h"
14 #include "toolbox.h"
15 #include "Options.h"
16 
18 #include "boost/filesystem.hpp"
19 #include "boost/algorithm/string.hpp"
20 #include "boost/property_tree/ptree.hpp"
21 #include "boost/property_tree/json_parser.hpp"
22 #include "boost/optional.hpp"
23 #include "fmt/printf.h"
24 
25 #include "TString.h"
26 #include "TColor.h"
27 
30 
31 using namespace std;
32 using namespace AllInOneConfig;
33 namespace fs = boost::filesystem;
34 namespace bc = boost::container;
35 
36 static const char *bold = "\e[1m", *normal = "\e[0m";
37 static const float defaultConvertScale = 1000.;
38 static const int startRun2016 = 272930;
39 static const int endRun2018 = 325175;
40 namespace pt = boost::property_tree;
41 
42 int trends(int argc, char *argv[]) {
43  // parse the command line
44 
46  options.helper(argc, argv);
47  options.parser(argc, argv);
48 
49  //Read in AllInOne json config
50  pt::ptree main_tree;
51  pt::read_json(options.config, main_tree);
52 
53  pt::ptree alignments = main_tree.get_child("alignments");
54  pt::ptree validation = main_tree.get_child("validation");
55  pt::ptree style = main_tree.get_child("style");
56 
57  //Read all configure variables and set default for missing keys
58  string outputdir = main_tree.get<string>("output");
59  bool doRMS = validation.count("doRMS") ? validation.get<bool>("doRMS") : true;
60  bool doUnitTest = validation.count("doUnitTest") ? validation.get<bool>("doUnitTest") : false;
61  int nWorkers = validation.count("nWorkers") ? validation.get<int>("nWorkers") : 20;
62  TString lumiInputFile = style.get_child("trends").count("lumiInputFile")
63  ? style.get_child("trends").get<string>("lumiInputFile")
64  : "Alignment/OfflineValidation/data/lumiPerRun_Run2.txt";
65 
66  fs::path lumiFile = lumiInputFile.Data();
67  edm::FileInPath fip = edm::FileInPath(lumiFile.string());
68  fs::path pathToLumiFile = "";
69  if (!fs::exists(lumiFile)) {
70  pathToLumiFile = fip.fullPath();
71  } else {
72  pathToLumiFile = lumiFile;
73  }
74  if (!fs::exists(pathToLumiFile)) {
75  cout << "ERROR: lumi-per-run file (" << lumiFile.string().data() << ") not found!" << endl
76  << "Please check!" << endl;
77  exit(EXIT_FAILURE);
78  } else {
79  cout << "Found lumi-per-run file: " << pathToLumiFile.string().data() << endl;
80  }
81 
82  string lumiAxisType = "recorded";
83  if (lumiInputFile.Contains("delivered"))
84  lumiAxisType = "delivered";
85 
86  std::cout << Form("NOTE: using %s luminosity!", lumiAxisType.data()) << std::endl;
87 
88  for (auto const &childTreeAlignments : alignments) {
89  fs::create_directory(childTreeAlignments.first.c_str());
90  for (auto const &childTreeIOV : validation.get_child("IOV")) {
91  int iov = childTreeIOV.second.get_value<int>();
92  TString file = childTreeAlignments.second.get<string>("file");
93  fs::path input_dir = Form("%s/PVValidation_%s_%d.root",
94  file.ReplaceAll("{}", to_string(iov)).Data(),
95  childTreeAlignments.first.c_str(),
96  iov);
97  fs::path link_dir = Form(
98  "./%s/PVValidation_%s_%d.root", childTreeAlignments.first.c_str(), childTreeAlignments.first.c_str(), iov);
99  if (!fs::exists(link_dir) && fs::exists(input_dir))
100  fs::create_symlink(input_dir, link_dir);
101  }
102  }
103 
104  string labels_to_add = "";
105  if (validation.count("labels")) {
106  for (auto const &label : validation.get_child("labels")) {
107  labels_to_add += "_";
108  labels_to_add += label.second.get_value<string>();
109  }
110  }
111 
112  fs::path pname = Form("%s/PVtrends%s.root", outputdir.data(), labels_to_add.data());
113 
114  PreparePVTrends prepareTrends(pname.c_str(), nWorkers, alignments);
115  prepareTrends.multiRunPVValidation(doRMS, pathToLumiFile.string().data(), doUnitTest);
116 
117  assert(fs::exists(pname));
118 
119  float convertUnit = style.get_child("trends").count("convertUnit")
120  ? style.get_child("trends").get<float>("convertUnit")
122  int firstRun = validation.count("firstRun") ? validation.get<int>("firstRun") : startRun2016;
123  int lastRun = validation.count("lastRun") ? validation.get<int>("lastRun") : endRun2018;
124 
125  const Run2Lumi GetLumi(pathToLumiFile.string().data(), firstRun, lastRun, convertUnit);
126 
127  vector<string> variables{"dxy_phi_vs_run", "dxy_eta_vs_run", "dz_phi_vs_run", "dz_eta_vs_run"};
128 
129  vector<string> titles{"of impact parameter in transverse plane as a function of azimuth angle",
130  "of impact parameter in transverse plane as a function of pseudorapidity",
131  "of impact parameter along z-axis as a function of azimuth angle",
132  "of impact parameter along z-axis as a function of pseudorapidity"};
133 
134  vector<string> ytitles{
135  "of d_{xy}(#phi) [#mum]", "of d_{xy}(#eta) [#mum]", "of d_{z}(#phi) [#mum]", "of d_{z}(#eta) [#mum]"};
136 
137  auto f = TFile::Open(pname.c_str());
138  for (size_t i = 0; i < variables.size(); i++) {
139  Trend mean(Form("mean_%s", variables[i].data()),
140  outputdir.data(),
141  Form("mean %s", titles[i].data()),
142  Form("mean %s", ytitles[i].data()),
143  -7.,
144  10.,
145  style,
146  GetLumi,
147  lumiAxisType.data());
148  Trend RMS(Form("RMS_%s", variables[i].data()),
149  outputdir.data(),
150  Form("RMS %s", titles[i].data()),
151  Form("RMS %s", ytitles[i].data()),
152  0.,
153  35.,
154  style,
155  GetLumi,
156  lumiAxisType.data());
157  mean.lgd.SetHeader("p_{T} (track) > 3 GeV");
158  RMS.lgd.SetHeader("p_{T} (track) > 3 GeV");
159 
160  for (auto const &alignment : alignments) {
161  bool fullRange = true;
162  if (style.get_child("trends").count("earlyStops")) {
163  for (auto const &earlyStop : style.get_child("trends.earlyStops")) {
164  if (earlyStop.second.get_value<string>() == alignment.first)
165  fullRange = false;
166  }
167  }
168 
169  TString gname = alignment.second.get<string>("title");
170  gname.ReplaceAll(" ", "_");
171 
172  auto gMean = Get<TGraph>(fmt::sprintf("mean_%s_%s", gname.Data(), variables[i].data()).c_str());
173  auto hRMS = Get<TH1>(fmt::sprintf("RMS_%s_%s", gname.Data(), variables[i].data()).c_str());
174  assert(gMean != nullptr);
175  assert(hRMS != nullptr);
176 
177  TString gtitle = alignment.second.get<string>("title");
178  gMean->SetTitle(gtitle); // for the legend
179  //gMean->SetTitle(""); // for the legend
180  gMean->SetMarkerSize(0.6);
181  int color = alignment.second.get<int>("color");
182  int style = floor(alignment.second.get<double>("style") / 100.);
183  gMean->SetMarkerColor(color);
184  gMean->SetMarkerStyle(style);
185 
186  hRMS->SetTitle(gtitle); // for the legend
187  //hRMS ->SetTitle(""); // for the legend
188  hRMS->SetMarkerSize(0.6);
189  hRMS->SetMarkerColor(color);
190  hRMS->SetMarkerStyle(style);
191 
192  mean(gMean, "PZ", "p", fullRange);
193  RMS(hRMS, "", "p", fullRange);
194  }
195  }
196 
197  f->Close();
198  cout << bold << "Done" << normal << endl;
199 
200  return EXIT_SUCCESS;
201 }
202 
203 #ifndef DOXYGEN_SHOULD_SKIP_THIS
204 int main(int argc, char *argv[]) { return exceptions<trends>(argc, argv); }
205 #endif
static const int endRun2018
Definition: PVtrends.cc:39
static const std::array< std::string, 5 > titles
std::string fullPath() const
Definition: FileInPath.cc:161
int main(int argc, char *argv[])
Definition: PVtrends.cc:204
void multiRunPVValidation(bool useRMS=true, TString lumiInputFile="", bool doUnitTest=false)
Definition: Trend.h:22
static const float defaultConvertScale
Definition: PVtrends.cc:37
assert(be >=bs)
int trends(int argc, char *argv[])
Definition: PVtrends.cc:42
static std::string to_string(const XMLCh *ch)
char const * label
Definition: Trend.h:78
static const char * normal
Definition: PVtrends.cc:36
Definition: style.py:1
double f[11][100]
static const char * bold
Definition: PVtrends.cc:36
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
firstRun
Definition: dataset.py:940
static const int startRun2016
Definition: PVtrends.cc:38
def exit(msg="")