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 
14 #include <TH1.h> //needed by the deltaR->Fill() call
15 
16 class PFCandidateMonitor : public Benchmark {
17 
18  public:
19 
20  PFCandidateMonitor( float dRMax = 0.3, bool matchCharge = true,
22 
23  virtual ~PFCandidateMonitor();
24 
26  void setParameters(float dRMax, bool matchCharge, Benchmark::Mode mode,
27  float ptmin, float ptmax, float etamin, float etamax,
28  float phimin, float phimax, bool refHistoFlag);
29 
32 
34  void setDirectory(TDirectory* dir);
35 
37  void setup();
38 
40  void setup(const edm::ParameterSet & parameterSet);
41 
43  template< class T, class C>
44  /*void fill(const T& candidateCollection,
45  const C& matchedCandCollection, float& minVal, float& maxVal);*/
46 
47  void fill(const T& candidateCollection,
48  const C& matchedCandCollection, float& minVal, float& maxVal,
49  const edm::ParameterSet & parameterSet);
50 
51  void fillOne(const reco::Candidate& cand);
52 
53  protected:
56 
57  TH1F* pt_gen_;
58  TH1F* eta_gen_;
59  TH1F* phi_gen_;
60 
61  TH1F* pt_ref_;
62  TH1F* eta_ref_;
63  TH1F* phi_ref_;
64 
65  TH1F* deltaR_;
66  float dRMax_;
70 
73 };
74 
76 template< class T, class C>
77  /*void PFCandidateMonitor::fill(const T& candCollection,
78  const C& matchedCandCollection, float& minVal, float& maxVal) {*/
79  void PFCandidateMonitor::fill(const T& candCollection,
80  const C& matchedCandCollection, float& minVal, float& maxVal,
82 
83  matching_done_ = false;
85  for( unsigned i=0; i<candCollection.size(); ++i) {
86  if( !isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi() ) ) continue;
87  fillOne(candCollection[i]); // fill pt_gen, eta_gen and phi_gen histos for UNMATCHED generated candidate
88 
89  for( unsigned j=0; j<matchedCandCollection.size(); ++j) // for DeltaR spectrum
90  if (deltaR_) deltaR_->Fill( reco::deltaR( candCollection[i], matchedCandCollection[j] ) ) ;
91  }
92  }
93 
94  std::vector<int> matchIndices;
95  PFB::match( candCollection, matchedCandCollection, matchIndices, matchCharge_, dRMax_ );
96  // now matchIndices[i] stores the j-th closest matched jet
97  matching_done_ = true;
98 
99  for (unsigned int i = 0; i < (candCollection).size(); i++) {
100  const reco::Candidate& cand = candCollection[i];
101 
102  if( !isInRange(cand.pt(), cand.eta(), cand.phi() ) ) continue;
103 
104  int iMatch = matchIndices[i];
105  assert(iMatch< static_cast<int>(matchedCandCollection.size()));
106 
107  if( iMatch!=-1 ) {
108  const reco::Candidate& matchedCand = matchedCandCollection[ iMatch ];
109  if(!isInRange(matchedCand.pt(),matchedCand.eta(),matchedCand.phi() ) ) continue;
110  //std::cout <<"PFJet pT " <<cand.pt() <<" eta " <<cand.eta() <<" phi " <<cand.phi() ;
111  //std::cout <<"\nmatched genJet pT " <<matchedCand.pt() <<" eta " <<matchedCand.eta() <<" phi " <<matchedCand.phi() <<"\n" <<std::endl ;
112  float ptRes = (cand.pt() - matchedCand.pt()) / matchedCand.pt();
113 
114  if (ptRes > maxVal) maxVal = ptRes;
115  if (ptRes < minVal) minVal = ptRes;
116 
117  if ( !createEfficiencyHistos_ ) {
118  candBench_.fillOne(cand); // fill pt, eta phi and charge histos for MATCHED candidate
119  //matchCandBench_.fillOne(cand, matchedCand); // fill delta_x_VS_y histos for matched couple
120  matchCandBench_.fillOne(cand, matchedCand, parameterSet); // fill delta_x_VS_y histos for matched couple
121  if (createReferenceHistos_) fillOne(matchedCand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED reference candidate
122  }
123  else {
124  candBench_.fillOne(matchedCand); // fill pt, eta phi and charge histos for MATCHED candidate
125  //matchCandBench_.fillOne(matchedCand, cand); // fill delta_x_VS_y histos for matched couple
126  matchCandBench_.fillOne(cand, matchedCand, parameterSet); // fill delta_x_VS_y histos for matched couple
127  if (createReferenceHistos_) fillOne(cand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED reference candidate
128  }
129 
130  }
131  }
132 }
133 #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:14
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:21
virtual float phi() const =0
momentum azimuthal angle
T eta() const
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 fill(const T &candidateCollection, const C &matchedCandCollection, float &minVal, float &maxVal, const edm::ParameterSet &parameterSet)
fill histograms with all particle
void setDirectory(TDirectory *dir)
set directory (to use in ROOT)
int j
Definition: DBlmapReader.cc:9
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
CandidateBenchmark candBench_
void fillOne(const reco::Candidate &cand)
double ptmin
Definition: HydjetWrapper.h:85
MatchCandidateBenchmark matchCandBench_
dbl *** dir
Definition: mlp_gen.cc:35
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:59
long double T
tuple size
Write out results.
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
PFCandidateMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)
Definition: DDAxes.h:10