Counts

PlotsVsDiLeptonRegion

PlotsVsKinematics

etaRegion

flavour

xAxis

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351
#ifndef Alignment_OfflineValidation_DiLeptonVertexHelpers_h
#define Alignment_OfflineValidation_DiLeptonVertexHelpers_h

#include <vector>
#include <string>
#include <fmt/printf.h>
#include "TH2F.h"
#include "TLorentzVector.h"
#include "CommonTools/UtilAlgos/interface/TFileService.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"

namespace DiLeptonHelp {

  //
  // Ancillary struct for counting
  //
  struct Counts {
    unsigned int eventsTotal;
    unsigned int eventsAfterMult;
    unsigned int eventsAfterPt;
    unsigned int eventsAfterEta;
    unsigned int eventsAfterVtx;
    unsigned int eventsAfterDist;

  public:
    void printCounts() {
      edm::LogInfo("DiLeptonHelpCounts") << " Total Events: " << eventsTotal << "\n"
                                         << " After multiplicity: " << eventsAfterMult << "\n"
                                         << " After pT cut: " << eventsAfterPt << "\n"
                                         << " After eta cut: " << eventsAfterEta << "\n"
                                         << " After Vtx: " << eventsAfterVtx << "\n"
                                         << " After VtxDist: " << eventsAfterDist << std::endl;
    }

    void zeroAll() {
      eventsTotal = 0;
      eventsAfterMult = 0;
      eventsAfterPt = 0;
      eventsAfterEta = 0;
      eventsAfterVtx = 0;
      eventsAfterDist = 0;
    }
  };

  enum flavour { MM = 0, EE = 1, UNDEF = -1 };

  enum etaRegion { BARBAR, BARFWD, BARBWD, FWDFWD, BWDBWD, FWDBWD, END };

  //
  // Ancillary class for plotting in different kinematics regions
  // of the two muon tracks
  //
  class PlotsVsDiLeptonRegion {
  public:
    PlotsVsDiLeptonRegion(const float etaBoundary) : m_etaBoundary(etaBoundary) {}
    ~PlotsVsDiLeptonRegion() = default;

    //________________________________________________________________________________//
    inline void bookSet(const TFileDirectory& fs, const TH1* histo) {
      const std::string name = histo->GetName();
      const std::string title = histo->GetTitle();
      const std::string xTitle = histo->GetXaxis()->GetTitle();
      const std::string yTitle = histo->GetYaxis()->GetTitle();
      std::string zTitle = "";
      if (((TObject*)histo)->InheritsFrom("TH2")) {
        zTitle = histo->GetZaxis()->GetTitle();
      }

      for (const auto& etaReg : m_etaRegions) {
        if (etaReg == etaRegion::END)
          continue;

        const auto& toSub = m_etaRegionNames[etaReg];

        if (((TObject*)histo)->InheritsFrom("TH2")) {
          m_h2_map[etaReg] = fs.make<TH2F>(
              (name + "_" + toSub).c_str(),
              (title + " (" + toSub + ");" + xTitle + " (" + toSub + ") ;" + yTitle + " (" + toSub + ");" + zTitle)
                  .c_str(),
              histo->GetNbinsX(),
              histo->GetXaxis()->GetXmin(),
              histo->GetXaxis()->GetXmax(),
              histo->GetNbinsY(),
              histo->GetYaxis()->GetXmin(),
              histo->GetYaxis()->GetXmax());
        } else {
          m_h1_map[etaReg] = fs.make<TH1F>(
              (name + "_" + toSub).c_str(),
              (title + " (" + toSub + ");" + xTitle + " (" + toSub + ") ;" + yTitle + " (" + toSub + ")").c_str(),
              histo->GetNbinsX(),
              histo->GetXaxis()->GetXmin(),
              histo->GetXaxis()->GetXmax());
        }
      }

      // flip the is booked bit
      m_isBooked = true;
    }

