ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/auterman/Demo/PATJetIDAnalyzer/src/PATJetIDAnalyzer.cc
Revision: 1.1.1.1 (vendor branch)
Committed: Mon Feb 25 15:54:04 2008 UTC (17 years, 2 months ago) by auterman
Content type: text/plain
Branch: tex, JetID, Demo
Changes since 1.1: +0 -0 lines
Log Message:
PAT Jet Analyzer

File Contents

# User Rev Content
1 auterman 1.1 // -*- C++ -*-
2     //
3     // Package: PATJetIDAnalyzer
4     // Class: PATJetIDAnalyzer
5     //
6     /**\class PATJetIDAnalyzer PATJetIDAnalyzer.cc Demo/PATJetIDAnalyzer/src/PATJetIDAnalyzer.cc
7    
8     Description: <one line class summary>
9    
10     Implementation:
11     <Notes on implementation>
12     */
13     //
14     // Original Author: Christian Autermann
15     // Created: Mon Feb 25 11:33:02 CET 2008
16     // $Id$
17     //
18     //
19    
20     //system includes
21     #include <fstream>
22     #include <iostream>
23     //framework includes
24     #include "FWCore/Utilities/interface/EDMException.h"
25     #include "FWCore/ServiceRegistry/interface/Service.h"
26     #include "PhysicsTools/UtilAlgos/interface/TFileService.h"
27     //root includes
28     #include "TH1F.h"
29     #include "TH2F.h"
30     #include "TString.h"
31     //User include files
32     #include "DataFormats/Math/interface/deltaR.h"
33     //AOD
34     #include "DataFormats/JetReco/interface/GenJet.h"
35     #include "DataFormats/JetReco/interface/GenJetfwd.h"
36     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
37     #include "DataFormats/METReco/interface/CaloMETCollection.h"
38     #include "DataFormats/METReco/interface/GenMETCollection.h"
39     #include "DataFormats/HepMCCandidate/interface/GenParticleCandidate.h"
40     //PAT
41     #include "DataFormats/PatCandidates/interface/Jet.h"
42     #include "DataFormats/PatCandidates/interface/MET.h"
43    
44     #include "Demo/PATJetIDAnalyzer/interface/NameScheme.h"
45     #include "Demo/PATJetIDAnalyzer/interface/PATJetIDAnalyzer.h"
46    
47    
48    
49     //
50     // constructors and destructor
51     //
52     PATJetIDAnalyzer::PATJetIDAnalyzer(const edm::ParameterSet& iConfig) :
53     _recJet ( iConfig.getParameter<edm::InputTag>( "recJet" ) ),
54     _genJet ( iConfig.getParameter<edm::InputTag>( "genJet" ) ),
55     _recMet ( iConfig.getParameter<edm::InputTag>( "recMet" ) ),
56     _genMet ( iConfig.getParameter<edm::InputTag>( "genMet" ) ),
57     _hist ( iConfig.getParameter<std::string>( "hist" ) ),
58     _jetminpt( iConfig.getParameter<double>( "MinJetPt" ) ),
59     _jetmaxeta(iConfig.getParameter<double>( "MaxJetEta" ) )
60     {
61     //now do what ever initialization is needed
62    
63     }
64    
65    
66     PATJetIDAnalyzer::~PATJetIDAnalyzer()
67     {
68    
69     // do anything here that needs to be done at desctruction time
70     // (e.g. close files, deallocate resources etc.)
71    
72     }
73    
74    
75     //
76     // member functions
77     //
78    
79     // ------------ method called to for each event ------------
80     void
81     PATJetIDAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
82     {
83     using namespace edm;
84     using namespace std;
85    
86     //CaloJets
87     edm::Handle<reco::CaloJetCollection> CJets;
88     iEvent.getByLabel( _recJet, CJets );
89     reco::CaloJetCollection CaloJets = *CJets;
90     std::sort( CaloJets.begin(), CaloJets.end(), PtGreater());
91    
92     //GenJets
93     edm::Handle<reco::GenJetCollection> GJets;
94     iEvent.getByLabel( _genJet, GJets );
95     reco::GenJetCollection GenJets = *GJets;
96     std::sort( GenJets.begin(), GenJets.end(), PtGreater());
97    
98     // matching maps for GenJets and CaloJets and vice versa
99     // fills some histograms inside (dR vs. dE, rec. eff., gen. vs. matched ...)
100     makeMatchingMaps(GJets,CJets);
101    
102     //CaloJets
103     unsigned multiplicity = 0;
104     double met=0.0, ht=0.0;
105     for (reco::CaloJetCollection::const_iterator it=CaloJets.begin();
106     it!=CaloJets.end(); ++it) {
107     if (it->pt() > _jetminpt && fabs( it->eta() )<_jetmaxeta){
108     if (multiplicity<_njets ) {
109     _pt_jet[ multiplicity]->Fill( it->pt() );
110     _eta_jet[ multiplicity]->Fill( it->eta() );
111     _phi_jet[ multiplicity]->Fill( it->phi() );
112     _emfrac_jet[ multiplicity]->Fill( it->emEnergyFraction() );
113     _hadfrac_jet[multiplicity]->Fill( it->energyFractionHadronic() );
114     _n60_jet[ multiplicity]->Fill( it->n60() );
115     _n90_jet[ multiplicity]->Fill( it->n90() );
116     _area_jet[ multiplicity]->Fill( it->towersArea() );
117     }
118     ++multiplicity;
119     ht += it->pt();
120     }
121     }
122     _jetmult->Fill( multiplicity );
123     _ht->Fill( ht );
124    
125     //MET
126     edm::Handle<reco::CaloMETCollection> recMet;
127     iEvent.getByLabel(_recMet,recMet);
128    
129     for ( reco::CaloMETCollection::const_iterator it=recMet->begin();
130     it!=recMet->end(); ++it) {
131     met = it->pt();
132     _met->Fill( met );
133     _metx->Fill(met*sin(it->phi()));
134     _mety->Fill(met*cos(it->phi()));
135     break;
136     }
137     _metmult->Fill( recMet->size() );
138     _htmet->Fill( ht+met );
139     _metsig->Fill( met / (0.97*sqrt(ht)) );
140     if (multiplicity>1 ) {
141     _nv->Fill( (met*met-(CaloJets[0].pt()-CaloJets[1].pt())*(CaloJets[0].pt()-CaloJets[1].pt())) / (CaloJets[1].pt()*CaloJets[1].pt()) );
142     _nv2->Fill( (met*met-(CaloJets[0].pt()-CaloJets[1].pt())*(CaloJets[0].pt()-CaloJets[1].pt())) / (CaloJets[0].pt()*CaloJets[0].pt()) );
143     _dijet->Fill(sqrt(pow(CaloJets[0].energy()+CaloJets[1].energy(),2)
144     -(pow(CaloJets[0].px()+CaloJets[1].px(),2)
145     +pow(CaloJets[0].py()+CaloJets[1].py(),2)
146     +pow(CaloJets[0].pz()+CaloJets[1].pz(),2) ) ) );
147     }
148    
149     //GenJets
150     multiplicity = 0;
151     met=0.0;
152     ht=0.0;
153     for (reco::GenJetCollection::const_iterator it=GenJets.begin();
154     it!=GenJets.end(); ++it) {
155     if (it->pt() > _jetminpt && fabs( it->eta() )<_jetmaxeta) {
156     if (multiplicity<_ngenjets ) {
157     _pt_genjet[ multiplicity]->Fill( it->pt() );
158     _eta_genjet[ multiplicity]->Fill( it->eta() );
159     _phi_genjet[ multiplicity]->Fill( it->phi() );
160     _emfrac_genjet[ multiplicity]->Fill( it->emEnergy() / it->energy() );
161     _hadfrac_genjet[ multiplicity]->Fill( it->hadEnergy() / it->energy() );
162     _invisible_genjet[multiplicity]->Fill( it->invisibleEnergy() / it->energy() );
163     _aux_genjet[ multiplicity]->Fill( it->auxiliaryEnergy() / it->energy() );
164     }
165     ++multiplicity;
166     ht += it->pt();
167     }
168     }
169     _genjetmult->Fill( multiplicity );
170    
171     //GenMET
172     edm::Handle<reco::GenMETCollection> genMet;
173     iEvent.getByLabel(_genMet,genMet);
174    
175     for ( reco::GenMETCollection::const_iterator it=genMet->begin();
176     it!=genMet->end(); ++it) {
177     met = it->pt();
178     _genmet->Fill( met );
179     _genmetx->Fill(met*sin(it->phi()));
180     _genmety->Fill(met*cos(it->phi()));
181     break;
182     }
183     _genmetmult->Fill( recMet->size() );
184     _genht->Fill( ht );
185     _genhtmet->Fill( ht+met );
186     if (multiplicity>1 ) {
187     _gendijet->Fill(sqrt(pow(GenJets[0].energy()+GenJets[1].energy(),2)
188     -(pow(GenJets[0].px()+GenJets[1].px(),2)
189     +pow(GenJets[0].py()+GenJets[1].py(),2)
190     +pow(GenJets[0].pz()+GenJets[1].pz(),2) ) ) );
191     }
192    
193     }
194    
195    
196     // ------------ method called once each job just before starting event loop ------------
197     void
198     PATJetIDAnalyzer::beginJob(const edm::EventSetup&)
199     {
200     if( _hist.empty() )
201     return;
202    
203     edm::Service<TFileService> fs;
204     if( !fs ){
205     throw edm::Exception( edm::errors::Configuration,
206     "TFile Service is not registered in cfg file" );
207     }
208    
209     NameScheme hsusy("hsusy");
210     ofstream hist(_hist.c_str(), std::ios::out);
211     //CaloJets
212     for (unsigned i=0; i<_njets; ++i) {
213     _pt_jet[i] = fs->make<TH1F>(hsusy.name(hist, "ptj", i ), hsusy.name("PtJet", i), 100, 0. , 1000.);
214     _eta_jet[i] = fs->make<TH1F>(hsusy.name(hist, "etaj", i ), hsusy.name("EtaJet", i), 100, -5. , 5.);
215     _phi_jet[i] = fs->make<TH1F>(hsusy.name(hist, "phij", i ), hsusy.name("PhiJet", i), 70, -3.5, 3.5);
216     _emfrac_jet[i] = fs->make<TH1F>(hsusy.name(hist, "emfj", i ), hsusy.name("EmfJet", i), 50, 0. , 2.);
217     _hadfrac_jet[i] = fs->make<TH1F>(hsusy.name(hist, "hadfj",i ), hsusy.name("HadfJet",i), 50, 0. , 2.);
218     _n60_jet[i] = fs->make<TH1F>(hsusy.name(hist, "n60j", i ), hsusy.name("n60Jet", i), 50, 0. , 50.);
219     _n90_jet[i] = fs->make<TH1F>(hsusy.name(hist, "n90j", i ), hsusy.name("n90Jet", i), 50, 0. , 50.);
220     _area_jet[i] = fs->make<TH1F>(hsusy.name(hist, "areaj",i ), hsusy.name("AreaJet",i), 50, 0. , 5.);
221     }
222     _jetmult = fs->make<TH1F>(hsusy.name(hist, "multj" ), hsusy.name("JetMultiplicity"), 21, -0.5, 20.5);
223     _ht = fs->make<TH1F>(hsusy.name(hist, "ht" ), hsusy.name("Ht"), 100, 0.0, 2000.);
224     _htmet = fs->make<TH1F>(hsusy.name(hist, "htmet" ), hsusy.name("HtMet"), 100, 0.0, 2000.);
225     _dijet = fs->make<TH1F>(hsusy.name(hist, "dijet" ), hsusy.name("DiJet"), 100, 0.0, 2000.);
226     //GenJets
227     for (unsigned i=0; i<_ngenjets; ++i) {
228     _pt_genjet[i] = fs->make<TH1F>(hsusy.name(hist, "ptgj", i ), hsusy.name("PtGenJet", i), 100, 0. , 1000.);
229     _eta_genjet[i] = fs->make<TH1F>(hsusy.name(hist, "etagj", i ), hsusy.name("EtaGenJet", i), 100, -5. , 5.);
230     _phi_genjet[i] = fs->make<TH1F>(hsusy.name(hist, "phigj", i ), hsusy.name("PhiGenJet", i), 70, -3.5, 3.5);
231     _emfrac_genjet[i] = fs->make<TH1F>(hsusy.name(hist, "emfgj", i ), hsusy.name("EmfGenJet", i), 50, 0. , 5.);
232     _hadfrac_genjet[i] = fs->make<TH1F>(hsusy.name(hist, "hadfgj",i ), hsusy.name("HadfGenJet",i), 50, 0. , 5.);
233     _invisible_genjet[i]=fs->make<TH1F>(hsusy.name(hist, "invisgj",i), hsusy.name("InvisiblefGenJet", i), 50, 0.,5.);
234     _aux_genjet[i] = fs->make<TH1F>(hsusy.name(hist, "auxgj", i ), hsusy.name("AuxfGenJet", i), 50, 0. , 5.);
235     }
236     _genjetmult = fs->make<TH1F>(hsusy.name(hist, "multgj" ), hsusy.name("GenJetMultiplicity"), 21, -0.5, 20.5);
237     _genht = fs->make<TH1F>(hsusy.name(hist, "htgj" ), hsusy.name("GenHt"), 100, 0.0, 2000.);
238     _genhtmet = fs->make<TH1F>(hsusy.name(hist, "htmetgj"), hsusy.name("GenHtMet"), 100, 0.0, 2000.);
239     _gendijet = fs->make<TH1F>(hsusy.name(hist, "dijetgj"), hsusy.name("GenDiJet"), 100, 0.0, 2000.);
240     //MET
241     _met = fs->make<TH1F>(hsusy.name(hist, "met" ), hsusy.name("Met"), 100, 0.0, 1000.);
242     _metmult = fs->make<TH1F>(hsusy.name(hist, "metmult"),hsusy.name("MetMult"), 11, -0.5, 10.5);
243     _metx = fs->make<TH1F>(hsusy.name(hist, "metx" ), hsusy.name("MetX"), 100, -1000.0, 1000.);
244     _mety = fs->make<TH1F>(hsusy.name(hist, "mety" ), hsusy.name("MetY"), 100, -1000.0, 1000.);
245     _metsig = fs->make<TH1F>(hsusy.name(hist, "metsig"), hsusy.name("MetSig"), 100, 0.0, 100.);
246     //GenMET
247     _genmet = fs->make<TH1F>(hsusy.name(hist, "metgj" ), hsusy.name("GenMet"), 100, 0.0, 1000.);
248     _genmetmult = fs->make<TH1F>(hsusy.name(hist, "metmultgj"),hsusy.name("GenMetMult"), 11, -0.5, 10.5);
249     _genmetx = fs->make<TH1F>(hsusy.name(hist, "metxgj" ), hsusy.name("GenMetX"), 100, -1000.0, 1000.);
250     _genmety = fs->make<TH1F>(hsusy.name(hist, "metygj" ), hsusy.name("GenMetY"), 100, -1000.0, 1000.);
251     //other
252     _nv = fs->make<TH1F>(hsusy.name(hist, "nv" ), hsusy.name("NV"), 100, -5.0, 5.);
253     _nv2 = fs->make<TH1F>(hsusy.name(hist, "nv2" ), hsusy.name("NV2"), 100, -5.0, 5.);
254     //matching
255     _GenOnCalo_match = fs->make<TH2F>(hsusy.name(hist,"gencalo"), hsusy.name("GenOnCalo"), 60, 0., 3., 120, -1, 5.);
256     _CaloOnGen_match = fs->make<TH2F>(hsusy.name(hist,"calogen"), hsusy.name("CaloOnGen"), 60, 0., 3., 120, -1, 5.);
257     _GenVsMatched_pt = fs->make<TH2F>(hsusy.name(hist,"genmatched_pt"), hsusy.name("GenVsMatched_Pt") , 100, 0. , 2000. , 100, 0. , 2000.);
258     _GenVsMatched_eta = fs->make<TH2F>(hsusy.name(hist,"genmatched_eta"), hsusy.name("GenVsMatched_Eta"), 100, -5. , 5. , 100, -5. , 5.);
259     _GenVsMatched_phi = fs->make<TH2F>(hsusy.name(hist,"genmatched_phi"), hsusy.name("GenVsMatched_Phi"), 70, -3.5, 3.5, 70, -3.5, 3.5);
260     _RecoEff_pt = fs->make<TH1F>(hsusy.name(hist, "recoeffpt") , hsusy.name("RecoEff_Pt") , 100, 0. , 2000. );
261     _RecoEff_eta = fs->make<TH1F>(hsusy.name(hist, "recoeffeta"), hsusy.name("RecoEff_Eta"), 100, -5. , 5. );
262     _RecoEff_phi = fs->make<TH1F>(hsusy.name(hist, "recoeffphi"), hsusy.name("RecoEff_Phi"), 70, -3.5, 3.5);
263     //helper histograms
264     _RecoTot_pt = new TH1F("recototpt" , "RecoTot_Pt" , 100, 0. , 2000. );
265     _RecoTot_eta= new TH1F("recototeta", "RecoTot_Eta", 100, -5. , 5. );
266     _RecoTot_phi= new TH1F("recototphi", "RecoTot_Phi", 70, -3.5, 3.5);
267     }
268    
269     // ------------ method called once each job just after ending the event loop ------------
270     void
271     PATJetIDAnalyzer::endJob() {
272     _RecoTot_pt->Sumw2();
273     _RecoEff_pt->Sumw2();
274     _RecoEff_pt->Divide(_RecoTot_pt);
275     _RecoTot_eta->Sumw2();
276     _RecoEff_eta->Sumw2();
277     _RecoEff_eta->Divide(_RecoTot_eta);
278     _RecoTot_phi->Sumw2();
279     _RecoEff_phi->Sumw2();
280     _RecoEff_phi->Divide(_RecoTot_phi);
281     }
282    
283     // ------------ redo matching
284     void
285     PATJetIDAnalyzer::makeMatchingMaps(edm::Handle<reco::GenJetCollection> GenJets, edm::Handle<reco::CaloJetCollection> CaloJets)
286     {
287     MatchingMapGen.clear();
288     MatchingMapJet.clear();
289    
290     // Loop over generator Jets
291     for(reco::GenJetCollection::const_iterator genjet = GenJets->begin(); genjet!=GenJets->end(); ++genjet) {
292     const reco::CaloJet* jetbest=0;
293     const reco::CaloJet* jetmatched=0;
294     double minDelR= 0.3;
295     double minDelE=-1/4.;
296     double maxDelE= 3.0;
297     double mindE=999.;
298     double mindR=999.;
299     double delR;
300     double delE;
301     double E1=(*genjet).energy();
302     bool matched = false;
303     // Loop over Calo jets
304     for(reco::CaloJetCollection::const_iterator calojet = CaloJets->begin(); calojet!=CaloJets->end(); ++calojet){
305     double E2=(*calojet).energy();
306     delR=deltaR(*calojet,*genjet);
307     delE=(E2-E1)/E1;
308     // This is the matching condition
309     if(delR<minDelR && delE>minDelE && delE<maxDelE){
310     matched=true;
311     minDelR=delR;
312     jetmatched=&(*calojet);
313     }
314     if(delR<mindR){
315     mindR=delR;
316     mindE=delE;
317     jetbest=&(*calojet);
318     }
319     }
320     if (matched) {
321     MatchingMapGen.insert(make_pair(&(*genjet),&(*jetmatched)));
322     _GenVsMatched_pt->Fill((*genjet).pt(),(*jetmatched).pt());
323     _GenVsMatched_eta->Fill((*genjet).eta(),(*jetmatched).eta());
324     _GenVsMatched_phi->Fill((*genjet).phi(),(*jetmatched).phi());
325     _RecoEff_pt->Fill((*genjet).pt());
326     _RecoEff_eta->Fill((*genjet).eta());
327     _RecoEff_phi->Fill((*genjet).phi());
328     }
329     _RecoTot_pt->Fill((*genjet).pt());
330     _RecoTot_eta->Fill((*genjet).eta());
331     _RecoTot_phi->Fill((*genjet).phi());
332     _GenOnCalo_match->Fill(mindR,mindE);
333     }
334    
335     // Loop over Calo Jets
336     for(reco::CaloJetCollection::const_iterator calojet = CaloJets->begin(); calojet!=CaloJets->end(); ++calojet) {
337     const reco::GenJet* jetbest=0;
338     const reco::GenJet* jetmatched=0;
339     double minDelR= 0.3;
340     double minDelE=-1/4.;
341     double maxDelE= 3.0;
342     double mindE=999.;
343     double mindR=999.;
344     double delR;
345     double delE;
346     double E1=(*calojet).energy();
347     bool matched = false;
348     // Loop over jets
349     for(reco::GenJetCollection::const_iterator genjet = GenJets->begin(); genjet!=GenJets->end(); ++genjet){
350     double E2=(*genjet).energy();
351     delR=deltaR(*genjet,*calojet);
352     delE=(E2-E1)/E1;
353     // This is the matching condition
354     if(delR<minDelR && delE>minDelE && delE<maxDelE){
355     matched=true;
356     minDelR=delR;
357     jetmatched=&(*genjet);
358     }
359     if(delR<mindR){
360     mindR=delR;
361     mindE=delE;
362     jetbest=&(*genjet);
363     }
364     }
365     if (matched) {
366     MatchingMapJet.insert(make_pair(&(*calojet),&(*jetmatched)));
367     }
368     _CaloOnGen_match->Fill(mindR,mindE);
369     }
370    
371     }
372    
373     //define this as a plug-in
374     DEFINE_FWK_MODULE(PATJetIDAnalyzer);