CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/PhysicsTools/PatExamples/bin/PatMcMatchFWLiteAnalyzer.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include <string>
00003 #include <vector>
00004 #include <sstream>
00005 #include <fstream>
00006 #include <iostream>
00007 
00008 #include <TH2F.h> // Use the correct histograms
00009 #include <TROOT.h>
00010 #include <TFile.h>
00011 #include <TSystem.h>
00012 
00013 #include "DataFormats/FWLite/interface/Handle.h"
00014 #include "DataFormats/PatCandidates/interface/Muon.h" // Include the examined data formats
00015 #include "DataFormats/PatCandidates/interface/Jet.h"  // Include the examined data formats
00016 #include "FWCore/FWLite/interface/AutoLibraryLoader.h"
00017 
00018 
00019 int main(int argc, char* argv[]) 
00020 {
00021   // ----------------------------------------------------------------------
00022   // First Part: 
00023   //
00024   //  * enable the AutoLibraryLoader 
00025   //  * book the histograms of interest 
00026   //  * open the input file
00027   // ----------------------------------------------------------------------
00028 
00029   // load framework libraries
00030   gSystem->Load( "libFWCoreFWLite" );
00031   AutoLibraryLoader::enable();
00032   
00033   // book a set of histograms
00034   TH2F* muonPt_  = new TH2F( "muonPt", "Muon Pt", 60, 0., 300., 60, 0., 300. ); // 2-D histo for muon Pt
00035   muonPt_->SetXTitle( "gen." );
00036   muonPt_->SetYTitle( "reco." );
00037   TH2F* jetPt_   = new TH2F( "jetPt", "Jet Pt", 100, 0., 500., 100, 0., 500. ); // 2-D histo for jet Pt
00038   jetPt_->SetXTitle( "gen." );
00039   jetPt_->SetYTitle( "reco." );
00040   
00041   // open input file (can be located on castor)
00042   TFile* inFile = TFile::Open( "file:edmPatMcMatch.root" ); // Adapt the input file name
00043 
00044   // ----------------------------------------------------------------------
00045   // Second Part: 
00046   //
00047   //  * loop the events in the input file 
00048   //  * receive the collections of interest via fwlite::Handle
00049   //  * fill the histograms
00050   //  * after the loop close the input file
00051   // ----------------------------------------------------------------------
00052 
00053   // loop the events
00054   unsigned int iEvent=0;
00055   fwlite::Event event(inFile);
00056   for(event.toBegin(); !event.atEnd(); ++event, ++iEvent){
00057     // break loop after end of file is reached 
00058     // or after 1000 events have been processed
00059     if( iEvent==1000 ) break;
00060     
00061     // simple event counter
00062     if(iEvent>0 && iEvent%25==0){ // Reduce print-out
00063       std::cout << "  processing event: " << iEvent << std::endl;
00064     }
00065 
00066     // fwlite::Handle to the muon collection
00067     fwlite::Handle<std::vector<pat::Muon> > muons; // Access the muon collection
00068     muons.getByLabel(event, "cleanLayer1Muons");
00069     // fwlite::Handle to the muon collection
00070     fwlite::Handle<std::vector<pat::Jet> > jets; // Access the jet collection
00071     jets.getByLabel(event, "cleanLayer1Jets");
00072     
00073     // loop muon collection and fill histograms
00074     for(unsigned i=0; i<muons->size(); ++i){
00075       const reco::GenParticle * genMuon = (*muons)[i].genParticle(); // Get the matched generator muon
00076       if ( genMuon ) {                                               // Check for valid reference
00077         muonPt_->Fill( genMuon->pt(), (*muons)[i].pt() );            // Fill 2-D histo
00078       }
00079     }
00080     // loop jet collection and fill histograms
00081     for(unsigned i=0; i<jets->size(); ++i){
00082       const reco::GenJet * genJet = (*jets)[i].genJet(); // Get the matched generator jet
00083       if ( genJet ) {                                    // Check for valid reference
00084         jetPt_->Fill( genJet->pt(), (*jets)[i].pt() );   // Fill 2-D histo
00085       }
00086     }
00087   }  
00088   // close input file
00089   inFile->Close();
00090 
00091   // ----------------------------------------------------------------------
00092   // Third Part: 
00093   //
00094   //  * open the output file 
00095   //  * write the histograms to the output file
00096   //  * close the output file
00097   // ----------------------------------------------------------------------
00098   
00099   //open output file
00100   TFile outFile( "rootPatMcMatch.root", "recreate" ); // Adapt the output file name
00101   outFile.mkdir("analyzeMcMatchPat");                 // Adapt output file according to modifications
00102   outFile.cd("analyzeMcMatchPat");
00103   muonPt_->Write( );
00104   jetPt_->Write( );
00105   outFile.Close();
00106   
00107   // ----------------------------------------------------------------------
00108   // Fourth Part: 
00109   //
00110   //  * never forgett to free the memory of the histograms
00111   // ----------------------------------------------------------------------
00112 
00113   // free allocated space
00114   delete muonPt_; // Delete the muon histo
00115   delete jetPt_;  // Delete the jet histo
00116   
00117   // that's it!
00118   return 0;
00119 }