CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/PhysicsTools/PatExamples/bin/PatMuonEDMAnalyzer.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 int main(int argc, char* argv[]) 
00023 {
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   // initialize command line parser
00037   optutl::CommandLineParser parser ("Analyze FWLite Histograms");
00038 
00039   // set defaults
00040   parser.integerValue ("maxEvents"  ) = 1000;
00041   parser.integerValue ("outputEvery") =   10;
00042   parser.stringValue  ("outputFile" ) = "analyzeEdmTuple.root";
00043 
00044   // parse arguments
00045   parser.parseArguments (argc, argv);
00046   int maxEvents_ = parser.integerValue("maxEvents");
00047   unsigned int outputEvery_ = parser.integerValue("outputEvery");
00048   std::string outputFile_ = parser.stringValue("outputFile");
00049   std::vector<std::string> inputFiles_ = parser.stringVector("inputFiles");
00050 
00051   // book a set of histograms
00052   fwlite::TFileService fs = fwlite::TFileService(outputFile_.c_str());
00053   TFileDirectory dir = fs.mkdir("analyzePatMuon");
00054   TH1F* muonPt_  = dir.make<TH1F>("muonPt"  , "pt"  ,   100,   0., 300.);
00055   TH1F* muonEta_ = dir.make<TH1F>("muonEta" , "eta" ,   100,  -3.,   3.);
00056   TH1F* muonPhi_ = dir.make<TH1F>("muonPhi" , "phi" ,   100,  -5.,   5.);  
00057 
00058   // loop the events
00059   int ievt=0;  
00060   for(unsigned int iFile=0; iFile<inputFiles_.size(); ++iFile){
00061     // open input file (can be located on castor)
00062     TFile* inFile = TFile::Open(inputFiles_[iFile].c_str());
00063     if( inFile ){
00064       // ----------------------------------------------------------------------
00065       // Second Part: 
00066       //
00067       //  * loop the events in the input file 
00068       //  * receive the collections of interest via fwlite::Handle
00069       //  * fill the histograms
00070       //  * after the loop close the input file
00071       // ----------------------------------------------------------------------      
00072       fwlite::Event ev(inFile);
00073       for(ev.toBegin(); !ev.atEnd(); ++ev, ++ievt){
00074         edm::EventBase const & event = ev;
00075         // break loop if maximal number of events is reached 
00076         if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00077         // simple event counter
00078         if(outputEvery_!=0 ? (ievt>0 && ievt%outputEvery_==0) : false) 
00079           std::cout << "  processing event: " << ievt << std::endl;
00080 
00081         // Handle to the muon pt
00082         edm::Handle<std::vector<float> > muonPt;
00083         event.getByLabel(std::string("patMuonAnalyzer:pt"), muonPt);
00084         // loop muon collection and fill histograms
00085         for(std::vector<float>::const_iterator mu1=muonPt->begin(); mu1!=muonPt->end(); ++mu1){
00086           muonPt_ ->Fill( *mu1 );
00087         }
00088         // Handle to the muon eta
00089         edm::Handle<std::vector<float> > muonEta;
00090         event.getByLabel(std::string("patMuonAnalyzer:eta"), muonEta);
00091         for(std::vector<float>::const_iterator mu1=muonEta->begin(); mu1!=muonEta->end(); ++mu1){
00092           muonEta_ ->Fill( *mu1 );
00093         }
00094         // Handle to the muon phi
00095         edm::Handle<std::vector<float> > muonPhi;
00096         event.getByLabel(std::string("patMuonAnalyzer:phi"), muonPhi);
00097         for(std::vector<float>::const_iterator mu1=muonPhi->begin(); mu1!=muonPhi->end(); ++mu1){
00098           muonPhi_ ->Fill( *mu1 );
00099         }
00100       }  
00101       // close input file
00102       inFile->Close();
00103     }
00104     // break loop if maximal number of events is reached:
00105     // this has to be done twice to stop the file loop as well
00106     if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
00107   }
00108   return 0;
00109 }
00110 
00111