CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/PhysicsTools/FWLite/examples/jetPt.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 
00003 // CMS includes
00004 #include "FWCore/Utilities/interface/InputTag.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "DataFormats/PatCandidates/interface/Jet.h"
00007 
00008 #include "PhysicsTools/FWLite/interface/EventContainer.h"
00009 #include "PhysicsTools/FWLite/interface/CommandLineParser.h" 
00010 
00011 // Root includes
00012 #include "TROOT.h"
00013 
00014 using namespace std;
00015 
00016 
00018 // ///////////////////// //
00019 // // Main Subroutine // //
00020 // ///////////////////// //
00022 
00023 int main (int argc, char* argv[]) 
00024 {
00026    // ////////////////////////// //
00027    // // Command Line Options // //
00028    // ////////////////////////// //
00030 
00031    // Tell people what this analysis code does and setup default options.
00032    optutl::CommandLineParser parser ("Plots Jet Pt");
00033 
00035    // Change any defaults or add any new command //
00036    //      line options you would like here.     //
00038    parser.stringValue ("outputFile") = "jetPt"; // .root added automatically
00039 
00040    // Parse the command line arguments
00041    parser.parseArguments (argc, argv);
00042 
00044    // //////////////////////////// //
00045    // // Create Event Container // //
00046    // //////////////////////////// //
00048 
00049    // This object 'event' is used both to get all information from the
00050    // event as well as to store histograms, etc.
00051    fwlite::EventContainer eventCont (parser);
00052 
00054    // ////////////////////////////////// //
00055    // //         Begin Run            // //
00056    // // (e.g., book histograms, etc) // //
00057    // ////////////////////////////////// //
00059 
00060    // Setup a style
00061    gROOT->SetStyle ("Plain");
00062 
00063    // Book those histograms!
00064    eventCont.add( new TH1F( "jetPt", "jetPt", 1000, 0, 1000) );
00065 
00067    // //////////////// //
00068    // // Event Loop // //
00069    // //////////////// //
00071 
00072    // create labels
00073    edm::InputTag jetLabel ("selectedLayer1Jets");
00074 
00075    for (eventCont.toBegin(); ! eventCont.atEnd(); ++eventCont) 
00076    {
00078       // Take What We Need From Event //
00080       edm::Handle< vector< pat::Jet > > jetHandle;
00081       eventCont.getByLabel (jetLabel, jetHandle);
00082       assert ( jetHandle.isValid() );
00083                                                 
00084       // Loop over the jets
00085       const vector< pat::Jet >::const_iterator kJetEnd = jetHandle->end();
00086       for (vector< pat::Jet >::const_iterator jetIter = jetHandle->begin();
00087            kJetEnd != jetIter; 
00088            ++jetIter) 
00089       {         
00090          eventCont.hist("jetPt")->Fill (jetIter->pt());
00091       } // for jetIter
00092    } // for eventCont
00093 
00094       
00096    // ////////////////// //
00097    // // Clean Up Job // //
00098    // ////////////////// //
00100 
00101    // Histograms will be automatically written to the root file
00102    // specificed by command line options.
00103 
00104    // All done!  Bye bye.
00105    return 0;
00106 }