CMS 3D CMS Logo

PatMuonEDMAnalyzer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 #include <vector>
4 #include <sstream>
5 #include <fstream>
6 #include <iostream>
7 
8 #include <TH1F.h>
9 #include <TROOT.h>
10 #include <TFile.h>
11 #include <TSystem.h>
12 
16 
21 
22 int main(int argc, char* argv[]) {
23  // ----------------------------------------------------------------------
24  // First Part:
25  //
26  // * enable FWLite
27  // * book the histograms of interest
28  // * open the input file
29  // ----------------------------------------------------------------------
30 
31  // load framework libraries
32  gSystem->Load("libFWCoreFWLite");
34 
35  // initialize command line parser
36  optutl::CommandLineParser parser("Analyze FWLite Histograms");
37 
38  // set defaults
39  parser.integerValue("maxEvents") = 1000;
40  parser.integerValue("outputEvery") = 10;
41  parser.stringValue("outputFile") = "analyzeEdmTuple.root";
42 
43  // parse arguments
44  parser.parseArguments(argc, argv);
45  int maxEvents_ = parser.integerValue("maxEvents");
46  unsigned int outputEvery_ = parser.integerValue("outputEvery");
47  std::string outputFile_ = parser.stringValue("outputFile");
48  std::vector<std::string> inputFiles_ = parser.stringVector("inputFiles");
49 
50  // book a set of histograms
52  TFileDirectory dir = fs.mkdir("analyzePatMuon");
53  TH1F* muonPt_ = dir.make<TH1F>("muonPt", "pt", 100, 0., 300.);
54  TH1F* muonEta_ = dir.make<TH1F>("muonEta", "eta", 100, -3., 3.);
55  TH1F* muonPhi_ = dir.make<TH1F>("muonPhi", "phi", 100, -5., 5.);
56 
57  // loop the events
58  int ievt = 0;
59  for (unsigned int iFile = 0; iFile < inputFiles_.size(); ++iFile) {
60  // open input file (can be located on castor)
61  TFile* inFile = TFile::Open(inputFiles_[iFile].c_str());
62  if (inFile) {
63  // ----------------------------------------------------------------------
64  // Second Part:
65  //
66  // * loop the events in the input file
67  // * receive the collections of interest via fwlite::Handle
68  // * fill the histograms
69  // * after the loop close the input file
70  // ----------------------------------------------------------------------
72  for (ev.toBegin(); !ev.atEnd(); ++ev, ++ievt) {
73  edm::EventBase const& event = ev;
74  // break loop if maximal number of events is reached
75  if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false)
76  break;
77  // simple event counter
78  if (outputEvery_ != 0 ? (ievt > 0 && ievt % outputEvery_ == 0) : false)
79  std::cout << " processing event: " << ievt << std::endl;
80 
81  // Handle to the muon pt
83  event.getByLabel(std::string("patMuonAnalyzer:pt"), muonPt);
84  // loop muon collection and fill histograms
85  for (std::vector<float>::const_iterator mu1 = muonPt->begin(); mu1 != muonPt->end(); ++mu1) {
86  muonPt_->Fill(*mu1);
87  }
88  // Handle to the muon eta
90  event.getByLabel(std::string("patMuonAnalyzer:eta"), muonEta);
91  for (std::vector<float>::const_iterator mu1 = muonEta->begin(); mu1 != muonEta->end(); ++mu1) {
92  muonEta_->Fill(*mu1);
93  }
94  // Handle to the muon phi
96  event.getByLabel(std::string("patMuonAnalyzer:phi"), muonPhi);
97  for (std::vector<float>::const_iterator mu1 = muonPhi->begin(); mu1 != muonPhi->end(); ++mu1) {
98  muonPhi_->Fill(*mu1);
99  }
100  }
101  // close input file
102  inFile->Close();
103  }
104  // break loop if maximal number of events is reached:
105  // this has to be done twice to stop the file loop as well
106  if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false)
107  break;
108  }
109  return 0;
110 }
int main(int argc, char *argv[])
static void enable()
enable automatic library loading