File indexing completed on 2024-04-06 12:09:45
0001 #include "DQMOffline/PFTau/interface/MatchCandidateBenchmark.h"
0002
0003 #include "DataFormats/Candidate/interface/Candidate.h"
0004
0005 #include <TFile.h>
0006 #include <TH1.h>
0007 #include <TH2.h>
0008 #include <TProfile.h>
0009 #include <TROOT.h>
0010
0011 using namespace std;
0012
0013 MatchCandidateBenchmark::MatchCandidateBenchmark(Mode mode) : Benchmark(mode) {
0014 delta_et_Over_et_VS_et_ = nullptr;
0015 delta_et_VS_et_ = nullptr;
0016 delta_eta_VS_et_ = nullptr;
0017 delta_phi_VS_et_ = nullptr;
0018
0019 BRdelta_et_Over_et_VS_et_ = nullptr;
0020 ERdelta_et_Over_et_VS_et_ = nullptr;
0021
0022
0023 histogramBooked_ = false;
0024 }
0025
0026 MatchCandidateBenchmark::~MatchCandidateBenchmark() {}
0027
0028 void MatchCandidateBenchmark::setup(DQMStore::IBooker &b) {
0029 if (!histogramBooked_) {
0030 PhaseSpace ptPS;
0031 PhaseSpace dptOvptPS;
0032 PhaseSpace dptPS;
0033 PhaseSpace detaPS;
0034 PhaseSpace dphiPS;
0035 switch (mode_) {
0036 case VALIDATION:
0037 ptPS = PhaseSpace(100, 0, 1000);
0038 dptOvptPS = PhaseSpace(200, -1, 1);
0039 dphiPS = PhaseSpace(200, -1, 1);
0040 detaPS = PhaseSpace(200, -1, 1);
0041 dptPS = PhaseSpace(100, -100, 100);
0042 break;
0043 case DQMOFFLINE:
0044 default:
0045 ptPS = PhaseSpace(50, 0, 100);
0046 dptOvptPS = PhaseSpace(50, -1, 1);
0047 dphiPS = PhaseSpace(50, -1, 1);
0048 detaPS = PhaseSpace(50, -1, 1);
0049 dptPS = PhaseSpace(50, -50, 50);
0050 break;
0051 }
0052 float ptBins[11] = {0, 1, 2, 5, 10, 20, 50, 100, 200, 400, 1000};
0053 int size = sizeof(ptBins) / sizeof(*ptBins);
0054
0055 delta_et_Over_et_VS_et_ = book2D(b,
0056 "delta_et_Over_et_VS_et_",
0057 ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
0058 size,
0059 ptBins,
0060 dptOvptPS.n,
0061 dptOvptPS.m,
0062 dptOvptPS.M);
0063
0064 BRdelta_et_Over_et_VS_et_ = book2D(b,
0065 "BRdelta_et_Over_et_VS_et_",
0066 ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
0067 size,
0068 ptBins,
0069 dptOvptPS.n,
0070 dptOvptPS.m,
0071 dptOvptPS.M);
0072 ERdelta_et_Over_et_VS_et_ = book2D(b,
0073 "ERdelta_et_Over_et_VS_et_",
0074 ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
0075 size,
0076 ptBins,
0077 dptOvptPS.n,
0078 dptOvptPS.m,
0079 dptOvptPS.M);
0080
0081 delta_et_VS_et_ =
0082 book2D(b, "delta_et_VS_et_", ";E_{T, true} (GeV);#DeltaE_{T}", size, ptBins, dptPS.n, dptPS.m, dptPS.M);
0083
0084 delta_eta_VS_et_ =
0085 book2D(b, "delta_eta_VS_et_", ";#E_{T, true} (GeV);#Delta#eta", size, ptBins, detaPS.n, detaPS.m, detaPS.M);
0086
0087 delta_phi_VS_et_ =
0088 book2D(b, "delta_phi_VS_et_", ";E_{T, true} (GeV);#Delta#phi", size, ptBins, dphiPS.n, dphiPS.m, dphiPS.M);
0089 pTRes_.resize(size);
0090 BRpTRes_.resize(size);
0091 ERpTRes_.resize(size);
0092 for (size_t i = 0; i < pTRes_.size(); i++) {
0093 pTRes_[i] = nullptr;
0094 BRpTRes_[i] = nullptr;
0095 ERpTRes_[i] = nullptr;
0096 }
0097
0098 histogramBooked_ = true;
0099 }
0100 }
0101
0102 void MatchCandidateBenchmark::computePtBins(const edm::ParameterSet &ps, const edm::ParameterSet &ptPS) {
0103 const std::vector<double> &ptBinsPS = ps.getParameter<std::vector<double> >("VariablePtBins");
0104 if (ptBinsPS.size() > 1) {
0105 ptBins_.reserve(ptBinsPS.size());
0106 for (size_t i = 0; i < ptBinsPS.size(); i++)
0107 ptBins_.push_back(ptBinsPS[i]);
0108 } else {
0109 Int_t nFixedBins = ptPS.getParameter<int32_t>("nBin");
0110 ptBins_.reserve(nFixedBins + 1);
0111 for (Int_t i = 0; i <= nFixedBins; i++)
0112 ptBins_.push_back(ptPS.getParameter<double>("xMin") +
0113 i * ((ptPS.getParameter<double>("xMax") - ptPS.getParameter<double>("xMin")) / nFixedBins));
0114 }
0115 }
0116 void MatchCandidateBenchmark::setup(DQMStore::IBooker &b, const edm::ParameterSet ¶meterSet) {
0117 if (!histogramBooked_) {
0118 edm::ParameterSet ptPS = parameterSet.getParameter<edm::ParameterSet>("PtHistoParameter");
0119 edm::ParameterSet dptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtHistoParameter");
0120 edm::ParameterSet dptOvptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtOvPtHistoParameter");
0121 edm::ParameterSet detaPS = parameterSet.getParameter<edm::ParameterSet>("DeltaEtaHistoParameter");
0122 edm::ParameterSet dphiPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPhiHistoParameter");
0123 computePtBins(parameterSet, ptPS);
0124 pTRes_.resize(ptBins_.size() - 1);
0125 BRpTRes_.resize(ptBins_.size() - 1);
0126 ERpTRes_.resize(ptBins_.size() - 1);
0127 if (!pTRes_.empty()) {
0128 for (size_t i = 0; i < pTRes_.size(); i++) {
0129 pTRes_[i] = nullptr;
0130 BRpTRes_[i] = nullptr;
0131 ERpTRes_[i] = nullptr;
0132 }
0133 }
0134
0135 if (dptOvptPS.getParameter<bool>("switchOn")) {
0136 delta_et_Over_et_VS_et_ = book2D(b,
0137 "delta_et_Over_et_VS_et_",
0138 ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
0139 ptBins_.size() - 1,
0140 &(ptBins_.front()),
0141 dptOvptPS.getParameter<int32_t>("nBin"),
0142 dptOvptPS.getParameter<double>("xMin"),
0143 dptOvptPS.getParameter<double>("xMax"));
0144 }
0145 if (dptOvptPS.getParameter<bool>("slicingOn")) {
0146 for (size_t i = 0; i < pTRes_.size(); i++) {
0147 pTRes_[i] = book1D(b,
0148 TString::Format("Pt%d_%d", (int)ptBins_[i], (int)ptBins_[i + 1]),
0149 ";#Deltap_{T}/p_{T};Entries",
0150 dptOvptPS.getParameter<int32_t>("nBin"),
0151 dptOvptPS.getParameter<double>("xMin"),
0152 dptOvptPS.getParameter<double>("xMax"));
0153 BRpTRes_[i] = book1D(b,
0154 TString::Format("BRPt%d_%d", (int)ptBins_[i], (int)ptBins_[i + 1]),
0155 ";#Deltap_{T}/p_{T};Entries",
0156 dptOvptPS.getParameter<int32_t>("nBin"),
0157 dptOvptPS.getParameter<double>("xMin"),
0158 dptOvptPS.getParameter<double>("xMax"));
0159 ERpTRes_[i] = book1D(b,
0160 TString::Format("ERPt%d_%d", (int)ptBins_[i], (int)ptBins_[i + 1]),
0161 ";#Deltap_{T}/p_{T};Entries",
0162 dptOvptPS.getParameter<int32_t>("nBin"),
0163 dptOvptPS.getParameter<double>("xMin"),
0164 dptOvptPS.getParameter<double>("xMax"));
0165 }
0166 }
0167 if (dptOvptPS.getParameter<bool>("BROn")) {
0168 BRdelta_et_Over_et_VS_et_ = book2D(b,
0169 "BRdelta_et_Over_et_VS_et_",
0170 ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
0171 ptBins_.size() - 1,
0172 &(ptBins_.front()),
0173 dptOvptPS.getParameter<int32_t>("nBin"),
0174 dptOvptPS.getParameter<double>("xMin"),
0175 dptOvptPS.getParameter<double>("xMax"));
0176 }
0177 if (dptOvptPS.getParameter<bool>("EROn")) {
0178 ERdelta_et_Over_et_VS_et_ = book2D(b,
0179 "ERdelta_et_Over_et_VS_et_",
0180 ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
0181 ptBins_.size() - 1,
0182 &(ptBins_.front()),
0183 dptOvptPS.getParameter<int32_t>("nBin"),
0184 dptOvptPS.getParameter<double>("xMin"),
0185 dptOvptPS.getParameter<double>("xMax"));
0186 }
0187
0188 if (dptPS.getParameter<bool>("switchOn")) {
0189 delta_et_VS_et_ = book2D(b,
0190 "delta_et_VS_et_",
0191 ";E_{T, true} (GeV);#DeltaE_{T}",
0192 ptBins_.size() - 1,
0193 &(ptBins_.front()),
0194 dptPS.getParameter<int32_t>("nBin"),
0195 dptPS.getParameter<double>("xMin"),
0196 dptPS.getParameter<double>("xMax"));
0197 }
0198
0199 if (detaPS.getParameter<bool>("switchOn")) {
0200 delta_eta_VS_et_ = book2D(b,
0201 "delta_eta_VS_et_",
0202 ";E_{T, true} (GeV);#Delta#eta",
0203 ptBins_.size() - 1,
0204 &(ptBins_.front()),
0205 detaPS.getParameter<int32_t>("nBin"),
0206 detaPS.getParameter<double>("xMin"),
0207 detaPS.getParameter<double>("xMax"));
0208 }
0209
0210 if (dphiPS.getParameter<bool>("switchOn")) {
0211 delta_phi_VS_et_ = book2D(b,
0212 "delta_phi_VS_et_",
0213 ";E_{T, true} (GeV);#Delta#phi",
0214 ptBins_.size() - 1,
0215 &(ptBins_.front()),
0216 dphiPS.getParameter<int32_t>("nBin"),
0217 dphiPS.getParameter<double>("xMin"),
0218 dphiPS.getParameter<double>("xMax"));
0219 }
0220 eta_min_barrel_ = dptOvptPS.getParameter<double>("BREtaMin");
0221 eta_max_barrel_ = dptOvptPS.getParameter<double>("BREtaMax");
0222 eta_min_endcap_ = dptOvptPS.getParameter<double>("EREtaMin");
0223 eta_max_endcap_ = dptOvptPS.getParameter<double>("EREtaMax");
0224 histogramBooked_ = true;
0225 }
0226 }
0227
0228 void MatchCandidateBenchmark::fillOne(const reco::Candidate &cand, const reco::Candidate &matchedCand) {
0229 if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
0230 return;
0231
0232 if (histogramBooked_) {
0233 if (delta_et_Over_et_VS_et_)
0234 delta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
0235 if (fabs(cand.eta()) <= 1.4)
0236 if (BRdelta_et_Over_et_VS_et_)
0237 BRdelta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
0238 if (fabs(cand.eta()) >= 1.6 && fabs(cand.eta()) <= 2.4)
0239 if (ERdelta_et_Over_et_VS_et_)
0240 ERdelta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
0241 if (delta_et_VS_et_)
0242 delta_et_VS_et_->Fill(matchedCand.pt(), cand.pt() - matchedCand.pt());
0243 if (delta_eta_VS_et_)
0244 delta_eta_VS_et_->Fill(matchedCand.pt(), cand.eta() - matchedCand.eta());
0245 if (delta_phi_VS_et_)
0246 delta_phi_VS_et_->Fill(matchedCand.pt(), cand.phi() - matchedCand.phi());
0247 }
0248 }
0249
0250 bool MatchCandidateBenchmark::inEtaRange(double value, bool inBarrel) {
0251 if (inBarrel) {
0252 return std::abs(value) >= eta_min_barrel_ && std::abs(value) <= eta_max_barrel_;
0253 }
0254 return std::abs(value) >= eta_min_endcap_ && std::abs(value) <= eta_max_endcap_;
0255 }
0256
0257 void MatchCandidateBenchmark::fillOne(const reco::Candidate &cand,
0258 const reco::Candidate &matchedCand,
0259 const edm::ParameterSet ¶meterSet) {
0260 if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
0261 return;
0262
0263 if (histogramBooked_) {
0264 edm::ParameterSet dptOvptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtOvPtHistoParameter");
0265
0266 if (matchedCand.pt() > ptBins_.at(0)) {
0267 if (delta_et_Over_et_VS_et_)
0268 delta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
0269 if (BRdelta_et_Over_et_VS_et_ and inBarrelRange(matchedCand.eta()))
0270 BRdelta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
0271 if (ERdelta_et_Over_et_VS_et_ and inEndcapRange(matchedCand.eta()))
0272 ERdelta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
0273 if (delta_et_VS_et_)
0274 delta_et_VS_et_->Fill(matchedCand.pt(), cand.pt() - matchedCand.pt());
0275 if (delta_eta_VS_et_)
0276 delta_eta_VS_et_->Fill(matchedCand.pt(), cand.eta() - matchedCand.eta());
0277 if (delta_phi_VS_et_)
0278 delta_phi_VS_et_->Fill(matchedCand.pt(), cand.phi() - matchedCand.phi());
0279 }
0280
0281 for (size_t i = 0; i < pTRes_.size(); i++) {
0282 if (matchedCand.pt() >= ptBins_.at(i) && matchedCand.pt() < ptBins_.at(i + 1)) {
0283 if (pTRes_[i])
0284 pTRes_[i]->Fill((cand.pt() - matchedCand.pt()) / matchedCand.pt());
0285 if (BRpTRes_[i])
0286 BRpTRes_[i]->Fill((cand.pt() - matchedCand.pt()) / matchedCand.pt());
0287 if (ERpTRes_[i])
0288 ERpTRes_[i]->Fill((cand.pt() - matchedCand.pt()) / matchedCand.pt());
0289 }
0290 }
0291 }
0292 }