CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/PhysicsTools/FWLite/bin/FWLiteWithSelectorUtils.cc

Go to the documentation of this file.
00001 #include <TH1F.h>
00002 #include <TROOT.h>
00003 #include <TFile.h>
00004 #include <TSystem.h>
00005 
00006 #include "DataFormats/FWLite/interface/Event.h"
00007 #include "DataFormats/Common/interface/Handle.h"
00008 #include "FWCore/FWLite/interface/AutoLibraryLoader.h"
00009 
00010 #include "DataFormats/FWLite/interface/InputSource.h"
00011 #include "DataFormats/FWLite/interface/OutputFiles.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "FWCore/PythonParameterSet/interface/MakeParameterSets.h"
00014 
00015 #include "TStopwatch.h"
00016 #include "PhysicsTools/FWLite/interface/WSelector.h"
00017 #include "PhysicsTools/FWLite/interface/TFileService.h"
00018 //#include "PhysicsTools/FWLite/interface/WSelectorFast.h"
00019 
00020 
00021 int main(int argc, char* argv[]) 
00022 {
00023   // define what muon you are using; this is necessary as FWLite is not 
00024   // capable of reading edm::Views
00025   using reco::Muon;
00026 
00027   // ----------------------------------------------------------------------
00028   // First Part: 
00029   //
00030   //  * enable the AutoLibraryLoader 
00031   //  * book the histograms of interest 
00032   //  * open the input file
00033   // ----------------------------------------------------------------------
00034 
00035   // load framework libraries
00036   gSystem->Load( "libFWCoreFWLite" );
00037   AutoLibraryLoader::enable();
00038 
00039   if ( argc < 2 ) {
00040     std::cout << "Usage : " << argv[0] << " [parameters.py]" << std::endl;
00041     return 0;
00042   }
00043 
00044   if( !edm::readPSetsFrom(argv[1])->existsAs<edm::ParameterSet>("process") ){
00045     std::cout << " ERROR: ParametersSet 'process' is missing in your configuration file" << std::endl; exit(0);
00046   }
00047   // get the python configuration
00048   const edm::ParameterSet& process = edm::readPSetsFrom(argv[1])->getParameter<edm::ParameterSet>("process");
00049   fwlite::InputSource inputHandler_(process); fwlite::OutputFiles outputHandler_(process);
00050 
00051   // initialize the W selector
00052   edm::ParameterSet selection = process.getParameter<edm::ParameterSet>("selection");
00053   WSelector wSelector( selection ); pat::strbitset wSelectorReturns = wSelector.getBitTemplate();
00054   
00055   // book a set of histograms
00056   fwlite::TFileService fs = fwlite::TFileService(outputHandler_.file().c_str());
00057   TFileDirectory theDir = fs.mkdir("analyzeBasicPat");
00058   TH1F* muonPt_  = theDir.make<TH1F>("muonPt", "pt",    100,  0.,300.);
00059   TH1F* muonEta_ = theDir.make<TH1F>("muonEta","eta",   100, -3.,  3.);
00060   TH1F* muonPhi_ = theDir.make<TH1F>("muonPhi","phi",   100, -5.,  5.);  
00061 
00062   // start a CPU timer
00063   TStopwatch timer; timer.Start();
00064 
00065   // loop the events
00066   int ievt=0;  
00067   unsigned int nEventsAnalyzed = 0;
00068   int maxEvents_( inputHandler_.maxEvents() );
00069   for(unsigned int iFile=0; iFile<inputHandler_.files().size(); ++iFile){
00070     // open input file (can be located on castor)
00071     TFile* inFile = TFile::Open(inputHandler_.files()[iFile].c_str());
00072     if( inFile ){
00073       // ----------------------------------------------------------------------
00074       // Second Part: 
00075       //
00076       //  * loop the events in the input file 
00077       //  * receive the collections of interest via fwlite::Handle
00078       //  * fill the histograms
00079       //  * after the loop close the input file
00080       // ----------------------------------------------------------------------
00081       fwlite::Event ev(inFile);
00082       for(ev.toBegin(); !ev.atEnd(); ++ev, ++ievt){
00083         edm::EventBase const & event = ev;
00084         // break loop if maximal number of events is reached 
00085         if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00086         // simple event counter
00087         if(inputHandler_.reportAfter()!=0 ? (ievt>0 && ievt%inputHandler_.reportAfter()==0) : false) 
00088           std::cout << "  processing event: " << ievt << std::endl;
00089     
00090         if ( wSelector(event, wSelectorReturns ) ) {
00091           pat::Muon const & wMuon = wSelector.wMuon();
00092           muonPt_ ->Fill( wMuon.pt()  );
00093           muonEta_->Fill( wMuon.eta() );
00094           muonPhi_->Fill( wMuon.phi() );
00095         } 
00096         ++nEventsAnalyzed;
00097       }  
00098       // close input file
00099       inFile->Close();
00100     }
00101     // break loop if maximal number of events is reached:
00102     // this has to be done twice to stop the file loop as well
00103     if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00104   }
00105   // stop CPU timer
00106   timer.Stop();
00107 
00108   // print selector
00109   wSelector.print(std::cout);
00110 
00111   // print some timing statistics
00112   double rtime = timer.RealTime();
00113   double ctime = timer.CpuTime ();
00114   // timing printouts
00115   printf("Analyzed events: %d \n",nEventsAnalyzed);
00116   printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
00117   printf("%4.2f events / RealTime second .\n", (double)nEventsAnalyzed/rtime);
00118   printf("%4.2f events / CpuTime second .\n", (double)nEventsAnalyzed/ctime);
00119   return 0;
00120 }