    //________________________________________________________________________________//
    // Determine the eta region based on eta values
    etaRegion getEtaRegion(const double eta1, const double eta2) {
      bool isEta1Barrel = std::abs(eta1) <= m_etaBoundary;
      bool isEta2Barrel = std::abs(eta2) <= m_etaBoundary;

      if (isEta1Barrel && isEta2Barrel) {
        return etaRegion::BARBAR;
      } else if ((isEta1Barrel && eta2 > m_etaBoundary) || (isEta2Barrel && eta1 > m_etaBoundary)) {
        return etaRegion::BARFWD;
      } else if ((isEta1Barrel && eta2 < -m_etaBoundary) || (isEta2Barrel && eta1 < -m_etaBoundary)) {
        return etaRegion::BARBWD;
      } else if (eta1 > m_etaBoundary && eta2 > m_etaBoundary) {
        return etaRegion::FWDFWD;
      } else if (eta1 < -m_etaBoundary && eta2 < -m_etaBoundary) {
        return etaRegion::BWDBWD;
      } else if ((eta1 > m_etaBoundary && eta2 < -m_etaBoundary) || (eta2 > m_etaBoundary && eta1 < -m_etaBoundary)) {
        return etaRegion::FWDBWD;
      }

      // Default case if none of the conditions match
      return etaRegion::END;  // Adjust the default based on your logic
    }

    //________________________________________________________________________________//
    inline void fillTH1Plots(const float val, const std::pair<TLorentzVector, TLorentzVector>& momenta) {
      if (!m_isBooked) {
        edm::LogError("PlotsVsDiLeptonRegion")
            << "In" << __FUNCTION__ << "," << __LINE__ << "trying to fill a plot not booked!" << std::endl;
        return;
      }

      etaRegion region = getEtaRegion(momenta.first.Eta(), momenta.second.Eta());
      if (region == etaRegion::END) {
        edm::LogError("PlotsVsDiLeptonRegion") << "undefined di-muon kinematics" << std::endl;
      }
      m_h1_map[region]->Fill(val);
    }

    //________________________________________________________________________________//
    inline void fillTH2Plots(const float valX,
                             const float valY,
                             const std::pair<TLorentzVector, TLorentzVector>& momenta) {
      if (!m_isBooked) {
        edm::LogError("PlotsVsDiLeptonRegion")
            << "In" << __FUNCTION__ << "," << __LINE__ << "trying to fill a plot not booked!" << std::endl;
        return;
      }

      etaRegion region = getEtaRegion(momenta.first.Eta(), momenta.second.Eta());
      if (region == etaRegion::END) {
        edm::LogError("PlotsVsDiLeptonRegion") << "undefined di-muon kinematics" << std::endl;
      }
      m_h2_map[region]->Fill(valX, valY);
    }

  private:
    const std::vector<etaRegion> m_etaRegions = {etaRegion::BARBAR,
                                                 etaRegion::BARFWD,
                                                 etaRegion::BARBWD,
                                                 etaRegion::FWDFWD,
                                                 etaRegion::BWDBWD,
                                                 etaRegion::FWDBWD};

    const std::vector<std::string> m_etaRegionNames = {"barrel-barrel",
                                                       "barrel-forward",
                                                       "barrel-backward",
                                                       "forward-forward",
                                                       "backward-backward",
                                                       "forward-backward"};
    const float m_etaBoundary;
    bool m_isBooked;
    std::map<etaRegion, TH1F*> m_h1_map;
    std::map<etaRegion, TH2F*> m_h2_map;
  };

  //
  // Ancillary class for plotting
  //
  class PlotsVsKinematics {
  public:
    PlotsVsKinematics(flavour FLAV) : m_name(""), m_title(""), m_ytitle(""), m_isBooked(false), m_flav(FLAV) {}

    //________________________________________________________________________________//
    // overloaded constructor
    PlotsVsKinematics(flavour FLAV, const std::string& name, const std::string& tt, const std::string& ytt)
        : m_name(name), m_title(tt), m_ytitle(ytt), m_isBooked(false), m_flav(FLAV) {
      if (m_flav < 0) {
        edm::LogError("PlotsVsKinematics") << "The initialization flavour is not correct!" << std::endl;
      }
    }

    ~PlotsVsKinematics() = default;

