CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenericBenchmark.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_Benchmark_GenericBenchmark_h
2 #define RecoParticleFlow_Benchmark_GenericBenchmark_h
3 
4 //COLIN: necessary?
6 
7 
11 
12 #include <string>
13 
14 #include <TH1F.h>
15 #include <TH2F.h>
16 #include <TFile.h>
17 
18 class DQMStore; // CMSSW_2_X_X
19 
20 class BenchmarkTree;
21 
22 
23 //COLIN: this class REALLY needs to be cleaned up and rationalized, on the model of PFCandidateBenchmark.
25 
26  public:
27 
29  virtual ~GenericBenchmark();
30 
31  void setup(DQMStore *DQM = NULL,
32  bool PlotAgainstReco_=true,
33  float minDeltaEt = -100., float maxDeltaEt = 50.,
34  float minDeltaPhi = -0.5, float maxDeltaPhi = 0.5,
35  bool doMetPlots=false);
36 
37  template< typename C>
38  void fill(const C *RecoCollection,
39  const C *GenCollection,
40  bool startFromGen=false,
41  bool PlotAgainstReco =true,
42  bool onlyTwoJets = false,
43  double recPt_cut = -1.,
44  double minEta_cut = -1.,
45  double maxEta_cut = -1.,
46  double deltaR_cut = -1.);
47 
48 
49  void write(std::string Filename);
50 
51  void setfile(TFile *file);
52 
53  private:
54 
55  bool accepted(const reco::Candidate* particle,
56  double ptCut,
57  double minEtaCut,
58  double maxEtaCut ) const;
59 
60  void fillHistos( const reco::Candidate* genParticle,
61  const reco::Candidate* recParticle,
62  double deltaR_cut,
63  bool plotAgainstReco);
64 
65  TFile *file_;
66 
67  TH1F *hDeltaEt;
68  TH1F *hDeltaEx;
69  TH1F *hDeltaEy;
70  TH2F *hDeltaEtvsEt;
78 
79  TH2F *hEtRecvsEt;
81 
82  TH1F *hDeltaEta;
85 
86  TH1F *hDeltaPhi;
89 
90  TH1F *hDeltaR;
91  TH2F *hDeltaRvsEt;
92  TH2F *hDeltaRvsEta;
93 
94  TH1F *hNRec;
95 
96 
97  TH1F *hEtGen;
98  TH1F *hEtaGen;
99  TH2F *hEtvsEtaGen;
100  TH2F *hEtvsPhiGen;
101  TH1F *hPhiGen;
102 
103  TH1F *hNGen;
104 
105  TH1F *hEtSeen;
106  TH1F *hEtaSeen;
107  TH1F *hPhiSeen;
110 
111  TH1F *hEtRec;
112  TH1F *hExRec;
113  TH1F *hEyRec;
114  TH1F *hPhiRec;
115 
116  TH1F *hSumEt;
117  TH1F *hTrueSumEt;
124 
126 
128 
129  protected:
130 
133 
134 };
135 
136 template< typename C>
137 void GenericBenchmark::fill(const C *RecoCollection,
138  const C *GenCollection,
139  bool startFromGen,
140  bool PlotAgainstReco,
141  bool onlyTwoJets,
142  double recPt_cut,
143  double minEta_cut,
144  double maxEta_cut,
145  double deltaR_cut) {
146 
147  //if (doMetPlots_)
148  //{
149  // const reco::MET* met1=static_cast<const reco::MET*>(&((*RecoCollection)[0]));
150  // if (met1!=NULL) std::cout << "FL: met1.sumEt() = " << (*met1).sumEt() << std::endl;
151  //}
152 
153  // loop over reco particles
154 
155  if( !startFromGen) {
156  int nRec = 0;
157  for (unsigned int i = 0; i < RecoCollection->size(); i++) {
158 
159  // generate histograms comparing the reco and truth candidate (truth = closest in delta-R)
160  const reco::Candidate *particle = &(*RecoCollection)[i];
161 
162  assert( particle!=NULL );
163  if( !accepted(particle, recPt_cut,
164  minEta_cut, maxEta_cut)) continue;
165 
166 
167  // Count the number of jets with a larger energy
168  if( onlyTwoJets ) {
169  unsigned highJets = 0;
170  for(unsigned j=0; j<RecoCollection->size(); j++) {
171  const reco::Candidate *otherParticle = &(*RecoCollection)[j];
172  if ( j != i && otherParticle->pt() > particle->pt() ) highJets++;
173  }
174  if ( highJets > 1 ) continue;
175  }
176  nRec++;
177 
178  const reco::Candidate *gen_particle = algo_->matchByDeltaR(particle,
179  GenCollection);
180  if(gen_particle==NULL) continue;
181 
182 
183 
184  // fill histograms
185  fillHistos( gen_particle, particle, deltaR_cut, PlotAgainstReco);
186  }
187 
188  hNRec->Fill(nRec);
189  }
190 
191  // loop over gen particles
192 
193  // std::cout<<"Reco size = "<<RecoCollection->size()<<", ";
194  // std::cout<<"Gen size = "<<GenCollection->size()<<std::endl;
195 
196  int nGen = 0;
197  for (unsigned int i = 0; i < GenCollection->size(); i++) {
198 
199  const reco::Candidate *gen_particle = &(*GenCollection)[i];
200 
201  if( !accepted(gen_particle, recPt_cut, minEta_cut, maxEta_cut)) {
202  continue;
203  }
204 
205  hEtaGen->Fill(gen_particle->eta() );
206  hPhiGen->Fill(gen_particle->phi() );
207  hEtGen->Fill(gen_particle->et() );
208  hEtvsEtaGen->Fill(gen_particle->eta(),gen_particle->et() );
209  hEtvsPhiGen->Fill(gen_particle->phi(),gen_particle->et() );
210 
211  const reco::Candidate *rec_particle = algo_->matchByDeltaR(gen_particle,
212  RecoCollection);
213  nGen++;
214  if(! rec_particle) continue; // no match
215  // must make a cut on delta R - so let's do the cut
216 
217  double deltaR = algo_->deltaR(rec_particle,gen_particle);
218  if (deltaR>deltaR_cut && deltaR_cut != -1.)
219  continue;
220 
221  hEtaSeen->Fill(gen_particle->eta() );
222  hPhiSeen->Fill(gen_particle->phi() );
223  hEtSeen->Fill(gen_particle->et() );
224  hEtvsEtaSeen->Fill(gen_particle->eta(),gen_particle->et() );
225  hEtvsPhiSeen->Fill(gen_particle->phi(),gen_particle->et() );
226 
227  hPhiRec->Fill(rec_particle->phi() );
228  hEtRec->Fill(rec_particle->et() );
229  hExRec->Fill(rec_particle->px() );
230  hEyRec->Fill(rec_particle->py() );
231 
232  hEtRecvsEt->Fill(gen_particle->et(),rec_particle->et());
233  if (gen_particle->et()!=0.0) hEtRecOverTrueEtvsTrueEt->Fill(gen_particle->et(),rec_particle->et()/gen_particle->et());
234 
235  if( startFromGen )
236  fillHistos( gen_particle, rec_particle, deltaR_cut, PlotAgainstReco);
237 
238 
239  }
240  hNGen->Fill(nGen);
241 
242 }
243 
244 
245 #endif // RecoParticleFlow_Benchmark_GenericBenchmark_h
void fillHistos(const reco::Candidate *genParticle, const reco::Candidate *recParticle, double deltaR_cut, bool plotAgainstReco)
void fill(const C *RecoCollection, const C *GenCollection, bool startFromGen=false, bool PlotAgainstReco=true, bool onlyTwoJets=false, double recPt_cut=-1., double minEta_cut=-1., double maxEta_cut=-1., double deltaR_cut=-1.)
int i
Definition: DBlmapReader.cc:9
virtual ~GenericBenchmark()
bool onlyTwoJets
virtual double et() const =0
transverse energy
virtual float eta() const =0
momentum pseudorapidity
void setfile(TFile *file)
PFBenchmarkAlgo * algo_
virtual float phi() const =0
momentum azimuthal angle
#define NULL
Definition: scimark2.h:8
bool accepted(const reco::Candidate *particle, double ptCut, double minEtaCut, double maxEtaCut) const
virtual float pt() const =0
transverse momentum
tuple DQM
Definition: DQM_cfg.py:6
TH2F * hRecSetOverTrueSetvsTrueSet
static double deltaR(const T *, const U *)
int j
Definition: DBlmapReader.cc:9
TH2F * hDeltaEtOverEtvsDeltaR
virtual double py() const =0
y coordinate of momentum vector
virtual double px() const =0
x coordinate of momentum vector
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
bool plotAgainstReco
static const Collection::value_type * matchByDeltaR(const T *, const Collection *)
TH2F * hEtRecOverTrueEtvsTrueEt
void write(std::string Filename)
BenchmarkTree * tree_
void setup(DQMStore *DQM=NULL, bool PlotAgainstReco_=true, float minDeltaEt=-100., float maxDeltaEt=50., float minDeltaPhi=-0.5, float maxDeltaPhi=0.5, bool doMetPlots=false)