00001 #ifndef RecoParticleFlow_Benchmark_GenericBenchmark_h
00002 #define RecoParticleFlow_Benchmark_GenericBenchmark_h
00003
00004
00005 #include "RecoParticleFlow/Benchmark/interface/PFBenchmarkAlgo.h"
00006
00007
00008 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00009 #include "DataFormats/Candidate/interface/Candidate.h"
00010 #include "DataFormats/METReco/interface/MET.h"
00011
00012 #include <string>
00013
00014 #include <TH1F.h>
00015 #include <TH2F.h>
00016 #include <TFile.h>
00017
00018 class DQMStore;
00019
00020 class BenchmarkTree;
00021
00022
00023
00024 class GenericBenchmark{
00025
00026 public:
00027
00028 GenericBenchmark();
00029 virtual ~GenericBenchmark();
00030
00031 void setup(DQMStore *DQM = NULL,
00032 bool PlotAgainstReco_=true,
00033 float minDeltaEt = -100., float maxDeltaEt = 50.,
00034 float minDeltaPhi = -0.5, float maxDeltaPhi = 0.5,
00035 bool doMetPlots=false);
00036
00037 template< typename C>
00038 void fill(const C *RecoCollection,
00039 const C *GenCollection,
00040 bool startFromGen=false,
00041 bool PlotAgainstReco =true,
00042 bool onlyTwoJets = false,
00043 double recPt_cut = -1.,
00044 double minEta_cut = -1.,
00045 double maxEta_cut = -1.,
00046 double deltaR_cut = -1.);
00047
00048
00049 void write(std::string Filename);
00050
00051 void setfile(TFile *file);
00052
00053 private:
00054
00055 bool accepted(const reco::Candidate* particle,
00056 double ptCut,
00057 double minEtaCut,
00058 double maxEtaCut ) const;
00059
00060 void fillHistos( const reco::Candidate* genParticle,
00061 const reco::Candidate* recParticle,
00062 double deltaR_cut,
00063 bool plotAgainstReco);
00064
00065 TFile *file_;
00066
00067 TH1F *hDeltaEt;
00068 TH1F *hDeltaEx;
00069 TH1F *hDeltaEy;
00070 TH2F *hDeltaEtvsEt;
00071 TH2F *hDeltaEtOverEtvsEt;
00072 TH2F *hDeltaEtvsEta;
00073 TH2F *hDeltaEtOverEtvsEta;
00074 TH2F *hDeltaEtvsPhi;
00075 TH2F *hDeltaEtOverEtvsPhi;
00076 TH2F *hDeltaEtvsDeltaR;
00077 TH2F *hDeltaEtOverEtvsDeltaR;
00078
00079 TH2F *hEtRecvsEt;
00080 TH2F *hEtRecOverTrueEtvsTrueEt;
00081
00082 TH1F *hDeltaEta;
00083 TH2F *hDeltaEtavsEt;
00084 TH2F *hDeltaEtavsEta;
00085
00086 TH1F *hDeltaPhi;
00087 TH2F *hDeltaPhivsEt;
00088 TH2F *hDeltaPhivsEta;
00089
00090 TH1F *hDeltaR;
00091 TH2F *hDeltaRvsEt;
00092 TH2F *hDeltaRvsEta;
00093
00094 TH1F *hNRec;
00095
00096
00097 TH1F *hEtGen;
00098 TH1F *hEtaGen;
00099 TH2F *hEtvsEtaGen;
00100 TH2F *hEtvsPhiGen;
00101 TH1F *hPhiGen;
00102
00103 TH1F *hNGen;
00104
00105 TH1F *hEtSeen;
00106 TH1F *hEtaSeen;
00107 TH1F *hPhiSeen;
00108 TH2F *hEtvsEtaSeen;
00109 TH2F *hEtvsPhiSeen;
00110
00111 TH1F *hEtRec;
00112 TH1F *hExRec;
00113 TH1F *hEyRec;
00114 TH1F *hPhiRec;
00115
00116 TH1F *hSumEt;
00117 TH1F *hTrueSumEt;
00118 TH2F *hDeltaSetvsSet;
00119 TH2F *hDeltaMexvsSet;
00120 TH2F *hDeltaSetOverSetvsSet;
00121 TH2F *hRecSetvsTrueSet;
00122 TH2F *hRecSetOverTrueSetvsTrueSet;
00123 TH2F *hTrueMexvsTrueSet;
00124
00125 BenchmarkTree* tree_;
00126
00127 bool doMetPlots_;
00128
00129 protected:
00130
00131 DQMStore *dbe_;
00132 PFBenchmarkAlgo *algo_;
00133
00134 };
00135
00136 template< typename C>
00137 void GenericBenchmark::fill(const C *RecoCollection,
00138 const C *GenCollection,
00139 bool startFromGen,
00140 bool PlotAgainstReco,
00141 bool onlyTwoJets,
00142 double recPt_cut,
00143 double minEta_cut,
00144 double maxEta_cut,
00145 double deltaR_cut) {
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155 if( !startFromGen) {
00156 int nRec = 0;
00157 for (unsigned int i = 0; i < RecoCollection->size(); i++) {
00158
00159
00160 const reco::Candidate *particle = &(*RecoCollection)[i];
00161
00162 assert( particle!=NULL );
00163 if( !accepted(particle, recPt_cut,
00164 minEta_cut, maxEta_cut)) continue;
00165
00166
00167
00168 if( onlyTwoJets ) {
00169 unsigned highJets = 0;
00170 for(unsigned j=0; j<RecoCollection->size(); j++) {
00171 const reco::Candidate *otherParticle = &(*RecoCollection)[j];
00172 if ( j != i && otherParticle->pt() > particle->pt() ) highJets++;
00173 }
00174 if ( highJets > 1 ) continue;
00175 }
00176 nRec++;
00177
00178 const reco::Candidate *gen_particle = algo_->matchByDeltaR(particle,
00179 GenCollection);
00180 if(gen_particle==NULL) continue;
00181
00182
00183
00184
00185 fillHistos( gen_particle, particle, deltaR_cut, PlotAgainstReco);
00186 }
00187
00188 hNRec->Fill(nRec);
00189 }
00190
00191
00192
00193
00194
00195
00196 int nGen = 0;
00197 for (unsigned int i = 0; i < GenCollection->size(); i++) {
00198
00199 const reco::Candidate *gen_particle = &(*GenCollection)[i];
00200
00201 if( !accepted(gen_particle, recPt_cut, minEta_cut, maxEta_cut)) {
00202 continue;
00203 }
00204
00205 hEtaGen->Fill(gen_particle->eta() );
00206 hPhiGen->Fill(gen_particle->phi() );
00207 hEtGen->Fill(gen_particle->et() );
00208 hEtvsEtaGen->Fill(gen_particle->eta(),gen_particle->et() );
00209 hEtvsPhiGen->Fill(gen_particle->phi(),gen_particle->et() );
00210
00211 const reco::Candidate *rec_particle = algo_->matchByDeltaR(gen_particle,
00212 RecoCollection);
00213 nGen++;
00214 if(! rec_particle) continue;
00215
00216
00217 double deltaR = algo_->deltaR(rec_particle,gen_particle);
00218 if (deltaR>deltaR_cut && deltaR_cut != -1.)
00219 continue;
00220
00221 hEtaSeen->Fill(gen_particle->eta() );
00222 hPhiSeen->Fill(gen_particle->phi() );
00223 hEtSeen->Fill(gen_particle->et() );
00224 hEtvsEtaSeen->Fill(gen_particle->eta(),gen_particle->et() );
00225 hEtvsPhiSeen->Fill(gen_particle->phi(),gen_particle->et() );
00226
00227 hPhiRec->Fill(rec_particle->phi() );
00228 hEtRec->Fill(rec_particle->et() );
00229 hExRec->Fill(rec_particle->px() );
00230 hEyRec->Fill(rec_particle->py() );
00231
00232 hEtRecvsEt->Fill(gen_particle->et(),rec_particle->et());
00233 if (gen_particle->et()!=0.0) hEtRecOverTrueEtvsTrueEt->Fill(gen_particle->et(),rec_particle->et()/gen_particle->et());
00234
00235 if( startFromGen )
00236 fillHistos( gen_particle, rec_particle, deltaR_cut, PlotAgainstReco);
00237
00238
00239 }
00240 hNGen->Fill(nGen);
00241
00242 }
00243
00244
00245 #endif // RecoParticleFlow_Benchmark_GenericBenchmark_h