CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateMonitor.h
Go to the documentation of this file.
1 #ifndef DQMOffline_PFTau_PFCandidateMonitor_h
2 #define DQMOffline_PFTau_PFCandidateMonitor_h
3 
8 
11 
12 #include <vector>
13 class PFCandidateMonitor : public Benchmark {
14 
15  public:
16 
17  PFCandidateMonitor( float dRMax = 0.3,
18  bool matchCharge = true,
20 
21  virtual ~PFCandidateMonitor();
22 
24  void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode,
25  float ptmin, float ptmax, float etamin, float etamax,
26  float phimin, float phimax, bool refHistoFlag);
27 
29  void setParameters( const edm::ParameterSet& parameterSet);
30 
32  void setDirectory(TDirectory* dir);
33 
35  void setup();
36 
38  void setup(const edm::ParameterSet & parameterSet);
39 
41  template< class T, class C>
42  void fill(const T& candidateCollection,
43  const C& matchedCandCollection, float& minVal, float& maxVal);
44 
45 
46  void fillOne(const reco::Candidate& cand);
47 
48  protected:
51 
52  TH1F* pt_ref_;
53  TH1F* eta_ref_;
54  TH1F* phi_ref_;
55 
56  float dRMax_;
60 
61 };
62 
64 template< class T, class C>
65 void PFCandidateMonitor::fill(const T& candCollection,
66  const C& matchedCandCollection, float& minVal, float& maxVal) {
67 
68 
69  std::vector<int> matchIndices;
70  PFB::match( candCollection, matchedCandCollection, matchIndices,
72 
73  for (unsigned int i = 0; i < (candCollection).size(); i++) {
74  const reco::Candidate& cand = candCollection[i];
75 
76  if( !isInRange(cand.pt(), cand.eta(), cand.phi() ) ) continue;
77 
78  int iMatch = matchIndices[i];
79  assert(iMatch< static_cast<int>(matchedCandCollection.size()));
80 
81  if( iMatch!=-1 ) {
82  const reco::Candidate& matchedCand = matchedCandCollection[ iMatch ];
83  if(!isInRange(matchedCand.pt(),matchedCand.eta(),matchedCand.phi() ) ) continue;
84  float ptRes = (cand.pt() - matchedCand.pt())/matchedCand.pt();
85 
86  if (ptRes > maxVal) maxVal = ptRes;
87  if (ptRes < minVal) minVal = ptRes;
88 
89  candBench_.fillOne(cand);
90  matchCandBench_.fillOne(cand, matchedCand);
91  if (createReferenceHistos_) fillOne(matchedCand);
92  }
93  }
94 }
95 #endif
int i
Definition: DBlmapReader.cc:9
void setup()
book histograms
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:13
virtual float eta() const =0
momentum pseudorapidity
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
abstract base class
Definition: Benchmark.h:20
virtual float phi() const =0
momentum azimuthal angle
void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode, float ptmin, float ptmax, float etamin, float etamax, float phimin, float phimax, bool refHistoFlag)
set the parameters locally
virtual float pt() const =0
transverse momentum
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
void fill(const T &candidateCollection, const C &matchedCandCollection, float &minVal, float &maxVal)
fill histograms with all particle
CandidateBenchmark candBench_
void fillOne(const reco::Candidate &cand)
double ptmin
Definition: HydjetWrapper.h:86
MatchCandidateBenchmark matchCandBench_
dbl *** dir
Definition: mlp_gen.cc:35
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
long double T
tuple size
Write out results.
PFCandidateMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)