CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/Calibration/IsolatedParticles/plugins/IsolatedParticlesGeneratedJets.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    IsolatedParticlesGeneratedJets
00004 // Class:      IsolatedParticlesGeneratedJets
00005 // 
00013 //
00014 // Original Author:  Seema Sharma
00015 //         Created:  Thu Mar  4 18:52:02 CST 2010
00016 //
00017 //
00018 
00019 #include "Calibration/IsolatedParticles/plugins/IsolatedParticlesGeneratedJets.h"
00020 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00021 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00022 
00023 IsolatedParticlesGeneratedJets::IsolatedParticlesGeneratedJets(const edm::ParameterSet& iConfig) {
00024 
00025   debug   = iConfig.getUntrackedParameter<bool>  ("Debug", false);
00026   jetSrc  = iConfig.getParameter<edm::InputTag>("JetSource");
00027   partSrc = iConfig.getParameter<edm::InputTag>("ParticleSource");
00028 }
00029 
00030 
00031 IsolatedParticlesGeneratedJets::~IsolatedParticlesGeneratedJets() {
00032 
00033 }
00034 
00035 void IsolatedParticlesGeneratedJets::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00036   
00037   //using namespace edm;
00038   clearTreeVectors();
00039 
00040   //=== genJet information
00041   edm::Handle<reco::GenJetCollection> genJets;
00042   iEvent.getByLabel(jetSrc, genJets);
00043 
00044   //=== genJet information
00045   edm::Handle<reco::GenParticleCollection> genParticles;
00046   iEvent.getByLabel(partSrc, genParticles);
00047 
00048   JetMatchingTools jetMatching (iEvent);
00049   std::vector <std::vector <const reco::GenParticle*> > genJetConstituents (genJets->size());
00050 
00051   int njets = 0;
00052   for (unsigned iGenJet = 0; iGenJet < genJets->size(); ++iGenJet) {
00053     const reco::GenJet& genJet = (*genJets) [iGenJet];
00054 
00055     double genJetE   = genJet.energy();
00056     double genJetPt  = genJet.pt();
00057     double genJetEta = genJet.eta();
00058     double genJetPhi = genJet.phi();
00059 
00060     if( genJetPt> 30.0 && std::abs(genJetEta)<3.0 ) {
00061 
00062       njets++;
00063 
00064       std::vector <const reco::GenParticle*> genJetConstituents = jetMatching.getGenParticles ((*genJets) [iGenJet]);
00065 
00066       std::vector<double> v_trkP, v_trkPt, v_trkEta, v_trkPhi, v_trkPdg, v_trkCharge;
00067     
00068       if(debug) std::cout<<"Jet(pt,Eta,Phi) "<<genJetPt<<" "<<genJetEta<<" "<<genJetPhi <<std::endl;
00069       for(unsigned int ic=0; ic<genJetConstituents.size(); ic++) {
00070 
00071         if(debug) {
00072           std::cout << "p,pt,eta,phi "<<genJetConstituents[ic]->p()<<" "<<genJetConstituents[ic]->pt()
00073                     <<" "<<genJetConstituents[ic]->eta()<<" "<<genJetConstituents[ic]->phi()
00074                     <<std::endl;
00075         }
00076 
00077         v_trkP.push_back(genJetConstituents[ic]->p());
00078         v_trkPt.push_back(genJetConstituents[ic]->pt());
00079         v_trkEta.push_back(genJetConstituents[ic]->eta());
00080         v_trkPhi.push_back(genJetConstituents[ic]->phi());
00081         v_trkPdg.push_back(genJetConstituents[ic]->pdgId());
00082         v_trkCharge.push_back(genJetConstituents[ic]->charge());
00083 
00084       } //loop over genjet constituents
00085 
00086       t_gjetE   ->push_back(genJetE  );
00087       t_gjetPt  ->push_back(genJetPt );
00088       t_gjetEta ->push_back(genJetEta);
00089       t_gjetPhi ->push_back(genJetPhi);
00090       
00091       t_jetTrkP   ->push_back(v_trkP  );
00092       t_jetTrkPt  ->push_back(v_trkPt );
00093       t_jetTrkEta ->push_back(v_trkEta);
00094       t_jetTrkPhi ->push_back(v_trkPhi);
00095       t_jetTrkPdg ->push_back(v_trkPdg);
00096       t_jetTrkCharge ->push_back(v_trkCharge);
00097 
00098     } // if jetPt>30
00099 
00100   } //loop over genjets
00101 
00102   t_gjetN->push_back(njets);
00103 
00104   unsigned int indx = 0;
00105   for(reco::GenParticleCollection::const_iterator ig = genParticles->begin(); ig!= genParticles->end(); ++ig,++indx) {
00106  
00107     if (debug)
00108       std::cout << "Track " << indx << " Status " << ig->status() << " charge "
00109                 << ig->charge() << " pdgId " << ig->pdgId() << " mass "
00110                 << ig->mass() << " P " << ig->momentum() << " E "
00111                 << ig->energy() << " Origin " << ig->vertex() << std::endl;
00112   }
00113 
00114 
00115   tree->Fill();
00116 }
00117 
00118 void IsolatedParticlesGeneratedJets::beginJob() {
00119 
00120   BookHistograms();
00121 
00122 }
00123 
00124 void IsolatedParticlesGeneratedJets::clearTreeVectors() {
00125   t_gjetN     ->clear();
00126   t_gjetE     ->clear();
00127   t_gjetPt    ->clear();
00128   t_gjetEta   ->clear();
00129   t_gjetPhi   ->clear();
00130 
00131   t_jetTrkP   ->clear();
00132   t_jetTrkPt  ->clear();
00133   t_jetTrkEta ->clear();
00134   t_jetTrkPhi ->clear();
00135   t_jetTrkPdg ->clear();
00136   t_jetTrkCharge ->clear();
00137 }
00138 
00139 void IsolatedParticlesGeneratedJets::BookHistograms(){
00140 
00141   tree = fs->make<TTree>("tree", "tree");
00142 
00143   t_gjetN     = new std::vector<int>   ();
00144   t_gjetE     = new std::vector<double>();
00145   t_gjetPt    = new std::vector<double>();
00146   t_gjetEta   = new std::vector<double>();
00147   t_gjetPhi   = new std::vector<double>();
00148 
00149   t_jetTrkP   = new std::vector<std::vector<double> >();
00150   t_jetTrkPt  = new std::vector<std::vector<double> >();
00151   t_jetTrkEta = new std::vector<std::vector<double> >();
00152   t_jetTrkPhi = new std::vector<std::vector<double> >();
00153   t_jetTrkPdg = new std::vector<std::vector<double> >();
00154   t_jetTrkCharge = new std::vector<std::vector<double> >();
00155 
00156   tree->Branch("t_gjetN",     "vector<int>",             &t_gjetN);
00157   tree->Branch("t_gjetE",     "vector<double>",          &t_gjetE);
00158   tree->Branch("t_gjetPt",    "vector<double>",          &t_gjetPt);
00159   tree->Branch("t_gjetEta",   "vector<double>",          &t_gjetEta);
00160   tree->Branch("t_gjetPhi",   "vector<double>",          &t_gjetPhi);
00161 
00162   tree->Branch("t_jetTrkP",   "vector<vector<double> >", &t_jetTrkP);
00163   tree->Branch("t_jetTrkPt",  "vector<vector<double> >", &t_jetTrkPt);
00164   tree->Branch("t_jetTrkEta", "vector<vector<double> >", &t_jetTrkEta);
00165   tree->Branch("t_jetTrkPhi", "vector<vector<double> >", &t_jetTrkPhi);
00166   tree->Branch("t_jetTrkPdg", "vector<vector<double> >", &t_jetTrkPdg);
00167   tree->Branch("t_jetTrkCharge", "vector<vector<double> >", &t_jetTrkCharge);
00168 
00169 }
00170 
00171 void IsolatedParticlesGeneratedJets::endJob() {
00172 }
00173 
00174 //define this as a plug-in
00175 DEFINE_FWK_MODULE(IsolatedParticlesGeneratedJets);