File indexing completed on 2025-03-02 23:54:07
0001
0002 #include <cmath>
0003 #include <string>
0004
0005
0006 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0007 #include "DQMServices/Core/interface/DQMStore.h"
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013 #include "FWCore/Framework/interface/LuminosityBlock.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018
0019
0020 #include "TF1.h"
0021 #include "TH1F.h"
0022
0023 class ElectronEfficiencyPlotter : public DQMEDHarvester {
0024 public:
0025
0026 ElectronEfficiencyPlotter(const edm::ParameterSet &ps);
0027
0028 ~ElectronEfficiencyPlotter() override = default;
0029
0030 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0031
0032 protected:
0033
0034 void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override;
0035
0036 private:
0037
0038 const int ptBin_;
0039 const double ptMin_;
0040 const double ptMax_;
0041 const std::string ID_;
0042 const std::string theFolder_;
0043 const std::string sourceFolder_;
0044
0045 MonitorElement *h_eff_pt_EB_doubleEG_HLT;
0046 MonitorElement *h_eff_pt_EE_doubleEG_HLT;
0047 MonitorElement *h_eff_pt_EB_singlePhoton_HLT;
0048 MonitorElement *h_eff_pt_EE_singlePhoton_HLT;
0049
0050 void calculateEfficiency(MonitorElement *Numerator, MonitorElement *Denominator, MonitorElement *Efficiency);
0051 };
0052
0053 using namespace edm;
0054 using namespace std;
0055
0056
0057 ElectronEfficiencyPlotter::ElectronEfficiencyPlotter(const edm::ParameterSet &ps)
0058 : ptBin_{ps.getParameter<int>("ptBin")},
0059 ptMin_{ps.getParameter<double>("ptMin")},
0060 ptMax_{ps.getParameter<double>("ptMax")},
0061 ID_{ps.getParameter<string>("sctElectronID")},
0062 theFolder_{ps.getParameter<string>("folder")},
0063 sourceFolder_{ps.getParameter<string>("srcFolder")} {}
0064
0065 void ElectronEfficiencyPlotter::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0066 edm::ParameterSetDescription desc;
0067 desc.add<int>("ptBin", 5);
0068 desc.add<double>("ptMin", 0);
0069 desc.add<double>("ptMax", 100);
0070 desc.add<string>("sctElectronID", {});
0071 desc.add<string>("folder", {});
0072 desc.add<string>("srcFolder", {});
0073 descriptions.addWithDefaultLabel(desc);
0074 }
0075
0076 void ElectronEfficiencyPlotter::dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0077 ibooker.setCurrentFolder(theFolder_);
0078
0079 h_eff_pt_EB_doubleEG_HLT =
0080 ibooker.book1D("Eff_pt_barrel_DSTdoubleEG", "DSTdoubleEG Eff. vs Pt (barrel)", ptBin_, ptMin_, ptMax_);
0081 h_eff_pt_EE_doubleEG_HLT =
0082 ibooker.book1D("Eff_pt_endcap_DSTdoubleEG", "DSTdoubleEG Eff. vs Pt (endcap)", ptBin_, ptMin_, ptMax_);
0083 h_eff_pt_EB_singlePhoton_HLT =
0084 ibooker.book1D("Eff_pt_barrel_DSTsinglePhoton", "DSTsinglePhoton Eff. vs Pt (barrel)", ptBin_, ptMin_, ptMax_);
0085 h_eff_pt_EE_singlePhoton_HLT =
0086 ibooker.book1D("Eff_pt_endcap_DSTsinglePhoton", "DSTsinglePhoton Eff. vs Pt (endcap)", ptBin_, ptMin_, ptMax_);
0087
0088
0089 h_eff_pt_EB_singlePhoton_HLT->setAxisTitle("p_{T} (GeV)", 1);
0090 h_eff_pt_EE_singlePhoton_HLT->setAxisTitle("p_{T} (GeV)", 1);
0091 h_eff_pt_EB_doubleEG_HLT->setAxisTitle("p_{T} (GeV)", 1);
0092 h_eff_pt_EE_doubleEG_HLT->setAxisTitle("p_{T} (GeV)", 1);
0093
0094 MonitorElement *Numerator_pt_barrel_doubleEG_hlt =
0095 igetter.get(sourceFolder_ +
0096 "/resonanceZ_Tag_pat_Probe_sctElectron_passDoubleEG_DST_"
0097 "fireTrigObj_Pt_Barrel");
0098 MonitorElement *Numerator_pt_endcap_doubleEG_hlt =
0099 igetter.get(sourceFolder_ +
0100 "/resonanceZ_Tag_pat_Probe_sctElectron_passDoubleEG_DST_"
0101 "fireTrigObj_Pt_Endcap");
0102 MonitorElement *Numerator_pt_barrel_singlePhoton_hlt =
0103 igetter.get(sourceFolder_ +
0104 "/resonanceZ_Tag_pat_Probe_sctElectron_passSinglePhoton_DST_"
0105 "fireTrigObj_Pt_Barrel");
0106 MonitorElement *Numerator_pt_endcap_singlePhoton_hlt =
0107 igetter.get(sourceFolder_ +
0108 "/resonanceZ_Tag_pat_Probe_sctElectron_passSinglePhoton_DST_"
0109 "fireTrigObj_Pt_Endcap");
0110 MonitorElement *Denominator_pt_barrel =
0111 igetter.get(sourceFolder_ + "/resonanceZ_Tag_pat_Probe_sctElectron_Pt_Barrel");
0112 MonitorElement *Denominator_pt_endcap =
0113 igetter.get(sourceFolder_ + "/resonanceZ_Tag_pat_Probe_sctElectron_Pt_Endcap");
0114
0115 if (Numerator_pt_barrel_doubleEG_hlt && Denominator_pt_barrel)
0116 calculateEfficiency(Numerator_pt_barrel_doubleEG_hlt, Denominator_pt_barrel, h_eff_pt_EB_doubleEG_HLT);
0117 if (Numerator_pt_endcap_doubleEG_hlt && Denominator_pt_endcap)
0118 calculateEfficiency(Numerator_pt_endcap_doubleEG_hlt, Denominator_pt_endcap, h_eff_pt_EE_doubleEG_HLT);
0119 if (Numerator_pt_barrel_singlePhoton_hlt && Denominator_pt_barrel)
0120 calculateEfficiency(Numerator_pt_barrel_singlePhoton_hlt, Denominator_pt_barrel, h_eff_pt_EB_singlePhoton_HLT);
0121 if (Numerator_pt_endcap_singlePhoton_hlt && Denominator_pt_endcap)
0122 calculateEfficiency(Numerator_pt_endcap_singlePhoton_hlt, Denominator_pt_endcap, h_eff_pt_EE_singlePhoton_HLT);
0123 }
0124
0125 void ElectronEfficiencyPlotter::calculateEfficiency(MonitorElement *Numerator,
0126 MonitorElement *Denominator,
0127 MonitorElement *Efficiency) {
0128 TH1F *h_numerator_pt = Numerator->getTH1F();
0129 TH1F *h_denominator_pt = Denominator->getTH1F();
0130 TH1F *h_eff_pt = Efficiency->getTH1F();
0131 if (h_eff_pt->GetSumw2N() == 0)
0132 h_eff_pt->Sumw2();
0133
0134
0135 int nBins = h_eff_pt->GetNbinsX();
0136 double *binEdges = new double[nBins + 1];
0137 for (int i = 0; i <= nBins; i++)
0138 binEdges[i] = h_eff_pt->GetBinLowEdge(i + 1);
0139
0140 TH1F *h_numerator_pt_rebin = (TH1F *)h_numerator_pt->Rebin(nBins, "num_pt_rebinned", binEdges);
0141 TH1F *h_denominator_pt_rebin = (TH1F *)h_denominator_pt->Rebin(nBins, "num_pt_rebinned", binEdges);
0142 h_eff_pt->Divide(h_numerator_pt_rebin, h_denominator_pt_rebin, 1., 1., "B");
0143 }
0144
0145 DEFINE_FWK_MODULE(ElectronEfficiencyPlotter);