CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
jetPt.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 
3 // CMS includes
7 
10 
11 // Root includes
12 #include "TROOT.h"
13 
14 using namespace std;
15 
16 
18 // ///////////////////// //
19 // // Main Subroutine // //
20 // ///////////////////// //
22 
23 int main (int argc, char* argv[])
24 {
26  // ////////////////////////// //
27  // // Command Line Options // //
28  // ////////////////////////// //
30 
31  // Tell people what this analysis code does and setup default options.
32  optutl::CommandLineParser parser ("Plots Jet Pt");
33 
35  // Change any defaults or add any new command //
36  // line options you would like here. //
38  parser.stringValue ("outputFile") = "jetPt"; // .root added automatically
39 
40  // Parse the command line arguments
41  parser.parseArguments (argc, argv);
42 
44  // //////////////////////////// //
45  // // Create Event Container // //
46  // //////////////////////////// //
48 
49  // This object 'event' is used both to get all information from the
50  // event as well as to store histograms, etc.
51  fwlite::EventContainer eventCont (parser);
52 
54  // ////////////////////////////////// //
55  // // Begin Run // //
56  // // (e.g., book histograms, etc) // //
57  // ////////////////////////////////// //
59 
60  // Setup a style
61  gROOT->SetStyle ("Plain");
62 
63  // Book those histograms!
64  eventCont.add( new TH1F( "jetPt", "jetPt", 1000, 0, 1000) );
65 
67  // //////////////// //
68  // // Event Loop // //
69  // //////////////// //
71 
72  // create labels
73  edm::InputTag jetLabel ("selectedLayer1Jets");
74 
75  for (eventCont.toBegin(); ! eventCont.atEnd(); ++eventCont)
76  {
78  // Take What We Need From Event //
81  eventCont.getByLabel (jetLabel, jetHandle);
82  assert ( jetHandle.isValid() );
83 
84  // Loop over the jets
85  const vector< pat::Jet >::const_iterator kJetEnd = jetHandle->end();
86  for (vector< pat::Jet >::const_iterator jetIter = jetHandle->begin();
87  kJetEnd != jetIter;
88  ++jetIter)
89  {
90  eventCont.hist("jetPt")->Fill (jetIter->pt());
91  } // for jetIter
92  } // for eventCont
93 
94 
96  // ////////////////// //
97  // // Clean Up Job // //
98  // ////////////////// //
100 
101  // Histograms will be automatically written to the root file
102  // specificed by command line options.
103 
104  // All done! Bye bye.
105  return 0;
106 }
std::string & stringValue(std::string key)
void parseArguments(int argc, char **argv, bool allowArgs=false)
assert(m_qm.get())
bool getByLabel(const std::type_info &iInfo, const char *iModuleLabel, const char *iProductInstanceLabel, const char *iProcessLabel, void *oData) const
void add(TH1 *histPtr, const std::string &directory="")
bool isValid() const
Definition: HandleBase.h:75
tuple argc
Definition: dir2webdir.py:38
const EventContainer & toBegin()
TH1 * hist(const std::string &name)