File indexing completed on 2024-04-06 12:09:14
0001 #ifndef Alignment_OfflineValidation_DiLeptonVertexHelpers_h
0002 #define Alignment_OfflineValidation_DiLeptonVertexHelpers_h
0003
0004 #include <vector>
0005 #include <string>
0006 #include <fmt/printf.h>
0007 #include "TH2F.h"
0008 #include "TLorentzVector.h"
0009 #include "DQMServices/Core/interface/DQMStore.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012
0013 namespace DiLepPlotHelp {
0014
0015 enum flavour { MM = 0, EE = 1, UNDEF = -1 };
0016
0017
0018
0019
0020 class PlotsVsKinematics {
0021 public:
0022 PlotsVsKinematics(flavour FLAV) : m_name(""), m_title(""), m_ytitle(""), m_isBooked(false), m_flav(FLAV) {}
0023
0024
0025
0026 PlotsVsKinematics(flavour FLAV, const std::string& name, const std::string& tt, const std::string& ytt)
0027 : m_name(name), m_title(tt), m_ytitle(ytt), m_isBooked(false), m_flav(FLAV) {
0028 if (m_flav < 0) {
0029 edm::LogError("PlotsVsKinematics") << "The initialization flavour is not correct!" << std::endl;
0030 }
0031 }
0032
0033 ~PlotsVsKinematics() = default;
0034
0035
0036 inline void bookFromPSet(dqm::reco::DQMStore::IBooker& iBooker, const edm::ParameterSet& hpar) {
0037 std::string namePostfix;
0038 std::string titlePostfix;
0039 float xmin, xmax;
0040
0041 std::string sed = (m_flav ? "e" : "#mu");
0042
0043 for (const auto& xAx : axisChoices) {
0044 switch (xAx) {
0045 case xAxis::Z_PHI:
0046 xmin = -M_PI;
0047 xmax = M_PI;
0048 namePostfix = m_flav ? "EEPhi" : "MuMuPhi";
0049 titlePostfix = fmt::sprintf("%s%s pair #phi;%s^{+}%s^{-} #phi", sed, sed, sed, sed);
0050 break;
0051 case xAxis::Z_ETA:
0052 xmin = -3.5;
0053 xmax = 3.5;
0054 namePostfix = m_flav ? "EEEta" : "MuMuEta";
0055 titlePostfix = fmt::sprintf("%s%s pair #eta;%s^{+}%s^{-} #eta", sed, sed, sed, sed);
0056 break;
0057 case xAxis::LP_PHI:
0058 xmin = -M_PI;
0059 xmax = M_PI;
0060 namePostfix = m_flav ? "EPlusPhi" : "MuPlusPhi";
0061 titlePostfix = fmt::sprintf("%s^{+} #phi;%s^{+} #phi [rad]", sed, sed);
0062 break;
0063 case xAxis::LP_ETA:
0064 xmin = -2.4;
0065 xmax = 2.4;
0066 namePostfix = m_flav ? "EPlusEta" : "MuPlusEta";
0067 titlePostfix = fmt::sprintf("%s^{+} #eta;%s^{+} #eta", sed, sed);
0068 break;
0069 case xAxis::LM_PHI:
0070 xmin = -M_PI;
0071 xmax = M_PI;
0072 namePostfix = m_flav ? "EMinusPhi" : "MuMinusPhi";
0073 titlePostfix = fmt::sprintf("%s^{-} #phi;%s^{-} #phi [rad]", sed, sed);
0074 break;
0075 case xAxis::LM_ETA:
0076 xmin = -2.4;
0077 xmax = 2.4;
0078 namePostfix = m_flav ? "EMinusEta" : "MuMinusEta";
0079 titlePostfix = fmt::sprintf("%s^{-} #eta;%s^{+} #eta", sed, sed);
0080 break;
0081 case xAxis::DELTA_ETA:
0082 xmin = -hpar.getParameter<double>("maxDeltaEta");
0083 xmax = hpar.getParameter<double>("maxDeltaEta");
0084 namePostfix = m_flav ? "EEDeltEta" : "MuMuDeltaEta";
0085 titlePostfix = fmt::sprintf("%s^{-}%s^{+} #Delta#eta;%s^{+}%s^{-} #Delta#eta", sed, sed, sed, sed);
0086 break;
0087 case xAxis::COSTHETACS:
0088 xmin = -1.;
0089 xmax = 1.;
0090 namePostfix = "CosThetaCS";
0091 titlePostfix =
0092 fmt::sprintf("%s^{+}%s^{-} cos(#theta_{CS});%s^{+}%s^{-} cos(#theta_{CS})", sed, sed, sed, sed);
0093 break;
0094 default:
0095 throw cms::Exception("LogicalError") << " there is not such Axis choice as " << xAx;
0096 }
0097
0098 const auto& h2name = fmt::sprintf("%sVs%s", hpar.getParameter<std::string>("name"), namePostfix);
0099 const auto& h2title = fmt::sprintf("%s vs %s;%s% s",
0100 hpar.getParameter<std::string>("title"),
0101 titlePostfix,
0102 hpar.getParameter<std::string>("title"),
0103 hpar.getParameter<std::string>("yUnits"));
0104
0105 m_h2_map[xAx] = iBooker.book2D(h2name.c_str(),
0106 h2title.c_str(),
0107 hpar.getParameter<int32_t>("NxBins"),
0108 xmin,
0109 xmax,
0110 hpar.getParameter<int32_t>("NyBins"),
0111 hpar.getParameter<double>("ymin"),
0112 hpar.getParameter<double>("ymax"));
0113 }
0114
0115
0116 m_isBooked = true;
0117 }
0118
0119
0120 inline void bookPlots(dqm::reco::DQMStore::IBooker& iBooker,
0121 const float valmin,
0122 const float valmax,
0123 const int nxbins,
0124 const int nybins) {
0125 if (m_name.empty() && m_title.empty() && m_ytitle.empty()) {
0126 edm::LogError("PlotsVsKinematics")
0127 << "In" << __FUNCTION__ << "," << __LINE__
0128 << "trying to book plots without the right constructor being called!" << std::endl;
0129 return;
0130 }
0131
0132 std::string dilep = (m_flav ? "e^{+}e^{-}" : "#mu^{+}#mu^{-}");
0133 std::string lep = (m_flav ? "e^{+}" : "#mu^{+}");
0134 std::string lem = (m_flav ? "e^{-}" : "#mu^{-}");
0135
0136 static constexpr float maxMuEta = 2.4;
0137 static constexpr float maxMuMuEta = 3.5;
0138 TH1F::SetDefaultSumw2(kTRUE);
0139
0140
0141 m_h2_map[xAxis::Z_ETA] = iBooker.book2D(fmt::sprintf("%sVsMuMuEta", m_name).c_str(),
0142 fmt::sprintf("%s vs %s pair #eta;%s #eta;%s", m_title, dilep, dilep, m_ytitle).c_str(),
0143 nxbins, -M_PI, M_PI,
0144 nybins, valmin, valmax);
0145
0146 m_h2_map[xAxis::Z_PHI] = iBooker.book2D(fmt::sprintf("%sVsMuMuPhi", m_name).c_str(),
0147 fmt::sprintf("%s vs %s pair #phi;%s #phi [rad];%s", m_title, dilep, dilep, m_ytitle).c_str(),
0148 nxbins, -maxMuMuEta, maxMuMuEta,
0149 nybins, valmin, valmax);
0150
0151 m_h2_map[xAxis::LP_ETA] = iBooker.book2D(fmt::sprintf("%sVsMuPlusEta", m_name).c_str(),
0152 fmt::sprintf("%s vs %s #eta;%s #eta;%s", m_title, lep, lep, m_ytitle).c_str(),
0153 nxbins, -maxMuEta, maxMuEta,
0154 nybins, valmin, valmax);
0155
0156 m_h2_map[xAxis::LP_PHI] = iBooker.book2D(fmt::sprintf("%sVsMuPlusPhi", m_name).c_str(),
0157 fmt::sprintf("%s vs %s #phi;%s #phi [rad];%s", m_title, lep, lep, m_ytitle).c_str(),
0158 nxbins, -M_PI, M_PI,
0159 nybins, valmin, valmax);
0160
0161 m_h2_map[xAxis::LM_ETA] = iBooker.book2D(fmt::sprintf("%sVsMuMinusEta", m_name).c_str(),
0162 fmt::sprintf("%s vs %s #eta;%s #eta;%s", m_title, lem, lem, m_ytitle).c_str(),
0163 nxbins, -maxMuEta, maxMuEta,
0164 nybins, valmin, valmax);
0165
0166 m_h2_map[xAxis::LM_PHI] = iBooker.book2D(fmt::sprintf("%sVsMuMinusPhi", m_name).c_str(),
0167 fmt::sprintf("%s vs %s #phi;%s #phi [rad];%s", m_title, lem, lem, m_ytitle).c_str(),
0168 nxbins, -M_PI, M_PI,
0169 nybins, valmin, valmax);
0170
0171 m_h2_map[xAxis::DELTA_ETA] = iBooker.book2D(fmt::sprintf("%sVsMuMuDeltaEta", m_name).c_str(),
0172 fmt::sprintf("%s vs %s #Delta#eta;%s #Delta#eta;%s", m_title, dilep, dilep, m_ytitle).c_str(),
0173 nxbins, -4., 4.,
0174 nybins, valmin, valmax);
0175
0176 m_h2_map[xAxis::COSTHETACS] = iBooker.book2D(fmt::sprintf("%sVsCosThetaCS", m_name).c_str(),
0177 fmt::sprintf("%s vs %s cos(#theta_{CS});%s cos(#theta_{CS});%s", m_title, dilep, dilep, m_ytitle).c_str(),
0178 nxbins, -1., 1.,
0179 nybins, valmin, valmax);
0180
0181
0182
0183 m_isBooked = true;
0184 }
0185
0186
0187 inline void fillPlots(const float val, const std::pair<TLorentzVector, TLorentzVector>& momenta) {
0188 if (!m_isBooked) {
0189 edm::LogError("PlotsVsKinematics")
0190 << "In" << __FUNCTION__ << "," << __LINE__ << "trying to fill a plot not booked!" << std::endl;
0191 return;
0192 }
0193
0194 m_h2_map[xAxis::Z_ETA]->Fill((momenta.first + momenta.second).Eta(), val);
0195 m_h2_map[xAxis::Z_PHI]->Fill((momenta.first + momenta.second).Phi(), val);
0196 m_h2_map[xAxis::LP_ETA]->Fill((momenta.first).Eta(), val);
0197 m_h2_map[xAxis::LP_PHI]->Fill((momenta.first).Phi(), val);
0198 m_h2_map[xAxis::LM_ETA]->Fill((momenta.second).Eta(), val);
0199 m_h2_map[xAxis::LM_PHI]->Fill((momenta.second).Phi(), val);
0200
0201
0202 double deltaEta = (momenta.first).Eta() - (momenta.second).Eta();
0203
0204 double muplus = 1. / sqrt(2.) * (momenta.first.E() + momenta.first.Z());
0205 double muminus = 1. / sqrt(2.) * (momenta.first.E() - momenta.first.Z());
0206 double mubarplus = 1. / sqrt(2.) * (momenta.second.E() + momenta.second.Z());
0207 double mubarminus = 1. / sqrt(2.) * (momenta.second.E() - momenta.second.Z());
0208
0209 const auto& mother = momenta.first + momenta.second;
0210 double cosThetaCS = 2. / mother.Mag() / sqrt(pow(mother.Mag(), 2) + pow(mother.Pt(), 2)) *
0211 (muplus * mubarminus - muminus * mubarplus);
0212
0213 m_h2_map[xAxis::DELTA_ETA]->Fill(deltaEta, val);
0214 m_h2_map[xAxis::COSTHETACS]->Fill(cosThetaCS, val);
0215 }
0216
0217 private:
0218 enum xAxis { Z_PHI, Z_ETA, LP_PHI, LP_ETA, LM_PHI, LM_ETA, DELTA_ETA, COSTHETACS };
0219 const std::vector<xAxis> axisChoices = {xAxis::Z_PHI,
0220 xAxis::Z_ETA,
0221 xAxis::LP_PHI,
0222 xAxis::LP_ETA,
0223 xAxis::LM_PHI,
0224 xAxis::LM_ETA,
0225 xAxis::DELTA_ETA,
0226 xAxis::COSTHETACS};
0227
0228 const std::string m_name;
0229 const std::string m_title;
0230 const std::string m_ytitle;
0231
0232 bool m_isBooked;
0233 flavour m_flav;
0234
0235 std::map<xAxis, dqm::reco::MonitorElement*> m_h2_map;
0236 };
0237 }
0238 #endif