CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/PhysicsTools/FWLite/bin/FWLiteWithPythonConfig.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 "DataFormats/MuonReco/interface/Muon.h"
00016 #include "DataFormats/PatCandidates/interface/Muon.h"
00017 #include "PhysicsTools/FWLite/interface/TFileService.h"
00018 
00019 
00020 int main(int argc, char* argv[]) 
00021 {
00022   // define what muon you are using; this is necessary as FWLite is not 
00023   // capable of reading edm::Views
00024   using reco::Muon;
00025 
00026   // ----------------------------------------------------------------------
00027   // First Part: 
00028   //
00029   //  * enable the AutoLibraryLoader 
00030   //  * book the histograms of interest 
00031   //  * open the input file
00032   // ----------------------------------------------------------------------
00033 
00034   // load framework libraries
00035   gSystem->Load( "libFWCoreFWLite" );
00036   AutoLibraryLoader::enable();
00037 
00038   // parse arguments
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 
00052   // now get each parameter
00053   const edm::ParameterSet& ana = process.getParameter<edm::ParameterSet>("muonAnalyzer");
00054   edm::InputTag muons_( ana.getParameter<edm::InputTag>("muons") );
00055 
00056   // book a set of histograms
00057   fwlite::TFileService fs = fwlite::TFileService(outputHandler_.file().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   int maxEvents_( inputHandler_.maxEvents() );
00067   for(unsigned int iFile=0; iFile<inputHandler_.files().size(); ++iFile){
00068     // open input file (can be located on castor)
00069     TFile* inFile = TFile::Open(inputHandler_.files()[iFile].c_str());
00070     if( inFile ){
00071       // ----------------------------------------------------------------------
00072       // Second Part: 
00073       //
00074       //  * loop the events in the input file 
00075       //  * receive the collections of interest via fwlite::Handle
00076       //  * fill the histograms
00077       //  * after the loop close the input file
00078       // ----------------------------------------------------------------------
00079       fwlite::Event ev(inFile);
00080       for(ev.toBegin(); !ev.atEnd(); ++ev, ++ievt){
00081         edm::EventBase const & event = ev;
00082         // break loop if maximal number of events is reached 
00083         if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00084         // simple event counter
00085         if(inputHandler_.reportAfter()!=0 ? (ievt>0 && ievt%inputHandler_.reportAfter()==0) : false) 
00086           std::cout << "  processing event: " << ievt << std::endl;
00087         
00088         // Handle to the muon collection
00089         edm::Handle<std::vector<Muon> > muons;
00090         event.getByLabel(muons_, muons);
00091         
00092         // loop muon collection and fill histograms
00093         for(std::vector<Muon>::const_iterator mu1=muons->begin(); mu1!=muons->end(); ++mu1){
00094           muonPt_ ->Fill( mu1->pt () );
00095           muonEta_->Fill( mu1->eta() );
00096           muonPhi_->Fill( mu1->phi() );   
00097           if( mu1->pt()>20 && fabs(mu1->eta())<2.1 ){
00098             for(std::vector<Muon>::const_iterator mu2=muons->begin(); mu2!=muons->end(); ++mu2){
00099               if(mu2>mu1){ // prevent double conting
00100                 if( mu1->charge()*mu2->charge()<0 ){ // check only muon pairs of unequal charge 
00101                   if( mu2->pt()>20 && fabs(mu2->eta())<2.1 ){
00102                     mumuMass_->Fill( (mu1->p4()+mu2->p4()).mass() );
00103                   }
00104                 }
00105               }
00106             }
00107           }
00108         }
00109       }  
00110       // close input file
00111       inFile->Close();
00112     }
00113     // break loop if maximal number of events is reached:
00114     // this has to be done twice to stop the file loop as well
00115     if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00116   }
00117   return 0;
00118 }