    //________________________________________________________________________________//
    inline void bookFromPSet(const TFileDirectory& fs, const edm::ParameterSet& hpar) {
      std::string namePostfix;
      std::string titlePostfix;
      float xmin, xmax;

      std::string sed = (m_flav ? "e" : "#mu");

      for (const auto& xAx : axisChoices) {
        switch (xAx) {
          case xAxis::Z_PHI:
            xmin = -M_PI;
            xmax = M_PI;
            namePostfix = m_flav ? "EEPhi" : "MMPhi";
            titlePostfix = fmt::sprintf("%s%s pair #phi;%s^{+}%s^{-} #phi", sed, sed, sed, sed);
            break;
          case xAxis::Z_ETA:
            xmin = -3.5;
            xmax = 3.5;
            namePostfix = m_flav ? "EEEta" : "MuMuEta";
            titlePostfix = fmt::sprintf("%s%s pair #eta;%s^{+}%s^{-} #eta", sed, sed, sed, sed);
            break;
          case xAxis::LP_PHI:
            xmin = -M_PI;
            xmax = M_PI;
            namePostfix = m_flav ? "EPlusPhi" : "MuPlusPhi";
            titlePostfix = fmt::sprintf("%s^{+} #phi;%s^{+} #phi [rad]", sed, sed);
            break;
          case xAxis::LP_ETA:
            xmin = -2.4;
            xmax = 2.4;
            namePostfix = m_flav ? "EPlusEta" : "MuPlusEta";
            titlePostfix = fmt::sprintf("%s^{+} #eta;%s^{+} #eta", sed, sed);
            break;
          case xAxis::LM_PHI:
            xmin = -M_PI;
            xmax = M_PI;
            namePostfix = m_flav ? "EMinusPhi" : "MuMinusPhi";
            titlePostfix = fmt::sprintf("%s^{-} #phi;%s^{-} #phi [rad]", sed, sed);
            break;
          case xAxis::LM_ETA:
            xmin = -2.4;
            xmax = 2.4;
            namePostfix = m_flav ? "EMinusEta" : "MuMinusEta";
            titlePostfix = fmt::sprintf("%s^{-} #eta;%s^{+} #eta", sed, sed);
            break;
          default:
            throw cms::Exception("LogicalError") << " there is not such Axis choice as " << xAx;
        }

        const auto& h2name = fmt::sprintf("%sVs%s", hpar.getParameter<std::string>("name"), namePostfix);
        const auto& h2title = fmt::sprintf("%s vs %s;%s% s",
                                           hpar.getParameter<std::string>("title"),
                                           titlePostfix,
                                           hpar.getParameter<std::string>("title"),
                                           hpar.getParameter<std::string>("yUnits"));

        m_h2_map[xAx] = fs.make<TH2F>(h2name.c_str(),
                                      h2title.c_str(),
                                      hpar.getParameter<int32_t>("NxBins"),
                                      xmin,
                                      xmax,
                                      hpar.getParameter<int32_t>("NyBins"),
                                      hpar.getParameter<double>("ymin"),
                                      hpar.getParameter<double>("ymax"));
      }

      // flip the is booked bit
      m_isBooked = true;
    }

