CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
24 class PFCandidateManager : public Benchmark {
25 
26  public:
27 
28  PFCandidateManager( float dRMax = 0.3,
29  bool matchCharge = true,
31  :
32  Benchmark(mode),
34  dRMax_(dRMax), matchCharge_(matchCharge) {}
35 
36  virtual ~PFCandidateManager();
37 
39  void setParameters( float dRMax = 0.3,
40  bool matchCharge = true,
42 
44  void setDirectory(TDirectory* dir);
45 
47  void setup(DQMStore::IBooker& b);
48 
50  template< class C>
51  void fill(const reco::PFCandidateCollection& candCollection,
52  const C& matchedCandCollection );
53 
54  protected:
58 
59  float dRMax_;
61 
62 };
63 
64 
66 
67 template< class C>
69  const C& matchCandCollection) {
70 
71 
72  std::vector<int> matchIndices;
73  PFB::match( candCollection, matchCandCollection, matchIndices,
75 
76  for (unsigned int i = 0; i < candCollection.size(); i++) {
77  const reco::PFCandidate& cand = candCollection[i];
78 
79  if( !isInRange(cand.pt(), cand.eta(), cand.phi() ) ) continue;
80 
81  int iMatch = matchIndices[i];
82 
83  assert(iMatch< static_cast<int>(matchCandCollection.size()));
84 
85  //COLIN how to handle efficiency plots?
86 
87  // filling the histograms in CandidateBenchmark only in case
88  // of a matching.
89  if( iMatch!=-1 ) {
90  candBench_.fillOne(cand);
91  pfCandBench_.fillOne(cand);
92  matchCandBench_.fillOne(cand, matchCandCollection[ iMatch ]);
93  }
94  }
95 }
96 
97 
98 
99 #endif
int i
Definition: DBlmapReader.cc:9
A benchmark managing several benchmarks.
To plot Candidate quantities.
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
assert(m_qm.get())
abstract base class
Definition: Benchmark.h:22
CandidateBenchmark candBench_
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
PFCandidateManager(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
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)
set directory (to use in ROOT)
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
MatchCandidateBenchmark matchCandBench_
double b
Definition: hdecay.h:120
void fillOne(const reco::PFCandidate &pfCand)
fill histograms with a given particle
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
void fill(const reco::PFCandidateCollection &candCollection, const C &matchedCandCollection)
fill histograms with all particle
dbl *** dir
Definition: mlp_gen.cc:35
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
virtual double phi() const
momentum azimuthal angle
PFCandidateBenchmark pfCandBench_