CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/PhysicsTools/PatExamples/bin/PatCleaningExercise.cc

Go to the documentation of this file.
00001 #include <vector>
00002 
00003 #include "TH1.h"
00004 #include "TFile.h"
00005 #include <TROOT.h>
00006 #include <TSystem.h>
00007 
00008 #include "DataFormats/Math/interface/deltaR.h"
00009 #include "DataFormats/FWLite/interface/Event.h"
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "DataFormats/PatCandidates/interface/Jet.h"
00012 #include "FWCore/ParameterSet/interface/ProcessDesc.h"
00013 #include "FWCore/FWLite/interface/AutoLibraryLoader.h"
00014 #include "PhysicsTools/FWLite/interface/TFileService.h"
00015 #include "FWCore/PythonParameterSet/interface/PythonProcessDesc.h"
00016 
00017 
00018 //using namespace std;
00019 //using namespace reco;
00020 //using namespace pat;
00021 
00022 
00023 int main(int argc, char *argv[]){
00024   // ----------------------------------------------------------------------
00025   // First Part: 
00026   //
00027   //  * enable the AutoLibraryLoader 
00028   //  * book the histograms of interest 
00029   //  * open the input file
00030   // ----------------------------------------------------------------------
00031 
00032   // load framework libraries
00033   gSystem->Load( "libFWCoreFWLite" );
00034   AutoLibraryLoader::enable();
00035 
00036   // only allow one argument for this simple example which should be the
00037   // the python cfg file
00038   if ( argc < 2 ) {
00039     std::cout << "Usage : " << argv[0] << " [parameters.py]" << std::endl;
00040     return 0;
00041   }
00042 
00043   // get the python configuration
00044   PythonProcessDesc builder(argv[1]);
00045   const edm::ParameterSet& fwliteParameters = builder.processDesc()->getProcessPSet()->getParameter<edm::ParameterSet>("FWLiteParams");
00046 
00047   // now get each parameter
00048   std::string   input_   ( fwliteParameters.getParameter<std::string  >("inputFile"  ) );
00049   std::string   output_  ( fwliteParameters.getParameter<std::string  >("outputFile" ) );
00050   std::string   overlaps_( fwliteParameters.getParameter<std::string  >("overlaps") );
00051   edm::InputTag jets_    ( fwliteParameters.getParameter<edm::InputTag>("jets"  ) );
00052 
00053   // book a set of histograms
00054   fwlite::TFileService fs = fwlite::TFileService(output_.c_str());
00055   TFileDirectory theDir = fs.mkdir("analyzePatCleaning");
00056   TH1F* emfAllJets_    = theDir.make<TH1F>("emfAllJets"    , "f_{emf}(All Jets)"    ,  20,  0.,  1.);
00057   TH1F* emfCleanJets_  = theDir.make<TH1F>("emfCleanJets"  , "f_{emf}(Clean Jets)"  ,  20,  0.,  1.);
00058   TH1F* emfOverlapJets_= theDir.make<TH1F>("emfOverlapJets", "f_{emf}(Overlap Jets)",  20,  0.,  1.);
00059   TH1F* deltaRElecJet_ = theDir.make<TH1F>("deltaRElecJet" , "#DeltaR (elec, jet)"  ,  10,  0., 0.5);
00060   TH1F* elecOverJet_   = theDir.make<TH1F>("elecOverJet"   , "E_{elec}/E_{jet}"     , 100,  0.,  2.);
00061   TH1F* nOverlaps_     = theDir.make<TH1F>("nOverlaps"     , "Number of overlaps"   ,   5,  0.,  5.);
00062   
00063   // open input file (can be located on castor)
00064   TFile* inFile = TFile::Open(input_.c_str());
00065 
00066   // ----------------------------------------------------------------------
00067   // Second Part: 
00068   //
00069   //  * loop the events in the input file 
00070   //  * receive the collections of interest via fwlite::Handle
00071   //  * fill the histograms
00072   //  * after the loop close the input file
00073   // ----------------------------------------------------------------------
00074 
00075   // loop the events
00076   unsigned int iEvent=0;
00077   fwlite::Event ev(inFile); 
00078   for(ev.toBegin(); !ev.atEnd(); ++ev, ++iEvent){
00079     edm::EventBase const & event = ev;
00080 
00081     // break loop after end of file is reached 
00082     // or after 1000 events have been processed
00083     if( iEvent==1000 ) break;
00084 
00085     // simple event counter
00086     if(iEvent>0 && iEvent%1==0){
00087       std::cout << "  processing event: " << iEvent << std::endl;
00088     }
00089 
00090     // handle to jet collection
00091     edm::Handle<std::vector<pat::Jet> > jets;
00092     event.getByLabel(jets_, jets);
00093 
00094     // loop over the jets in the event
00095     for( std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); jet++ ){
00096       if(jet->pt()>20 && jet==jets->begin()){
00097         emfAllJets_->Fill( jet->emEnergyFraction() );
00098         if(! jet->hasOverlaps(overlaps_)){
00099           emfCleanJets_->Fill( jet->emEnergyFraction() );
00100         }
00101         else{
00102           //get all overlaps
00103           const reco::CandidatePtrVector overlaps = jet->overlaps(overlaps_);
00104           nOverlaps_->Fill( overlaps.size() );
00105           emfOverlapJets_->Fill( jet->emEnergyFraction() );
00106           //loop over the overlaps
00107           for( reco::CandidatePtrVector::const_iterator overlap = overlaps.begin(); overlap != overlaps.end(); overlap++){ 
00108             float deltaR = reco::deltaR( (*overlap)->eta(), (*overlap)->phi(), jet->eta(), jet->phi() );
00109             deltaRElecJet_->Fill( deltaR );
00110             elecOverJet_->Fill( (*overlap)->energy()/jet->energy() );
00111           }
00112         }
00113       }
00114     }
00115   }
00116   inFile->Close();
00117   return 0;
00118 }