ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/ForwardAnalysis/Utilities/plugins/SimpleJetResponseAnalyzer.cc
Revision: 1.1
Committed: Thu Aug 25 19:11:15 2011 UTC (13 years, 8 months ago) by antoniov
Content type: text/plain
Branch: MAIN
CVS Tags: V01-01-01, V01-01-00, antoniov-forwardAnalysis-09Jul2012-v1, antoniov-forwardAnalysis-29Jun2012-v1, V01-00-00, antoniov-utilities-11Jun2012-v1, antoniov-forwardAnalysis-Oct072011-v1, sfonseca_10_04_2011, antoniov-forwardAnalysis-Sep182011-v1, antoniov-forwardAnalysis-Sep102011-v1, eliza_09_02_2011, sfonseca_08_26_2011, HEAD
Log Message:
moving plugins to ForwardAnalysis

File Contents

# User Rev Content
1 antoniov 1.1
2     #include "FWCore/Framework/interface/EDAnalyzer.h"
3     #include "FWCore/Utilities/interface/InputTag.h"
4     #include "FWCore/Framework/interface/Frameworkfwd.h"
5    
6     class TTree;
7    
8     class SimpleJetResponseAnalyzer: public edm::EDAnalyzer
9     {
10     public:
11     explicit SimpleJetResponseAnalyzer(const edm::ParameterSet&);
12     ~SimpleJetResponseAnalyzer();
13    
14     virtual void beginJob();
15     virtual void analyze(const edm::Event&, const edm::EventSetup&);
16    
17    
18     private:
19     //edm::InputTag genJetTag_;
20     edm::InputTag genToRecoJetMatchTag_;
21    
22     TTree* data_;
23    
24     struct {
25     double genJetPt_;
26     double genJetEta_;
27     double genJetPhi_;
28     double recJetPt_;
29     double recJetEta_;
30     double recJetPhi_;
31     } eventData_;
32     };
33    
34    
35     #include "FWCore/Framework/interface/Event.h"
36     #include "FWCore/ParameterSet/interface/ParameterSet.h"
37     #include "FWCore/MessageLogger/interface/MessageLogger.h"
38    
39     #include "FWCore/Framework/interface/MakerMacros.h"
40    
41     //#include "DataFormats/JetReco/interface/Jet.h"
42     //#include "DataFormats/JetReco/interface/PFJet.h"
43     #include "DataFormats/Candidate/interface/CandidateFwd.h"
44     #include "DataFormats/Candidate/interface/CandMatchMap.h"
45    
46     #include "FWCore/ServiceRegistry/interface/Service.h"
47     #include "CommonTools/UtilAlgos/interface/TFileService.h"
48    
49     #include "TTree.h"
50    
51     using namespace reco;
52    
53     SimpleJetResponseAnalyzer::SimpleJetResponseAnalyzer(const edm::ParameterSet& pset):
54     genToRecoJetMatchTag_(pset.getParameter<edm::InputTag>("GenToRecoJetMatchTag"))
55     {
56     }
57    
58     SimpleJetResponseAnalyzer::~SimpleJetResponseAnalyzer() {}
59    
60     void SimpleJetResponseAnalyzer::beginJob(){
61     edm::Service<TFileService> fs;
62    
63     data_ = fs->make<TTree>("data","data");
64     data_->Branch("genJetPt",&eventData_.genJetPt_,"genJetPt/D");
65     data_->Branch("genJetEta",&eventData_.genJetEta_,"genJetEta/D");
66     data_->Branch("genJetPhi",&eventData_.genJetPhi_,"genJetPhi/D");
67     data_->Branch("recJetPt",&eventData_.recJetPt_,"recJetPt/D");
68     data_->Branch("recJetEta",&eventData_.recJetEta_,"recJetEta/D");
69     data_->Branch("recJetPhi",&eventData_.recJetPhi_,"recJetPhi/D");
70     }
71    
72     void SimpleJetResponseAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup){
73     edm::Handle<CandViewMatchMap> genToRecoJetAssociationH;
74     event.getByLabel(genToRecoJetMatchTag_,genToRecoJetAssociationH);
75     const CandViewMatchMap& genToRecoJetAssociation = *genToRecoJetAssociationH;
76    
77     CandViewMatchMap::const_iterator it_key = genToRecoJetAssociation.begin();
78     CandViewMatchMap::const_iterator match_end = genToRecoJetAssociation.end();
79     for(; it_key != match_end; ++it_key){
80     const CandidateBaseRef& genJet = it_key->key;
81     const CandidateBaseRef& recJet = it_key->val;
82    
83     double response = recJet->pt()/genJet->pt();
84     double phiDiff = recJet->phi() - genJet->phi();
85     double etaDiff = recJet->eta() - genJet->eta();
86     double ptDiff = recJet->pt() - genJet->pt();
87    
88     LogTrace("") << " gen. jet pt,eta,phi= " << genJet->pt() << ", " << genJet->eta() << ", " << genJet->phi()
89     << "\n rec. jet pt,eta,phi= " << recJet->pt() << ", " << recJet->eta() << ", " << recJet->phi()
90     << "\n response,phiDiff,etaDiff,ptDiff= " << response << ", " << phiDiff << ", " << etaDiff << ", " << ptDiff;
91    
92     eventData_.genJetPt_ = genJet->pt();
93     eventData_.genJetEta_ = genJet->eta();
94     eventData_.genJetPhi_ = genJet->phi();
95     eventData_.recJetPt_ = recJet->pt();
96     eventData_.recJetEta_ = recJet->eta();
97     eventData_.recJetPhi_ = recJet->phi();
98    
99     data_->Fill();
100     }
101     }
102    
103     DEFINE_FWK_MODULE(SimpleJetResponseAnalyzer);