CMS 3D CMS Logo

PFCandidateManager.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_Benchmark_BenchmarkManager_h
2 #define RecoParticleFlow_Benchmark_BenchmarkManager_h
3 
8 
11 
12 #include <vector>
13 
26 class PFCandidateManager : public Benchmark {
27 public:
29  : Benchmark(mode),
33  dRMax_(dRMax),
35 
36  ~PFCandidateManager() override;
37 
39  void setParameters(float dRMax = 0.3, bool matchCharge = true, Benchmark::Mode mode = Benchmark::DEFAULT);
40 
42  void setDirectory(TDirectory *dir) override;
43 
45  void setup(DQMStore::IBooker &b);
46 
48  template <class C>
49  void fill(const reco::PFCandidateCollection &candCollection, const C &matchedCandCollection);
50 
51 protected:
55 
56  float dRMax_;
58 };
59 
61 
62 template <class C>
63 void PFCandidateManager::fill(const reco::PFCandidateCollection &candCollection, const C &matchCandCollection) {
64  std::vector<int> matchIndices;
65  PFB::match(candCollection, matchCandCollection, matchIndices, matchCharge_, dRMax_);
66 
67  for (unsigned int i = 0; i < candCollection.size(); i++) {
68  const reco::PFCandidate &cand = candCollection[i];
69 
70  if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
71  continue;
72 
73  int iMatch = matchIndices[i];
74 
75  assert(iMatch < static_cast<int>(matchCandCollection.size()));
76 
77  // COLIN how to handle efficiency plots?
78 
79  // filling the histograms in CandidateBenchmark only in case
80  // of a matching.
81  if (iMatch != -1) {
82  candBench_.fillOne(cand);
83  pfCandBench_.fillOne(cand);
84  matchCandBench_.fillOne(cand, matchCandCollection[iMatch]);
85  }
86  }
87 }
88 
89 #endif
A benchmark managing several benchmarks.
To plot Candidate quantities.
double eta() const final
momentum pseudorapidity
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
To plot Candidate quantities.
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
abstract base class
Definition: Benchmark.h:19
double pt() const final
transverse momentum
CandidateBenchmark candBench_
PFCandidateManager(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
~PFCandidateManager() override
void setParameters(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
set the benchmark parameters
void setup(DQMStore::IBooker &b)
book histograms
void setDirectory(TDirectory *dir) override
set directory (to use in ROOT)
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
MatchCandidateBenchmark matchCandBench_
double b
Definition: hdecay.h:118
void fillOne(const reco::PFCandidate &pfCand)
fill histograms with a given particle
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
void fill(const reco::PFCandidateCollection &candCollection, const C &matchedCandCollection)
fill histograms with all particle
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:50
double phi() const final
momentum azimuthal angle
PFCandidateBenchmark pfCandBench_