CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/PhysicsTools/FWLite/examples/playingWithJets.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 #include "DataFormats/Math/interface/deltaPhi.h"
00008 
00009 #include "PhysicsTools/FWLite/interface/EventContainer.h"
00010 #include "PhysicsTools/FWLite/interface/CommandLineParser.h" 
00011 
00012 // Root includes
00013 #include "TROOT.h"
00014 
00015 using namespace std;
00016 
00017 
00019 // ///////////////////// //
00020 // // Main Subroutine // //
00021 // ///////////////////// //
00023 
00024 int main (int argc, char* argv[]) 
00025 {
00027    // ////////////////////////// //
00028    // // Command Line Options // //
00029    // ////////////////////////// //
00031 
00032    // Tell people what this analysis code does and setup default options.
00033    optutl::CommandLineParser parser ("Playing with jets");
00034 
00036    // Change any defaults or add any new command //
00037    //      line options you would like here.     //
00039    // change default output filename
00040    parser.stringValue ("outputFile") = "jetInfo.root";
00041    
00042    // Parse the command line arguments
00043    parser.parseArguments (argc, argv);
00044 
00046    // //////////////////////////// //
00047    // // Create Event Container // //
00048    // //////////////////////////// //
00050 
00051    // This object 'event' is used both to get all information from the
00052    // event as well as to store histograms, etc.
00053    fwlite::EventContainer eventCont (parser);
00054 
00056    // ////////////////////////////////// //
00057    // //         Begin Run            // //
00058    // // (e.g., book histograms, etc) // //
00059    // ////////////////////////////////// //
00061 
00062    // Setup a style
00063    gROOT->SetStyle ("Plain");
00064 
00065    // Book those histograms!
00066    eventCont.add( new TH1F( "jetpt",        "Jet p_{T} using standard absolute p_{T} calibration", 100, 0, 60) );
00067    eventCont.add( new TH1F( "jeteta",       "Jet eta using standard absolute p_{T} calibration",   100, 0, 10) );
00068    eventCont.add( new TH1F( "reljetpt",     "Jet p_{T} using relative inter eta calibration",      100, 0, 60) );
00069    eventCont.add( new TH1F( "reljeteta",    "Jet eta using relative inter eta calibration",        100, 0, 10) );
00070    eventCont.add( new TH1F( "phijet1jet2",  "Phi between Jet 1 and Jet 2",                        100, 0, 3.5) );
00071    eventCont.add( new TH1F( "invarMass",    "Invariant Mass of the 4-vector sum of Two Jets",     100, 0, 200) );
00072 
00074    // //////////////// //
00075    // // Event Loop // //
00076    // //////////////// //
00078 
00079    // create labels
00080    edm::InputTag jetLabel ("selectedLayer1Jets");
00081 
00082    for (eventCont.toBegin(); ! eventCont.atEnd(); ++eventCont) 
00083    {
00085       // Take What We Need From Event //
00087 
00088       // Get jets
00089       edm::Handle< vector< pat::Jet > > jetHandle;
00090       eventCont.getByLabel (jetLabel, jetHandle);
00091       assert ( jetHandle.isValid() );
00092 
00093       const vector< pat::Jet >::const_iterator kJetEnd = jetHandle->end();
00094       for (vector< pat::Jet >::const_iterator jetIter = jetHandle->begin(); 
00095            jetIter != kJetEnd; 
00096            ++jetIter) 
00097       {    
00098          // A few correctedJet options:
00099          // ABS - absolute pt calibration (automatic)
00100          // REL - relative inter eta calibration  
00101          // EMF - calibration as a function of the jet EMF
00102          eventCont.hist("reljetpt") ->Fill( jetIter->correctedJet("REL").pt() );
00103          eventCont.hist("reljeteta")->Fill( jetIter->correctedJet("REL").eta() );
00104          eventCont.hist("jetpt")    ->Fill( jetIter->correctedJet("ABS").pt() );
00105          // Automatically ABS
00106          eventCont.hist("jeteta")   ->Fill( jetIter->eta() );
00107       } // for jetIter
00108 
00109       // Do we have at least two jets?
00110       if (jetHandle->size() < 2)
00111       {
00112          // Nothing to do here
00113          continue;
00114       }
00115 
00116       // Store invariant mass and delta phi between two leading jets.
00117       eventCont.hist("invarMass")->Fill( (jetHandle->at(0).p4() + jetHandle->at(1).p4()).M() );
00118       eventCont.hist("phijet1jet2")->Fill( deltaPhi( jetHandle->at(0).phi(), jetHandle->at(1).phi() ) );
00119    } // for eventCont
00120 
00121       
00123    // ////////////////// //
00124    // // Clean Up Job // //
00125    // ////////////////// //
00127 
00128    // Histograms will be automatically written to the root file
00129    // specificed by command line options.
00130 
00131    // All done!  Bye bye.
00132    return 0;
00133 }