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(DQMStore::IBooker& b);
38  void setup(DQMStore::IBooker& b, 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  void fill(const T& candidateCollection,
45  const C& matchedCandCollection, float& minVal, float& maxVal,
46  const edm::ParameterSet & parameterSet) ;
47  template< class T, class C, class M>
48  void fill(const T& candidateCollection,
49  const C& matchedCandCollection, float& minVal, float& maxVal,
50  const edm::ParameterSet & parameterSet, const M& muonMatchedCandCollection ) ;
51 
52  void fillOne(const reco::Candidate& cand);
53 
54  protected:
57 
58  TH1F* pt_gen_;
59  TH1F* eta_gen_;
60  TH1F* phi_gen_;
61 
62  TH1F* pt_ref_;
63  TH1F* eta_ref_;
64  TH1F* phi_ref_;
65 
66  TH1F* deltaR_;
67  float dRMax_;
71 
74 };
75 
77 template< class T, class C>
78  void PFCandidateMonitor::fill(const T& candCollection,
79  const C& matchedCandCollection, float& minVal, float& maxVal,
81 
82  matching_done_ = false;
84  for( unsigned i=0; i<candCollection.size(); ++i) {
85  if( !isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi() ) ) continue;
86  fillOne(candCollection[i]); // fill pt_gen, eta_gen and phi_gen histos for UNMATCHED generated candidate
87 
88  for( unsigned j=0; j<matchedCandCollection.size(); ++j) // for DeltaR spectrum
89  if (deltaR_) deltaR_->Fill( reco::deltaR( candCollection[i], matchedCandCollection[j] ) ) ;
90  }
91  }
92 
93  std::vector<int> matchIndices;
94  PFB::match( candCollection, matchedCandCollection, matchIndices, matchCharge_, dRMax_ );
95  //PFB::match( candCollection, matchedCandCollection, matchIndices, parameterSet, 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 
134 template< class T, class C, class M>
135  /*void PFCandidateMonitor::fill(const T& candCollection,
136  const C& matchedCandCollection, float& minVal, float& maxVal) {*/
137  void PFCandidateMonitor::fill(const T& candCollection,
138  const C& matchedCandCollection, float& minVal, float& maxVal,
139  const edm::ParameterSet & parameterSet, const M& muonMatchedCandCollection) {
140 
141 
142  matching_done_ = false;
143  if ( createEfficiencyHistos_ ) {
144  for( unsigned i=0; i<candCollection.size(); ++i) {
145  if( !isInRange(candCollection[i].pt(), candCollection[i].eta(), candCollection[i].phi() ) ) continue;
146  fillOne(candCollection[i]); // fill pt_gen, eta_gen and phi_gen histos for UNMATCHED generated candidate
147 
148  for( unsigned j=0; j<matchedCandCollection.size(); ++j) // for DeltaR spectrum
149  if (deltaR_) deltaR_->Fill( reco::deltaR( candCollection[i], matchedCandCollection[j] ) ) ;
150  }
151  }
152 
153  std::vector<int> matchIndices;
154  //PFB::match( candCollection, matchedCandCollection, matchIndices, matchCharge_, dRMax_ );
155  //PFB::match( candCollection, matchedCandCollection, matchIndices, parameterSet, matchCharge_, dRMax_ );
156  PFB::match( candCollection, matchedCandCollection, matchIndices, parameterSet, muonMatchedCandCollection, matchCharge_, dRMax_ );
157  // now matchIndices[i] stores the j-th closest matched jet
158  matching_done_ = true;
159 
160  for (unsigned int i = 0; i < (candCollection).size(); i++) {
161  const reco::Candidate& cand = candCollection[i];
162 
163  if( !isInRange(cand.pt(), cand.eta(), cand.phi() ) ) continue;
164 
165  int iMatch = matchIndices[i];
166  assert(iMatch< static_cast<int>(matchedCandCollection.size()));
167 
168  if( iMatch!=-1 ) {
169  const reco::Candidate& matchedCand = matchedCandCollection[ iMatch ];
170  if(!isInRange(matchedCand.pt(),matchedCand.eta(),matchedCand.phi() ) ) continue;
171  //std::cout <<"PFJet pT " <<cand.pt() <<" eta " <<cand.eta() <<" phi " <<cand.phi() ;
172  //std::cout <<"\nmatched genJet pT " <<matchedCand.pt() <<" eta " <<matchedCand.eta() <<" phi " <<matchedCand.phi() <<"\n" <<std::endl ;
173  float ptRes = (cand.pt() - matchedCand.pt()) / matchedCand.pt();
174 
175  if (ptRes > maxVal) maxVal = ptRes;
176  if (ptRes < minVal) minVal = ptRes;
177 
178  if ( !createEfficiencyHistos_ ) {
179  candBench_.fillOne(cand); // fill pt, eta phi and charge histos for MATCHED candidate
180  matchCandBench_.fillOne(cand, matchedCand, parameterSet); // fill delta_x_VS_y histos for matched couple
181  if (createReferenceHistos_) fillOne(matchedCand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED reference candidate
182  }
183  else {
184  candBench_.fillOne(matchedCand); // fill pt, eta phi and charge histos for MATCHED candidate
185  matchCandBench_.fillOne(cand, matchedCand, parameterSet); // fill delta_x_VS_y histos for matched couple
186  if (createReferenceHistos_) fillOne(cand); // fill pt_ref, eta_ref and phi_ref histos for MATCHED reference candidate
187  }
188 
189  }
190  }
191 }
192 #endif
int i
Definition: DBlmapReader.cc:9
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
virtual double pt() const =0
transverse momentum
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
assert(m_qm.get())
abstract base class
Definition: Benchmark.h:22
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
void setup(DQMStore::IBooker &b)
book histograms
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
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
CandidateBenchmark candBench_
void fillOne(const reco::Candidate &cand)
double b
Definition: hdecay.h:120
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:58
long double T
tuple size
Write out results.
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
PFCandidateMonitor(float dRMax=0.3, bool matchCharge=true, Benchmark::Mode mode=Benchmark::DEFAULT)