Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:00

0001 #include <memory>
0002 #include <string>
0003 #include <vector>
0004 #include <sstream>
0005 #include <fstream>
0006 #include <iostream>
0007 
0008 #include <TH1F.h>
0009 #include <TROOT.h>
0010 #include <TFile.h>
0011 #include <TSystem.h>
0012 
0013 #include "DataFormats/FWLite/interface/Event.h"
0014 #include "DataFormats/Common/interface/Handle.h"
0015 #include "FWCore/FWLite/interface/FWLiteEnabler.h"
0016 
0017 #include "DataFormats/MuonReco/interface/Muon.h"
0018 #include "DataFormats/PatCandidates/interface/Muon.h"
0019 #include "PhysicsTools/FWLite/interface/TFileService.h"
0020 #include "PhysicsTools/FWLite/interface/CommandLineParser.h"
0021 
0022 int main(int argc, char* argv[]) {
0023   // ----------------------------------------------------------------------
0024   // First Part:
0025   //
0026   //  * enable FWLite
0027   //  * book the histograms of interest
0028   //  * open the input file
0029   // ----------------------------------------------------------------------
0030 
0031   // load framework libraries
0032   gSystem->Load("libFWCoreFWLite");
0033   FWLiteEnabler::enable();
0034 
0035   // initialize command line parser
0036   optutl::CommandLineParser parser("Analyze FWLite Histograms");
0037 
0038   // set defaults
0039   parser.integerValue("maxEvents") = 1000;
0040   parser.integerValue("outputEvery") = 10;
0041   parser.stringValue("outputFile") = "analyzeEdmTuple.root";
0042 
0043   // parse arguments
0044   parser.parseArguments(argc, argv);
0045   int maxEvents_ = parser.integerValue("maxEvents");
0046   unsigned int outputEvery_ = parser.integerValue("outputEvery");
0047   std::string outputFile_ = parser.stringValue("outputFile");
0048   std::vector<std::string> inputFiles_ = parser.stringVector("inputFiles");
0049 
0050   // book a set of histograms
0051   fwlite::TFileService fs = fwlite::TFileService(outputFile_);
0052   TFileDirectory dir = fs.mkdir("analyzePatMuon");
0053   TH1F* muonPt_ = dir.make<TH1F>("muonPt", "pt", 100, 0., 300.);
0054   TH1F* muonEta_ = dir.make<TH1F>("muonEta", "eta", 100, -3., 3.);
0055   TH1F* muonPhi_ = dir.make<TH1F>("muonPhi", "phi", 100, -5., 5.);
0056 
0057   // loop the events
0058   int ievt = 0;
0059   for (unsigned int iFile = 0; iFile < inputFiles_.size(); ++iFile) {
0060     // open input file (can be located on castor)
0061     TFile* inFile = TFile::Open(inputFiles_[iFile].c_str());
0062     if (inFile) {
0063       // ----------------------------------------------------------------------
0064       // Second Part:
0065       //
0066       //  * loop the events in the input file
0067       //  * receive the collections of interest via fwlite::Handle
0068       //  * fill the histograms
0069       //  * after the loop close the input file
0070       // ----------------------------------------------------------------------
0071       fwlite::Event ev(inFile);
0072       for (ev.toBegin(); !ev.atEnd(); ++ev, ++ievt) {
0073         edm::EventBase const& event = ev;
0074         // break loop if maximal number of events is reached
0075         if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false)
0076           break;
0077         // simple event counter
0078         if (outputEvery_ != 0 ? (ievt > 0 && ievt % outputEvery_ == 0) : false)
0079           std::cout << "  processing event: " << ievt << std::endl;
0080 
0081         // Handle to the muon pt
0082         edm::Handle<std::vector<float> > muonPt;
0083         event.getByLabel(std::string("patMuonAnalyzer:pt"), muonPt);
0084         // loop muon collection and fill histograms
0085         for (std::vector<float>::const_iterator mu1 = muonPt->begin(); mu1 != muonPt->end(); ++mu1) {
0086           muonPt_->Fill(*mu1);
0087         }
0088         // Handle to the muon eta
0089         edm::Handle<std::vector<float> > muonEta;
0090         event.getByLabel(std::string("patMuonAnalyzer:eta"), muonEta);
0091         for (std::vector<float>::const_iterator mu1 = muonEta->begin(); mu1 != muonEta->end(); ++mu1) {
0092           muonEta_->Fill(*mu1);
0093         }
0094         // Handle to the muon phi
0095         edm::Handle<std::vector<float> > muonPhi;
0096         event.getByLabel(std::string("patMuonAnalyzer:phi"), muonPhi);
0097         for (std::vector<float>::const_iterator mu1 = muonPhi->begin(); mu1 != muonPhi->end(); ++mu1) {
0098           muonPhi_->Fill(*mu1);
0099         }
0100       }
0101       // close input file
0102       inFile->Close();
0103     }
0104     // break loop if maximal number of events is reached:
0105     // this has to be done twice to stop the file loop as well
0106     if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false)
0107       break;
0108   }
0109   return 0;
0110 }