CMS 3D CMS Logo

Functions
mostPatObjects.cc File Reference
#include "FWCore/Utilities/interface/InputTag.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/PatCandidates/interface/Jet.h"
#include "DataFormats/PatCandidates/interface/MET.h"
#include "DataFormats/PatCandidates/interface/Photon.h"
#include "DataFormats/PatCandidates/interface/Electron.h"
#include "DataFormats/PatCandidates/interface/Tau.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "PhysicsTools/FWLite/interface/EventContainer.h"
#include "PhysicsTools/FWLite/interface/CommandLineParser.h"
#include "TROOT.h"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

int main ( int  argc,
char *  argv[] 
)

Definition at line 28 of file mostPatObjects.cc.

References fwlite::EventContainer::add(), fwlite::EventContainer::atEnd(), BoostedTopProducer_cfi::electronLabel, fwlite::EventContainer::getByLabel(), fwlite::EventContainer::hist(), edm::HandleBase::isValid(), BoostedTopProducer_cfi::jetLabel, GlobalHaloData_cfi::metLabel, optutl::CommandLineParser::parseArguments(), writedatasetfile::parser, photonId_cfi::photonLabel, optutl::VariableMapCont::stringValue(), fwlite::EventContainer::toBegin(), and MuonErrorMatrixAdjuster_cfi::trackLabel.

29 {
31  // ////////////////////////// //
32  // // Command Line Options // //
33  // ////////////////////////// //
35 
36  // Tell people what this analysis code does and setup default options.
37  optutl::CommandLineParser parser ("Accessing many different PAT objects");
38 
40  // Change any defaults or add any new command //
41  // line options you would like here. //
43  parser.stringValue ("outputFile") = "mostPat"; // .root added automatically
44 
45  // Parse the command line arguments
46  parser.parseArguments (argc, argv);
47 
49  // //////////////////////////// //
50  // // Create Event Container // //
51  // //////////////////////////// //
53 
54  // This object 'event' is used both to get all information from the
55  // event as well as to store histograms, etc.
56  fwlite::EventContainer eventCont (parser);
57 
59  // ////////////////////////////////// //
60  // // Begin Run // //
61  // // (e.g., book histograms, etc) // //
62  // ////////////////////////////////// //
64 
65  // Setup a style
66  gROOT->SetStyle ("Plain");
67 
68  // Book those histograms!
69  eventCont.add( new TH1F( "jetpt", "Jet Pt", 100, 0, 200) );
70  eventCont.add( new TH1F( "jetnum", "Jet Size", 100, 0, 50) );
71  eventCont.add( new TH1F( "metpt", "MET Pt", 100, 0, 200) );
72  eventCont.add( new TH1F( "photonpt", "Photon Pt", 100, 0, 200) );
73  eventCont.add( new TH1F( "trackpt", "Track Pt", 100, 0, 200) );
74  eventCont.add( new TH1F( "electronpt", "Electron Pt", 100, 0, 200) );
75  eventCont.add( new TH1F( "taupt", "Tau Pt", 100, 0, 200) );
76 
78  // //////////////// //
79  // // Event Loop // //
80  // //////////////// //
82  edm::InputTag jetLabel ("selectedLayer1Jets");
83  edm::InputTag metLabel ("selectedLayer1METs");
84  edm::InputTag photonLabel ("selectedLayer1Photons");
85  edm::InputTag trackLabel ("generalTracks");
86  edm::InputTag electronLabel ("selectedLayer1Electrons");
87  edm::InputTag tauLabel ("selectedLayer1Taus");
88 
89  for (eventCont.toBegin(); ! eventCont.atEnd(); ++eventCont)
90  {
92  // Take What We Need From Event //
94 
95  // Get jets, METs, photons, tracks, electrons, and taus
97  eventCont.getByLabel (jetLabel, h_jet);
98  assert( h_jet.isValid() );
99 
101  eventCont.getByLabel (metLabel, h_met);
102  assert( h_met.isValid() );
103 
105  eventCont.getByLabel (photonLabel, h_photon);
106  assert( h_photon.isValid() );
107 
109  eventCont.getByLabel (trackLabel, h_track);
110  assert( h_track.isValid() );
111 
113  eventCont.getByLabel (electronLabel, h_electron);
114  assert( h_electron.isValid() );
115 
117  eventCont.getByLabel (tauLabel, h_tau);
118  assert( h_tau.isValid() );
119 
120  // Fill, baby, fill!
121 
122  eventCont.hist("jetnum")->Fill( h_jet->size() );
123 
124  const vector< pat::Jet >::const_iterator kJetEnd = h_jet->end();
125  for (vector< pat::Jet >::const_iterator jetIter = h_jet->begin();
126  jetIter != kJetEnd;
127  ++jetIter)
128  {
129  eventCont.hist("jetpt")->Fill( jetIter->pt() );
130  }
131 
132  const vector< pat::MET >::const_iterator kMetEnd = h_met->end();
133  for (vector< pat::MET >::const_iterator metIter = h_met->begin();
134  metIter != kMetEnd;
135  ++metIter)
136  {
137  eventCont.hist("metpt")->Fill( metIter->pt() );
138  }
139 
140  const vector< pat::Photon >::const_iterator kPhotonEnd = h_photon->end();
141  for (vector< pat::Photon >::const_iterator photonIter = h_photon->begin();
142  photonIter != kPhotonEnd;
143  ++photonIter)
144  {
145  eventCont.hist("photonpt")->Fill( photonIter->pt() );
146  }
147 
148  const vector< reco::Track >::const_iterator kTrackEnd = h_track->end();
149  for (vector< reco::Track >::const_iterator trackIter = h_track->begin();
150  trackIter != kTrackEnd;
151  ++trackIter)
152  {
153  eventCont.hist("trackpt")->Fill( trackIter->pt() );
154  }
155 
156  const vector< pat::Electron >::const_iterator kElectronEnd = h_electron->end();
157  for (vector< pat::Electron >::const_iterator electronIter = h_electron->begin();
158  electronIter != kElectronEnd;
159  ++electronIter)
160  {
161  eventCont.hist("electronpt")->Fill( electronIter->pt() );
162  }
163 
164  const vector< pat::Tau >::const_iterator kTauEnd = h_tau->end();
165  for (vector< pat::Tau >::const_iterator tauIter = h_tau->begin();
166  tauIter != kTauEnd;
167  ++tauIter)
168  {
169  eventCont.hist ("taupt")->Fill (tauIter->pt() );
170  }
171  } // for eventCont
172 
173 
175  // ////////////////// //
176  // // Clean Up Job // //
177  // ////////////////// //
179 
180  // Histograms will be automatically written to the root file
181  // specificed by command line options.
182 
183  // All done! Bye bye.
184  return 0;
185 }
bool isValid() const
Definition: HandleBase.h:74