CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/JetMETCorrections/MCJet/plugins/PFMCTruthTreeProducer.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <sstream>
00003 #include <istream>
00004 #include <fstream>
00005 #include <iomanip>
00006 #include <string>
00007 #include <cmath>
00008 #include <functional>
00009 
00010 #include "JetMETCorrections/MCJet/plugins/PFMCTruthTreeProducer.h" 
00011 #include "JetMETCorrections/MCJet/plugins/JetUtilMC.h"
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "DataFormats/Common/interface/Handle.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00018 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00019 #include "DataFormats/JetReco/interface/PFJet.h"
00020 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00021 #include "DataFormats/JetReco/interface/GenJet.h"
00022 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00023 using namespace edm;
00024 using namespace reco;
00025 using namespace std;
00026 //namespace cms
00027 //{
00028 
00029 PFMCTruthTreeProducer::PFMCTruthTreeProducer(edm::ParameterSet const& cfg) 
00030 {
00031   jets_          = cfg.getParameter<std::string> ("jets");
00032   genjets_       = cfg.getParameter<std::string> ("genjets");
00033   histogramFile_ = cfg.getParameter<std::string> ("histogramFile");
00034 }
00036 void PFMCTruthTreeProducer::beginJob() 
00037 {
00038   file_          = new TFile(histogramFile_.c_str(),"RECREATE");
00039   mcTruthTree_   = new TTree("mcTruthTree","mcTruthTree");
00040 
00041   mcTruthTree_->Branch("ptJet",      &ptJet_,      "ptJet_/F");
00042   mcTruthTree_->Branch("ptGen",      &ptGen_,      "ptGen_/F");
00043   mcTruthTree_->Branch("ptHat",      &ptHat_,      "ptHat_/F");
00044   mcTruthTree_->Branch("chfJet",     &chfJet_,     "chfJet_/F");
00045   mcTruthTree_->Branch("nhfJet",     &nhfJet_,     "nhfJet_/F");
00046   mcTruthTree_->Branch("cemfJet",    &cemfJet_,    "cemfJet_/F");
00047   mcTruthTree_->Branch("nemfJet",    &nemfJet_,    "nemfJet_/F");
00048   mcTruthTree_->Branch("cmultiJet",  &cmultiJet_,  "cmultiJet_/I");
00049   mcTruthTree_->Branch("nmultiJet",  &nmultiJet_,  "nmultiJet_/I");
00050   mcTruthTree_->Branch("etaJet",     &etaJet_,    "etaJet_/F");
00051   mcTruthTree_->Branch("etaGen",     &etaGen_,     "etaGen_/F");
00052   mcTruthTree_->Branch("phiJet",     &phiJet_,     "phiJet_/F");
00053   mcTruthTree_->Branch("phiGen",     &phiGen_,     "phiGen_/F");
00054   mcTruthTree_->Branch("dR",         &dR_,         "dR_/F");
00055   mcTruthTree_->Branch("rank",       &rank_,       "rank_/I");
00056 }
00058 void PFMCTruthTreeProducer::endJob() 
00059 {
00060   if (file_ !=0) 
00061     {
00062       file_->cd();
00063       mcTruthTree_->Write();
00064     }
00065   file_ = 0;
00066 }
00068 void PFMCTruthTreeProducer::analyze(edm::Event const& event, edm::EventSetup const& iSetup) 
00069 { 
00070   edm::Handle<GenJetCollection> genjets;
00071   edm::Handle<PFJetCollection> jets;
00072   edm::Handle<GenEventInfoProduct> hEventInfo;
00073   PFJetCollection::const_iterator i_jet,i_matched;
00074   GenJetCollection::const_iterator i_genjet;
00075   event.getByLabel (genjets_,genjets);
00076   event.getByLabel (jets_,jets);
00077   event.getByLabel("generator",hEventInfo);
00078   ptHat_ = hEventInfo->binningValues()[0];
00079   float rr;  
00080   int njet(0);
00081   if (jets->size()>0 && genjets->size()>0)
00082     {
00083       for (i_genjet = genjets->begin(); i_genjet != genjets->end(); i_genjet++)
00084        {    
00085          float rmin(99);
00086          for(i_jet = jets->begin();i_jet != jets->end(); i_jet++)
00087            {
00088              rr=radius(i_genjet,i_jet);
00089              if (rr<rmin)
00090                {
00091                  rmin = rr;
00092                  i_matched = i_jet;
00093                }
00094            }
00095          ptGen_     = i_genjet->pt();
00096          etaGen_    = i_genjet->eta();
00097          phiGen_    = i_genjet->phi();
00098          ptJet_     = i_matched->pt();
00099          etaJet_    = i_matched->eta();
00100          phiJet_    = i_matched->phi();
00101          chfJet_    = i_matched->chargedHadronEnergyFraction();
00102          nhfJet_    = i_matched->neutralHadronEnergyFraction();
00103          cemfJet_   = i_matched->chargedEmEnergyFraction();
00104          nemfJet_   = i_matched->neutralEmEnergyFraction();
00105          cmultiJet_ = i_matched->chargedMultiplicity();
00106          nmultiJet_ = i_matched->neutralMultiplicity();
00107          dR_        = rmin;
00108          rank_      = njet;     
00109          mcTruthTree_->Fill();
00110          njet++;
00111        }  
00112     }      
00113 }
00115 PFMCTruthTreeProducer::~PFMCTruthTreeProducer() 
00116 {
00117   delete file_;
00118   delete mcTruthTree_;
00119 }
00120 //}