    //________________________________________________________________________________//
    inline void bookPlots(
        TFileDirectory& fs, const float valmin, const float valmax, const int nxbins, const int nybins) {
      if (m_name.empty() && m_title.empty() && m_ytitle.empty()) {
        edm::LogError("PlotsVsKinematics")
            << "In" << __FUNCTION__ << "," << __LINE__
            << "trying to book plots without the right constructor being called!" << std::endl;
        return;
      }

      std::string dilep = (m_flav ? "e^{+}e^{-}" : "#mu^{+}#mu^{-}");
      std::string lep = (m_flav ? "e^{+}" : "#mu^{+}");
      std::string lem = (m_flav ? "e^{-}" : "#mu^{-}");

      static constexpr float maxMuEta = 2.4;
      static constexpr float maxMuMuEta = 3.5;
      TH1F::SetDefaultSumw2(kTRUE);

      // clang-format off
      m_h2_map[xAxis::Z_ETA] = fs.make<TH2F>(fmt::sprintf("%sVsMuMuEta", m_name).c_str(),
					     fmt::sprintf("%s vs %s pair #eta;%s #eta;%s", m_title, dilep, dilep, m_ytitle).c_str(),
					     nxbins, -M_PI, M_PI,
					     nybins, valmin, valmax);
      
      m_h2_map[xAxis::Z_PHI] = fs.make<TH2F>(fmt::sprintf("%sVsMuMuPhi", m_name).c_str(),
					     fmt::sprintf("%s vs %s pair #phi;%s #phi [rad];%s", m_title, dilep, dilep, m_ytitle).c_str(),
					     nxbins, -maxMuMuEta, maxMuMuEta,
					     nybins, valmin, valmax);
      
      m_h2_map[xAxis::LP_ETA] = fs.make<TH2F>(fmt::sprintf("%sVsMuPlusEta", m_name).c_str(),
					      fmt::sprintf("%s vs %s #eta;%s #eta;%s", m_title, lep, lep, m_ytitle).c_str(),
					      nxbins, -maxMuEta, maxMuEta,
					      nybins, valmin, valmax);
      
      m_h2_map[xAxis::LP_PHI] = fs.make<TH2F>(fmt::sprintf("%sVsMuPlusPhi", m_name).c_str(),
					      fmt::sprintf("%s vs %s #phi;%s #phi [rad];%s", m_title, lep, lep, m_ytitle).c_str(),
					      nxbins, -M_PI, M_PI,
					      nybins, valmin, valmax);
      
      m_h2_map[xAxis::LM_ETA] = fs.make<TH2F>(fmt::sprintf("%sVsMuMinusEta", m_name).c_str(),
					      fmt::sprintf("%s vs %s #eta;%s #eta;%s", m_title, lem, lem, m_ytitle).c_str(),
					      nxbins, -maxMuEta, maxMuEta,
					      nybins, valmin, valmax);
      
      m_h2_map[xAxis::LM_PHI] = fs.make<TH2F>(fmt::sprintf("%sVsMuMinusPhi", m_name).c_str(),
					      fmt::sprintf("%s vs %s #phi;%s #phi [rad];%s", m_title, lem, lem,  m_ytitle).c_str(),
					      nxbins, -M_PI, M_PI,
					      nybins, valmin, valmax);
      // clang-format on

      // flip the is booked bit
      m_isBooked = true;
    }

    //________________________________________________________________________________//
    inline void fillPlots(const float val, const std::pair<TLorentzVector, TLorentzVector>& momenta) {
      if (!m_isBooked) {
        edm::LogError("PlotsVsKinematics")
            << "In" << __FUNCTION__ << "," << __LINE__ << "trying to fill a plot not booked!" << std::endl;
        return;
      }

      m_h2_map[xAxis::Z_ETA]->Fill((momenta.first + momenta.second).Eta(), val);
      m_h2_map[xAxis::Z_PHI]->Fill((momenta.first + momenta.second).Phi(), val);
      m_h2_map[xAxis::LP_ETA]->Fill((momenta.first).Eta(), val);
      m_h2_map[xAxis::LP_PHI]->Fill((momenta.first).Phi(), val);
      m_h2_map[xAxis::LM_ETA]->Fill((momenta.second).Eta(), val);
      m_h2_map[xAxis::LM_PHI]->Fill((momenta.second).Phi(), val);
    }

  private:
    enum xAxis { Z_PHI, Z_ETA, LP_PHI, LP_ETA, LM_PHI, LM_ETA };
    const std::vector<xAxis> axisChoices = {
        xAxis::Z_PHI, xAxis::Z_ETA, xAxis::LP_PHI, xAxis::LP_ETA, xAxis::LM_PHI, xAxis::LM_ETA};

    const std::string m_name;
    const std::string m_title;
    const std::string m_ytitle;

    bool m_isBooked;
    flavour m_flav;

    std::map<xAxis, TH2F*> m_h2_map;
  };
}  // namespace DiLeptonHelp
#endif