CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/PhysicsTools/FWLite/bin/FWLiteHistograms.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 <TH1F.h>
00009 #include <TROOT.h>
00010 #include <TFile.h>
00011 #include <TSystem.h>
00012 
00013 #include "DataFormats/FWLite/interface/Event.h"
00014 #include "DataFormats/Common/interface/Handle.h"
00015 #include "FWCore/FWLite/interface/AutoLibraryLoader.h"
00016 
00017 #include "DataFormats/MuonReco/interface/Muon.h"
00018 #include "DataFormats/PatCandidates/interface/Muon.h"
00019 #include "PhysicsTools/FWLite/interface/TFileService.h"
00020 #include "PhysicsTools/FWLite/interface/CommandLineParser.h"
00021 
00022 
00023 int main(int argc, char* argv[]) 
00024 {
00025   // define what muon you are using; this is necessary as FWLite is not 
00026   // capable of reading edm::Views
00027   using reco::Muon;
00028 
00029   // ----------------------------------------------------------------------
00030   // First Part: 
00031   //
00032   //  * enable the AutoLibraryLoader 
00033   //  * book the histograms of interest 
00034   //  * open the input file
00035   // ----------------------------------------------------------------------
00036 
00037   // load framework libraries
00038   gSystem->Load( "libFWCoreFWLite" );
00039   AutoLibraryLoader::enable();
00040 
00041   // initialize command line parser
00042   optutl::CommandLineParser parser ("Analyze FWLite Histograms");
00043 
00044   // set defaults
00045   parser.integerValue ("maxEvents"  ) = 1000;
00046   parser.integerValue ("outputEvery") =   10;
00047   parser.stringValue  ("outputFile" ) = "analyzeFWLiteHistograms.root";
00048 
00049   // parse arguments
00050   parser.parseArguments (argc, argv);
00051   int maxEvents_ = parser.integerValue("maxEvents");
00052   unsigned int outputEvery_ = parser.integerValue("outputEvery");
00053   std::string outputFile_ = parser.stringValue("outputFile");
00054   std::vector<std::string> inputFiles_ = parser.stringVector("inputFiles");
00055 
00056   // book a set of histograms
00057   fwlite::TFileService fs = fwlite::TFileService(outputFile_.c_str());
00058   TFileDirectory dir = fs.mkdir("analyzeBasicPat");
00059   TH1F* muonPt_  = dir.make<TH1F>("muonPt"  , "pt"  ,   100,   0., 300.);
00060   TH1F* muonEta_ = dir.make<TH1F>("muonEta" , "eta" ,   100,  -3.,   3.);
00061   TH1F* muonPhi_ = dir.make<TH1F>("muonPhi" , "phi" ,   100,  -5.,   5.);  
00062   TH1F* mumuMass_= dir.make<TH1F>("mumuMass", "mass",    90,  30.,  120.);
00063 
00064   // loop the events
00065   int ievt=0;  
00066   for(unsigned int iFile=0; iFile<inputFiles_.size(); ++iFile){
00067     // open input file (can be located on castor)
00068     TFile* inFile = TFile::Open(inputFiles_[iFile].c_str());
00069     if( inFile ){
00070       // ----------------------------------------------------------------------
00071       // Second Part: 
00072       //
00073       //  * loop the events in the input file 
00074       //  * receive the collections of interest via fwlite::Handle
00075       //  * fill the histograms
00076       //  * after the loop close the input file
00077       // ----------------------------------------------------------------------      
00078       fwlite::Event ev(inFile);
00079       for(ev.toBegin(); !ev.atEnd(); ++ev, ++ievt){
00080         edm::EventBase const & event = ev;
00081         // break loop if maximal number of events is reached 
00082         if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00083         // simple event counter
00084         if(outputEvery_!=0 ? (ievt>0 && ievt%outputEvery_==0) : false) 
00085           std::cout << "  processing event: " << ievt << std::endl;
00086 
00087         // Handle to the muon collection
00088         edm::Handle<std::vector<Muon> > muons;
00089         event.getByLabel(std::string("muons"), muons);
00090         
00091         // loop muon collection and fill histograms
00092         for(std::vector<Muon>::const_iterator mu1=muons->begin(); mu1!=muons->end(); ++mu1){
00093           muonPt_ ->Fill( mu1->pt () );
00094           muonEta_->Fill( mu1->eta() );
00095           muonPhi_->Fill( mu1->phi() );   
00096           if( mu1->pt()>20 && fabs(mu1->eta())<2.1 ){
00097             for(std::vector<Muon>::const_iterator mu2=muons->begin(); mu2!=muons->end(); ++mu2){
00098               if(mu2>mu1){ // prevent double conting
00099                 if( mu1->charge()*mu2->charge()<0 ){ // check only muon pairs of unequal charge 
00100                   if( mu2->pt()>20 && fabs(mu2->eta())<2.1 ){
00101                     mumuMass_->Fill( (mu1->p4()+mu2->p4()).mass() );
00102                   }
00103                 }
00104               }
00105             }
00106           }
00107         }
00108       }  
00109       // close input file
00110       inFile->Close();
00111     }
00112     // break loop if maximal number of events is reached:
00113     // this has to be done twice to stop the file loop as well
00114     if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00115   }
00116   return 0;
00117 }