CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DiLeptonVertexHelpers.h
Go to the documentation of this file.
1 #ifndef Alignment_OfflineValidation_DiLeptonVertexHelpers_h
2 #define Alignment_OfflineValidation_DiLeptonVertexHelpers_h
3 
4 #include <vector>
5 #include <string>
6 #include <fmt/printf.h>
7 #include "TH2F.h"
8 #include "TLorentzVector.h"
12 
13 namespace DiLeptonHelp {
14 
15  //
16  // Ancillary struct for counting
17  //
18  struct Counts {
19  unsigned int eventsTotal;
20  unsigned int eventsAfterMult;
21  unsigned int eventsAfterPt;
22  unsigned int eventsAfterEta;
23  unsigned int eventsAfterVtx;
24  unsigned int eventsAfterDist;
25 
26  public:
27  void printCounts() {
28  edm::LogInfo("DiLeptonHelpCounts") << " Total Events: " << eventsTotal << "\n"
29  << " After multiplicity: " << eventsAfterMult << "\n"
30  << " After pT cut: " << eventsAfterPt << "\n"
31  << " After eta cut: " << eventsAfterEta << "\n"
32  << " After Vtx: " << eventsAfterVtx << "\n"
33  << " After VtxDist: " << eventsAfterDist << std::endl;
34  }
35 
36  void zeroAll() {
37  eventsTotal = 0;
38  eventsAfterMult = 0;
39  eventsAfterPt = 0;
40  eventsAfterEta = 0;
41  eventsAfterVtx = 0;
42  eventsAfterDist = 0;
43  }
44  };
45 
46  enum flavour { MM = 0, EE = 1, UNDEF = -1 };
47 
48  //
49  // Ancillary class for plotting
50  //
52  public:
54 
55  //________________________________________________________________________________//
56  // overloaded constructor
58  : m_name(name), m_title(tt), m_ytitle(ytt), m_isBooked(false), m_flav(FLAV) {
59  if (m_flav < 0) {
60  edm::LogError("PlotsVsKinematics") << "The initialization flavour is not correct!" << std::endl;
61  }
62  }
63 
64  ~PlotsVsKinematics() = default;
65 
66  //________________________________________________________________________________//
67  inline void bookFromPSet(const TFileDirectory& fs, const edm::ParameterSet& hpar) {
68  std::string namePostfix;
69  std::string titlePostfix;
70  float xmin, xmax;
71 
72  std::string sed = (m_flav ? "e" : "#mu");
73 
74  for (const auto& xAx : axisChoices) {
75  switch (xAx) {
76  case xAxis::Z_PHI:
77  xmin = -M_PI;
78  xmax = M_PI;
79  namePostfix = m_flav ? "EEPhi" : "MMPhi";
80  titlePostfix = fmt::sprintf("%s%s pair #phi;%s^{+}%s^{-} #phi", sed, sed, sed, sed);
81  break;
82  case xAxis::Z_ETA:
83  xmin = -3.5;
84  xmax = 3.5;
85  namePostfix = m_flav ? "EEEta" : "MuMuEta";
86  titlePostfix = fmt::sprintf("%s%s pair #eta;%s^{+}%s^{-} #eta", sed, sed, sed, sed);
87  break;
88  case xAxis::LP_PHI:
89  xmin = -M_PI;
90  xmax = M_PI;
91  namePostfix = m_flav ? "EPlusPhi" : "MuPlusPhi";
92  titlePostfix = fmt::sprintf("%s^{+} #phi;%s^{+} #phi [rad]", sed, sed);
93  break;
94  case xAxis::LP_ETA:
95  xmin = -2.4;
96  xmax = 2.4;
97  namePostfix = m_flav ? "EPlusEta" : "MuPlusEta";
98  titlePostfix = fmt::sprintf("%s^{+} #eta;%s^{+} #eta", sed, sed);
99  break;
100  case xAxis::LM_PHI:
101  xmin = -M_PI;
102  xmax = M_PI;
103  namePostfix = m_flav ? "EMinusPhi" : "MuMinusPhi";
104  titlePostfix = fmt::sprintf("%s^{-} #phi;%s^{-} #phi [rad]", sed, sed);
105  break;
106  case xAxis::LM_ETA:
107  xmin = -2.4;
108  xmax = 2.4;
109  namePostfix = m_flav ? "EMinusEta" : "MuMinusEta";
110  titlePostfix = fmt::sprintf("%s^{-} #eta;%s^{+} #eta", sed, sed);
111  break;
112  default:
113  throw cms::Exception("LogicalError") << " there is not such Axis choice as " << xAx;
114  }
115 
116  const auto& h2name = fmt::sprintf("%sVs%s", hpar.getParameter<std::string>("name"), namePostfix);
117  const auto& h2title = fmt::sprintf("%s vs %s;%s% s",
118  hpar.getParameter<std::string>("title"),
119  titlePostfix,
120  hpar.getParameter<std::string>("title"),
121  hpar.getParameter<std::string>("yUnits"));
122 
123  m_h2_map[xAx] = fs.make<TH2F>(h2name.c_str(),
124  h2title.c_str(),
125  hpar.getParameter<int32_t>("NxBins"),
126  xmin,
127  xmax,
128  hpar.getParameter<int32_t>("NyBins"),
129  hpar.getParameter<double>("ymin"),
130  hpar.getParameter<double>("ymax"));
131  }
132 
133  // flip the is booked bit
134  m_isBooked = true;
135  }
136 
137  //________________________________________________________________________________//
138  inline void bookPlots(
139  TFileDirectory& fs, const float valmin, const float valmax, const int nxbins, const int nybins) {
140  if (m_name.empty() && m_title.empty() && m_ytitle.empty()) {
141  edm::LogError("PlotsVsKinematics")
142  << "In" << __FUNCTION__ << "," << __LINE__
143  << "trying to book plots without the right constructor being called!" << std::endl;
144  return;
145  }
146 
147  std::string dilep = (m_flav ? "e^{+}e^{-}" : "#mu^{+}#mu^{-}");
148  std::string lep = (m_flav ? "e^{+}" : "#mu^{+}");
149  std::string lem = (m_flav ? "e^{-}" : "#mu^{-}");
150 
151  static constexpr float maxMuEta = 2.4;
152  static constexpr float maxMuMuEta = 3.5;
153  TH1F::SetDefaultSumw2(kTRUE);
154 
155  // clang-format off
156  m_h2_map[xAxis::Z_ETA] = fs.make<TH2F>(fmt::sprintf("%sVsMuMuEta", m_name).c_str(),
157  fmt::sprintf("%s vs %s pair #eta;%s #eta;%s", m_title, dilep, dilep, m_ytitle).c_str(),
158  nxbins, -M_PI, M_PI,
159  nybins, valmin, valmax);
160 
161  m_h2_map[xAxis::Z_PHI] = fs.make<TH2F>(fmt::sprintf("%sVsMuMuPhi", m_name).c_str(),
162  fmt::sprintf("%s vs %s pair #phi;%s #phi [rad];%s", m_title, dilep, dilep, m_ytitle).c_str(),
163  nxbins, -maxMuMuEta, maxMuMuEta,
164  nybins, valmin, valmax);
165 
166  m_h2_map[xAxis::LP_ETA] = fs.make<TH2F>(fmt::sprintf("%sVsMuPlusEta", m_name).c_str(),
167  fmt::sprintf("%s vs %s #eta;%s #eta;%s", m_title, lep, lep, m_ytitle).c_str(),
168  nxbins, -maxMuEta, maxMuEta,
169  nybins, valmin, valmax);
170 
171  m_h2_map[xAxis::LP_PHI] = fs.make<TH2F>(fmt::sprintf("%sVsMuPlusPhi", m_name).c_str(),
172  fmt::sprintf("%s vs %s #phi;%s #phi [rad];%s", m_title, lep, lep, m_ytitle).c_str(),
173  nxbins, -M_PI, M_PI,
174  nybins, valmin, valmax);
175 
176  m_h2_map[xAxis::LM_ETA] = fs.make<TH2F>(fmt::sprintf("%sVsMuMinusEta", m_name).c_str(),
177  fmt::sprintf("%s vs %s #eta;%s #eta;%s", m_title, lem, lem, m_ytitle).c_str(),
178  nxbins, -maxMuEta, maxMuEta,
179  nybins, valmin, valmax);
180 
181  m_h2_map[xAxis::LM_PHI] = fs.make<TH2F>(fmt::sprintf("%sVsMuMinusPhi", m_name).c_str(),
182  fmt::sprintf("%s vs %s #phi;%s #phi [rad];%s", m_title, lem, lem, m_ytitle).c_str(),
183  nxbins, -M_PI, M_PI,
184  nybins, valmin, valmax);
185  // clang-format on
186 
187  // flip the is booked bit
188  m_isBooked = true;
189  }
190 
191  //________________________________________________________________________________//
192  inline void fillPlots(const float val, const std::pair<TLorentzVector, TLorentzVector>& momenta) {
193  if (!m_isBooked) {
194  edm::LogError("PlotsVsKinematics")
195  << "In" << __FUNCTION__ << "," << __LINE__ << "trying to fill a plot not booked!" << std::endl;
196  return;
197  }
198 
199  m_h2_map[xAxis::Z_ETA]->Fill((momenta.first + momenta.second).Eta(), val);
200  m_h2_map[xAxis::Z_PHI]->Fill((momenta.first + momenta.second).Phi(), val);
201  m_h2_map[xAxis::LP_ETA]->Fill((momenta.first).Eta(), val);
202  m_h2_map[xAxis::LP_PHI]->Fill((momenta.first).Phi(), val);
203  m_h2_map[xAxis::LM_ETA]->Fill((momenta.second).Eta(), val);
204  m_h2_map[xAxis::LM_PHI]->Fill((momenta.second).Phi(), val);
205  }
206 
207  private:
209  const std::vector<xAxis> axisChoices = {
210  xAxis::Z_PHI, xAxis::Z_ETA, xAxis::LP_PHI, xAxis::LP_ETA, xAxis::LM_PHI, xAxis::LM_ETA};
211 
215 
218 
219  std::map<xAxis, TH2F*> m_h2_map;
220  };
221 } // namespace DiLeptonHelp
222 #endif
Log< level::Error, false > LogError
std::map< xAxis, TH2F * > m_h2_map
const std::vector< xAxis > axisChoices
PlotsVsKinematics(flavour FLAV, const std::string &name, const std::string &tt, const std::string &ytt)
T * make(const Args &...args) const
make new ROOT object
#define M_PI
Log< level::Info, false > LogInfo
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void bookPlots(TFileDirectory &fs, const float valmin, const float valmax, const int nxbins, const int nybins)
void fillPlots(const float val, const std::pair< TLorentzVector, TLorentzVector > &momenta)
void bookFromPSet(const TFileDirectory &fs, const edm::ParameterSet &hpar)