File indexing completed on 2024-04-06 12:09:44
0001 #ifndef DQMOffline_PFTau_PFCandidateMonitor_h
0002 #define DQMOffline_PFTau_PFCandidateMonitor_h
0003
0004 #include "DQMOffline/PFTau/interface/Benchmark.h"
0005 #include "DQMOffline/PFTau/interface/CandidateBenchmark.h"
0006 #include "DQMOffline/PFTau/interface/MatchCandidateBenchmark.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008
0009 #include "DataFormats/Candidate/interface/Candidate.h"
0010 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0011
0012 #include <vector>
0013
0014 #include <TH1.h> //needed by the deltaR->Fill() call
0015
0016 class PFCandidateMonitor : public Benchmark {
0017 public:
0018 PFCandidateMonitor(float dRMax = 0.3, bool matchCharge = true, Benchmark::Mode mode = Benchmark::DEFAULT);
0019
0020 ~PFCandidateMonitor() override;
0021
0022
0023 void setParameters(float dRMax,
0024 bool matchCharge,
0025 Benchmark::Mode mode,
0026 float ptmin,
0027 float ptmax,
0028 float etamin,
0029 float etamax,
0030 float phimin,
0031 float phimax,
0032 bool refHistoFlag);
0033
0034
0035 void setParameters(const edm::ParameterSet ¶meterSet);
0036
0037
0038 void setDirectory(TDirectory *dir) override;
0039
0040
0041 void setup(DQMStore::IBooker &b);
0042 void setup(DQMStore::IBooker &b, const edm::ParameterSet ¶meterSet);
0043
0044
0045 template <class T, class C>
0046
0047
0048 void fill(const T &candidateCollection,
0049 const C &matchedCandCollection,
0050 float &minVal,
0051 float &maxVal,
0052 const edm::ParameterSet ¶meterSet);
0053 template <class T, class C, class M>
0054 void fill(const T &candidateCollection,
0055 const C &matchedCandCollection,
0056 float &minVal,
0057 float &maxVal,
0058 const edm::ParameterSet ¶meterSet,
0059 const M &muonMatchedCandCollection);
0060
0061 void fillOne(const reco::Candidate &cand);
0062
0063 protected:
0064 CandidateBenchmark candBench_;
0065 MatchCandidateBenchmark matchCandBench_;
0066
0067 TH1F *pt_gen_;
0068 TH1F *eta_gen_;
0069 TH1F *phi_gen_;
0070
0071 TH1F *pt_ref_;
0072 TH1F *eta_ref_;
0073 TH1F *phi_ref_;
0074
0075 TH1F *deltaR_;
0076 float dRMax_;
0077 bool matchCharge_;
0078 bool createReferenceHistos_;
0079 bool histogramBooked_;
0080
0081 bool matching_done_;
0082 bool createEfficiencyHistos_;
0083 };
0084
0085 #include "DQMOffline/PFTau/interface/Matchers.h"
0086 template <class T, class C>
0087 void PFCandidateMonitor::fill(const T &candCollection,
0088 const C &matchedCandCollection,
0089 float &minVal,
0090 float &maxVal,
0091 const edm::ParameterSet ¶meterSet) {
0092 matching_done_ = false;
0093 if (createEfficiencyHistos_) {
0094 for (unsigned i = 0; i < candCollection.size(); ++i) {
0095 if (!isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi()))
0096 continue;
0097 fillOne(candCollection[i]);
0098
0099
0100 for (unsigned j = 0; j < matchedCandCollection.size(); ++j)
0101 if (deltaR_)
0102 deltaR_->Fill(reco::deltaR(candCollection[i], matchedCandCollection[j]));
0103 }
0104 }
0105
0106 std::vector<int> matchIndices;
0107 PFB::match(candCollection, matchedCandCollection, matchIndices, matchCharge_, dRMax_);
0108
0109
0110
0111 matching_done_ = true;
0112
0113 for (unsigned int i = 0; i < (candCollection).size(); i++) {
0114 const reco::Candidate &cand = candCollection[i];
0115
0116 if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
0117 continue;
0118
0119 int iMatch = matchIndices[i];
0120 assert(iMatch < static_cast<int>(matchedCandCollection.size()));
0121
0122 if (iMatch != -1) {
0123 const reco::Candidate &matchedCand = matchedCandCollection[iMatch];
0124 if (!isInRange(matchedCand.pt(), matchedCand.eta(), matchedCand.phi()))
0125 continue;
0126
0127
0128
0129
0130 float ptRes = (cand.pt() - matchedCand.pt()) / matchedCand.pt();
0131
0132 if (ptRes > maxVal)
0133 maxVal = ptRes;
0134 if (ptRes < minVal)
0135 minVal = ptRes;
0136
0137 if (!createEfficiencyHistos_) {
0138 candBench_.fillOne(cand);
0139
0140
0141 matchCandBench_.fillOne(cand, matchedCand,
0142 parameterSet);
0143 if (createReferenceHistos_)
0144 fillOne(matchedCand);
0145
0146 } else {
0147 candBench_.fillOne(matchedCand);
0148
0149
0150
0151 matchCandBench_.fillOne(cand, matchedCand,
0152 parameterSet);
0153 if (createReferenceHistos_)
0154 fillOne(cand);
0155
0156 }
0157 }
0158 }
0159 }
0160
0161 template <class T, class C, class M>
0162
0163
0164 void PFCandidateMonitor::fill(const T &candCollection,
0165 const C &matchedCandCollection,
0166 float &minVal,
0167 float &maxVal,
0168 const edm::ParameterSet ¶meterSet,
0169 const M &muonMatchedCandCollection) {
0170 matching_done_ = false;
0171 if (createEfficiencyHistos_) {
0172 for (unsigned i = 0; i < candCollection.size(); ++i) {
0173 if (!isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi()))
0174 continue;
0175 fillOne(candCollection[i]);
0176
0177
0178 for (unsigned j = 0; j < matchedCandCollection.size(); ++j)
0179 if (deltaR_)
0180 deltaR_->Fill(reco::deltaR(candCollection[i], matchedCandCollection[j]));
0181 }
0182 }
0183
0184 std::vector<int> matchIndices;
0185
0186
0187
0188 PFB::match(candCollection,
0189 matchedCandCollection,
0190 matchIndices,
0191 parameterSet,
0192 muonMatchedCandCollection,
0193 matchCharge_,
0194 dRMax_);
0195
0196 matching_done_ = true;
0197
0198 for (unsigned int i = 0; i < (candCollection).size(); i++) {
0199 const reco::Candidate &cand = candCollection[i];
0200
0201 if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
0202 continue;
0203
0204 int iMatch = matchIndices[i];
0205 assert(iMatch < static_cast<int>(matchedCandCollection.size()));
0206
0207 if (iMatch != -1) {
0208 const reco::Candidate &matchedCand = matchedCandCollection[iMatch];
0209 if (!isInRange(matchedCand.pt(), matchedCand.eta(), matchedCand.phi()))
0210 continue;
0211
0212
0213
0214
0215 float ptRes = (cand.pt() - matchedCand.pt()) / matchedCand.pt();
0216
0217 if (ptRes > maxVal)
0218 maxVal = ptRes;
0219 if (ptRes < minVal)
0220 minVal = ptRes;
0221
0222 if (!createEfficiencyHistos_) {
0223 candBench_.fillOne(cand);
0224 matchCandBench_.fillOne(cand, matchedCand,
0225 parameterSet);
0226 if (createReferenceHistos_)
0227 fillOne(matchedCand);
0228
0229 } else {
0230 candBench_.fillOne(matchedCand);
0231
0232 matchCandBench_.fillOne(cand, matchedCand,
0233 parameterSet);
0234 if (createReferenceHistos_)
0235 fillOne(cand);
0236
0237 }
0238 }
0239 }
0240 }
0241 #endif