CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/RecoJets/JetAnalyzers/src/JetAnaPythia.cc

Go to the documentation of this file.
00001 // Name: JetAnaPythia
00002 // Description:  Example of analysis of Pythia produced partons & jets
00003 //               Based on Kostas Kousouris' templated JetPlotsExample.
00004 //               Plots are tailored to needs of dijet mass and ratio analysis.
00005 // Author: R. Harris
00006 // Date:  28 - Oct - 2008
00007 #include "RecoJets/JetAnalyzers/interface/JetAnaPythia.h"
00008 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00009 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00010 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00011 #include "DataFormats/JetReco/interface/CaloJet.h"
00012 #include "DataFormats/JetReco/interface/PFJet.h"
00013 #include "DataFormats/JetReco/interface/GenJet.h"
00014 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00015 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00016 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00017 #include "FWCore/Framework/interface/Event.h"
00018 #include "FWCore/Framework/interface/Run.h"
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include <TFile.h>
00021 #include <cmath>
00022 using namespace edm;
00023 using namespace reco;
00024 using namespace std;
00026 template<class Jet>
00027 JetAnaPythia<Jet>::JetAnaPythia(edm::ParameterSet const& cfg)
00028 {
00029   JetAlgorithm  = cfg.getParameter<std::string> ("JetAlgorithm"); 
00030   HistoFileName = cfg.getParameter<std::string> ("HistoFileName");
00031   NJets         = cfg.getParameter<int> ("NJets");
00032   debug         = cfg.getParameter<bool> ("debug");
00033   eventsGen     = cfg.getParameter<int> ("eventsGen");
00034   anaLevel      = cfg.getParameter<std::string> ("anaLevel");
00035   xsecGen       = cfg.getParameter< vector<double> > ("xsecGen");
00036   ptHatEdges    = cfg.getParameter< vector<double> > ("ptHatEdges");
00037 
00038 }
00040 template<class Jet>
00041 void JetAnaPythia<Jet>::beginJob() 
00042 {
00043   TString hname;
00044   m_file = new TFile(HistoFileName.c_str(),"RECREATE"); 
00046   const int nMassBins = 103;
00047   double massBoundaries[nMassBins+1] = {1, 3, 6, 10, 16, 23, 31, 40, 50, 61, 74, 88, 103, 119, 137, 156, 176, 197, 220, 244, 270, 296, 325, 354, 386, 419, 453, 489, 526, 565, 606, 649, 693, 740, 788, 838, 890, 944, 1000, 1058, 1118, 1181, 1246, 1313, 1383, 1455, 1530, 1607, 1687, 1770, 1856, 1945, 2037, 2132, 2231, 2332, 2438, 2546, 2659, 2775, 2895, 3019, 3147, 3279, 3416, 3558, 3704, 3854, 4010, 4171, 4337, 4509, 4686, 4869, 5058, 5253, 5455, 5663, 5877, 6099, 6328, 6564, 6808, 7060, 7320, 7589, 7866, 8152, 8447, 8752, 9067, 9391, 9726, 10072, 10430, 10798, 11179, 11571, 11977, 12395, 12827, 13272, 13732, 14000};  
00048 
00049   hname = "JetPt";
00050   m_HistNames1D[hname] = new TH1F(hname,hname,500,0,5000);
00051 
00052   hname = "JetEta";
00053   m_HistNames1D[hname] = new TH1F(hname,hname,120,-6,6);
00054 
00055   hname = "JetPhi";
00056   m_HistNames1D[hname] = new TH1F(hname,hname,100,-M_PI,M_PI);
00057 
00058   hname = "NumberOfJets";
00059   m_HistNames1D[hname] = new TH1F(hname,hname,100,0,100);
00060 
00061   hname = "DijetMass";
00062   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00063 
00064   hname = "DijetMassWt";
00065   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00066   m_HistNames1D.find(hname)->second->Sumw2();
00067 
00068   hname = "DijetMassIn";
00069   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00070 
00071   hname = "DijetMassInWt";
00072   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00073   m_HistNames1D.find(hname)->second->Sumw2();
00074 
00075   hname = "DijetMassOut";
00076   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00077 
00078   hname = "DijetMassOutWt";
00079   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00080   m_HistNames1D.find(hname)->second->Sumw2();
00081 
00082   hname = "ResonanceMass";
00083   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00084 
00085   hname = "DipartonMass";
00086   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00087 
00088   hname = "DipartonMassWt";
00089   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00090   m_HistNames1D.find(hname)->second->Sumw2();
00091 
00092   hname = "DipartonMassIn";
00093   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00094 
00095   hname = "DipartonMassInWt";
00096   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00097   m_HistNames1D.find(hname)->second->Sumw2();
00098 
00099   hname = "DipartonMassOut";
00100   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00101 
00102   hname = "DipartonMassOutWt";
00103   m_HistNames1D[hname] = new TH1F(hname,hname,nMassBins,massBoundaries);
00104   m_HistNames1D.find(hname)->second->Sumw2();
00105 
00106   hname = "PtHat";
00107   m_HistNames1D[hname] = new TH1F(hname,hname,1000,0,5000);
00108 
00109   hname = "PtHatWt";
00110   m_HistNames1D[hname] = new TH1F(hname,hname,1000,0,5000);
00111   m_HistNames1D.find(hname)->second->Sumw2();
00112 
00113   hname = "PtHatFine";
00114   m_HistNames1D[hname] = new TH1F(hname,hname,5000,0,5000);
00115 
00116   hname = "PtHatFineWt";
00117   m_HistNames1D[hname] = new TH1F(hname,hname,5000,0,5000);
00118   m_HistNames1D.find(hname)->second->Sumw2();
00119 
00120   mcTruthTree_   = new TTree("mcTruthTree","mcTruthTree");
00121   mcTruthTree_->Branch("xsec",     &xsec,      "xsec/F");
00122   mcTruthTree_->Branch("weight",     &weight,      "weight/F");
00123   mcTruthTree_->Branch("pt_hat",     &pt_hat,      "pt_hat/F");
00124   mcTruthTree_->Branch("nJets",     &nJets,      "nJets/I");
00125   mcTruthTree_->Branch("etaJet1",     &etaJet1,      "etaJet1/F");
00126   mcTruthTree_->Branch("etaJet2",     &etaJet2,      "etaJet2/F");
00127   mcTruthTree_->Branch("ptJet1",     &ptJet1,      "ptJet1/F");
00128   mcTruthTree_->Branch("ptJet2",     &ptJet2,      "ptJet2/F");
00129   mcTruthTree_->Branch("diJetMass",     &diJetMass,      "diJetMass/F");
00130   mcTruthTree_->Branch("etaPart1",     &etaPart1,      "etaPart1/F");
00131   mcTruthTree_->Branch("etaPart2",     &etaPart2,      "etaPart2/F");
00132   mcTruthTree_->Branch("ptPart1",     &ptPart1,      "ptPart1/F");
00133   mcTruthTree_->Branch("ptPart2",     &ptPart2,      "ptPart2/F");
00134   mcTruthTree_->Branch("diPartMass",     &diPartMass,      "diPartMass/F");
00135 
00136   
00137 }
00139 template<class Jet>
00140 void JetAnaPythia<Jet>::analyze(edm::Event const& evt, edm::EventSetup const& iSetup) 
00141 {
00142   int notDone=1;
00143   while(notDone)
00144   {  //while loop to allow us to tailor the analysis level for faster running.
00145     TString hname; 
00146       
00147     // Process Info
00148   
00149     //edm::Handle< double > genEventScale;
00150     //evt.getByLabel("genEventScale", genEventScale );
00151     //pt_hat = *genEventScale;
00152 
00153     edm::Handle<edm::HepMCProduct> MCevt;
00154     evt.getByLabel("generator", MCevt);
00155     HepMC::GenEvent * myGenEvent = new HepMC::GenEvent(*(MCevt->GetEvent()));
00156    
00157     double pthat = myGenEvent->event_scale();
00158     pt_hat = float(pthat);
00159 
00160     delete myGenEvent;
00161 
00162     if(anaLevel != "generating"){  //We are not generating events, so xsec is there
00163       //edm::Handle< GenRunInfoProduct > genInfoProduct;
00165       //xsec = (double)genInfoProduct->externalXSecLO();
00166        xsec=0.0;
00167        if( ptHatEdges.size()>xsecGen.size() ){
00168          for(unsigned int i_pthat = 0; i_pthat < xsecGen.size(); ++i_pthat){
00169             if( pthat >= ptHatEdges[i_pthat] &&  pthat < ptHatEdges[i_pthat+1])xsec=float(xsecGen[i_pthat]);
00170          }
00171        }
00172        else 
00173        {
00174         std::cout << "Number of PtHat bin edges too small. Xsec set to zero" << std::endl;  
00175        }
00176     }
00177     else
00178     {                        
00179       xsec = xsecGen[0];   //Generating events, no xsec in event, get xsec from user input
00180     } 
00181     if(debug)std::cout << "cross section=" <<xsec << " pb" << std::endl;           
00182     weight =  xsec/eventsGen;
00183 
00184     if(debug)std::cout << "pt_hat=" <<pt_hat  <<  std::endl;
00185     hname = "PtHat";
00186     FillHist1D(hname, pt_hat, 1.0); 
00187     hname = "PtHatFine";
00188     FillHist1D(hname, pt_hat, 1.0); 
00189     hname = "PtHatWt";
00190     FillHist1D(hname, pt_hat, weight); 
00191     hname = "PtHatFineWt";
00192     FillHist1D(hname, pt_hat, weight); 
00193     if(anaLevel=="PtHatOnly")break;  //ptHatOnly should be very fast
00194 
00195     // Jet Info
00196     math::XYZTLorentzVector p4jet[2];
00197     float etajet[2];
00199     Handle<JetCollection> jets;
00200     evt.getByLabel(JetAlgorithm,jets);
00201     typename JetCollection::const_iterator i_jet;
00202     int index = 0;
00203 
00205     hname = "NumberOfJets";
00206     nJets = jets->size();
00207     FillHist1D(hname,nJets,1.0); 
00208     
00209   
00210     // Two Leading Jet Info
00211     for(i_jet = jets->begin(); i_jet != jets->end() && index < 2; ++i_jet) 
00212       {
00213         hname = "JetPt";
00214         FillHist1D(hname,i_jet->pt(),1.0);   
00215         hname = "JetEta";
00216         FillHist1D(hname,i_jet->eta(),1.0);
00217         hname = "JetPhi";
00218         FillHist1D(hname,i_jet->phi(),1.0);
00219         p4jet[index] = i_jet->p4();
00220         etajet[index] = i_jet->eta();
00221         if(debug)std::cout << "jet " << index+1 <<": pt=" <<i_jet->pt() << ", eta=" <<etajet[index] <<  std::endl;
00222         index++;
00223       }
00224 
00225       // TTree variables //
00226       etaJet1 = etajet[0];
00227       etaJet2 = etajet[1];
00228       ptJet1 = p4jet[0].pt();
00229       ptJet2 = p4jet[1].pt();
00230       diJetMass = (p4jet[0]+p4jet[1]).mass();
00231 
00233       if(index==2&&abs(etaJet1)<1.3&&abs(etaJet2)<1.3){
00234        hname = "DijetMass";
00235        FillHist1D(hname,diJetMass ,1.0); 
00236        hname = "DijetMassWt";
00237        FillHist1D(hname,diJetMass ,weight); 
00238      }
00239 
00241       if(index==2&&abs(etaJet1)<0.7&&abs(etaJet2)<0.7){
00242        hname = "DijetMassIn";
00243        FillHist1D(hname,diJetMass ,1.0); 
00244        hname = "DijetMassInWt";
00245        FillHist1D(hname,diJetMass ,weight); 
00246      }
00248       if(index==2 && (abs(etaJet1)>0.7&&abs(etaJet1)<1.3) 
00249                   && (abs(etaJet2)>0.7&&abs(etaJet2)<1.3) ){
00250        hname = "DijetMassOut";
00251        FillHist1D(hname, diJetMass ,1.0); 
00252        hname = "DijetMassOutWt";
00253        FillHist1D(hname,diJetMass ,weight); 
00254      }
00255      if(anaLevel=="Jets")break;  //Jets level for samples without genParticles
00256   
00257 
00258      // Parton Info
00259      edm::Handle<std::vector<reco::GenParticle> > genParticlesHandle_;
00260      evt.getByLabel("genParticles",genParticlesHandle_);
00261      if(debug)for( size_t i = 0; i < genParticlesHandle_->size(); ++ i ) {
00262        const reco::GenParticle & p = (*genParticlesHandle_)[i];
00263        int id = p.pdgId();
00264        int st = p.status();
00265        math::XYZTLorentzVector genP4 = p.p4();
00266        if(i>=2&&i<=8)std::cout << "particle " << i << ": id=" << id << ", status=" << st << ", mass=" << genP4.mass() << ", pt=" <<  genP4.pt() << ", eta=" << genP4.eta() << std::endl; 
00267      }
00268      // Examine the 7th particle in pythia.
00269      // It should be either a resonance (abs(id)>=32) or the first outgoing parton
00270      // for the processes we will consider: dijet resonances, QCD, or QCD +contact interactions.
00271      const reco::GenParticle & p = (*genParticlesHandle_)[6];
00272      int id = p.pdgId();
00273      math::XYZTLorentzVector resonance_p, parton1_p, parton2_p;
00274      if(abs(id)>=32){
00276         resonance_p = p.p4();      
00277         hname = "ResonanceMass";
00278         FillHist1D(hname,resonance_p.mass() ,1.0); 
00279         const reco::GenParticle & q = (*genParticlesHandle_)[7];
00280         parton1_p = q.p4();
00281         const reco::GenParticle & r = (*genParticlesHandle_)[8];
00282         parton2_p = r.p4();
00283         if(debug)std::cout << "Resonance mass=" << resonance_p.mass() << ", parton 1 pt=" << parton1_p.pt()  << ", parton 2 pt=" << parton2_p.pt() << ", diparton mass=" << (parton1_p+parton2_p).mass() << std::endl;
00284      }
00285      else
00286      {
00288         parton1_p = p.p4();
00289         const reco::GenParticle & q = (*genParticlesHandle_)[7];
00290         parton2_p = q.p4();
00291         if(debug)std::cout <<  "parton 1 pt=" << parton1_p.pt()  << ", parton 2 pt=" << parton2_p.pt() << ", diparton mass=" << (parton1_p+parton2_p).mass() << std::endl;
00292      }
00293 
00294       etaPart1 = parton1_p.eta();
00295       etaPart2 = parton2_p.eta();
00296       ptPart1 = parton1_p.pt();
00297       ptPart2 = parton2_p.pt();  
00298       diPartMass = (parton1_p+parton2_p).mass();  
00300      if(abs(etaPart1)<1.3&&abs(etaPart2)<1.3){
00301        hname = "DipartonMass";
00302        FillHist1D(hname,diPartMass ,1.0); 
00303        hname = "DipartonMassWt";
00304        FillHist1D(hname,diPartMass ,weight); 
00305      }
00307      if(abs(etaPart1)<0.7&&abs(etaPart2)<0.7){
00308        hname = "DipartonMassIn";
00309        FillHist1D(hname,diPartMass ,1.0); 
00310        hname = "DipartonMassInWt";
00311        FillHist1D(hname,diPartMass ,weight); 
00312      }
00314      if(    (abs(etaPart1)>0.7&&abs(etaPart1)<1.3)
00315          && (abs(etaPart2)>0.7&&abs(etaPart2)<1.3) ){
00316        hname = "DipartonMassOut";
00317        FillHist1D(hname,diPartMass ,1.0); 
00318        hname = "DipartonMassOutWt";
00319        FillHist1D(hname,diPartMass ,weight); 
00320      }
00321 
00322      // Fill the TTree //
00323      mcTruthTree_->Fill();
00324    
00325      notDone=0;  //We are done, exit the while loop
00326    }//end of while
00327 
00328 }
00330 template<class Jet>
00331 void JetAnaPythia<Jet>::endJob() 
00332 {
00334   if (m_file !=0) 
00335     {
00336       m_file->cd();
00337       mcTruthTree_->Write(); 
00338       for (std::map<TString, TH1*>::iterator hid = m_HistNames1D.begin(); hid != m_HistNames1D.end(); hid++)
00339         hid->second->Write();
00340       delete m_file;
00341       m_file = 0;      
00342     }
00343 }
00345 template<class Jet>
00346 void JetAnaPythia<Jet>::FillHist1D(const TString& histName,const Double_t& value, const Double_t& wt) 
00347 {
00348   std::map<TString, TH1*>::iterator hid=m_HistNames1D.find(histName);
00349   if (hid==m_HistNames1D.end())
00350     std::cout << "%fillHist -- Could not find histogram with name: " << histName << std::endl;
00351   else
00352     hid->second->Fill(value,wt);
00353 }
00355 #include "FWCore/Framework/interface/MakerMacros.h"
00357 typedef JetAnaPythia<CaloJet> CaloJetAnaPythia;
00358 DEFINE_FWK_MODULE(CaloJetAnaPythia);
00360 typedef JetAnaPythia<GenJet> GenJetAnaPythia;
00361 DEFINE_FWK_MODULE(GenJetAnaPythia);
00363 typedef JetAnaPythia<PFJet> PFJetAnaPythia;
00364 DEFINE_FWK_MODULE(PFJetAnaPythia);