CMS 3D CMS Logo

CorJetsExample Class Reference

#include <RecoJets/JetAnalyzers/interface/CorJetsExample.h>

Inheritance diagram for CorJetsExample:

edm::EDAnalyzer

List of all members.

Public Member Functions

 CorJetsExample (const edm::ParameterSet &)

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &)
void beginJob (const edm::EventSetup &)
void endJob ()

Private Attributes

std::string CaloJetAlgorithm
std::string CorJetAlgorithm
std::string GenJetAlgorithm
TH1F h_ptCal
TH1F h_ptCor
TH1F h_ptCorOnFly
TH1F h_ptGen
std::string JetCorrectionService
TFile * m_file


Detailed Description

Definition at line 15 of file CorJetsExample.h.


Constructor & Destructor Documentation

CorJetsExample::CorJetsExample ( const edm::ParameterSet cfg  ) 

Definition at line 24 of file CorJetsExample.cc.

00024                                                          :
00025   CaloJetAlgorithm( cfg.getParameter<string>( "CaloJetAlgorithm" ) ), 
00026   CorJetAlgorithm( cfg.getParameter<string>( "CorJetAlgorithm" ) ), 
00027   JetCorrectionService( cfg.getParameter<string>( "JetCorrectionService" ) ), 
00028   GenJetAlgorithm( cfg.getParameter<string>( "GenJetAlgorithm" ) )
00029   {
00030 }


Member Function Documentation

void CorJetsExample::analyze ( const edm::Event evt,
const edm::EventSetup es 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 42 of file CorJetsExample.cc.

References CaloJetAlgorithm, CorJetAlgorithm, JetCorrector::correction(), GenJetAlgorithm, edm::Event::getByLabel(), JetCorrector::getJetCorrector(), h_ptCal, h_ptCor, h_ptCorOnFly, h_ptGen, JetCorrectionService, and scale.

00042                                                                      {
00043 
00044   //Get the CaloJet collection
00045   Handle<CaloJetCollection> caloJets;
00046   evt.getByLabel( CaloJetAlgorithm, caloJets );
00047 
00048   //Loop over the two leading CaloJets and fill a histogram
00049   int jetInd = 0;
00050   for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end() && jetInd<2; ++ cal ) {
00051     h_ptCal.Fill( cal->pt() );   
00052     jetInd++;
00053   }
00054 
00055   //Get the GenJet collection
00056   Handle<GenJetCollection> genJets;
00057   evt.getByLabel( GenJetAlgorithm, genJets );
00058 
00059   //Loop over the two leading GenJets and fill a histogram
00060   jetInd = 0;
00061   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end() && jetInd<2; ++ gen ) {
00062     h_ptGen.Fill( gen->pt() );   
00063     jetInd++;
00064   }
00065 
00066   //Get the Corrected Jet collection
00067   Handle<CaloJetCollection> corJets;
00068   evt.getByLabel( CorJetAlgorithm, corJets );
00069 
00070   //Loop over all the CorJets, save the two highest Pt, and fill the histogram.
00071   //Corrected jets have the original Pt order as the uncorrected CaloJets for releases before 1.5.0,
00072   // and hence need to be re-ordered to find the two leading jets.
00073   jetInd = 0;
00074   double highestPt=0.0;
00075   double nextPt=0.0;
00076   for( CaloJetCollection::const_iterator cor = corJets->begin(); cor != corJets->end(); ++ cor ) {
00077     double corPt=cor->pt();
00078     //std::cout << "cor=" << cor->pt() << std::endl;
00079     if(corPt>highestPt){
00080       nextPt=highestPt;
00081       highestPt=corPt;
00082     }
00083     else if(corPt>nextPt)nextPt=corPt;
00084   }
00085   h_ptCor.Fill( highestPt );   
00086   h_ptCor.Fill( nextPt );
00087 
00088   const JetCorrector* corrector = 
00089                  JetCorrector::getJetCorrector (JetCorrectionService, es);
00090   //Loop over the CaloJets, correct them on the fly, save the two highest Pt and fill the histogram.
00091   highestPt=0.0;
00092   nextPt=0.0;
00093   for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) {
00094     double scale = corrector->correction (*cal);
00095     double corPt=scale*cal->pt();
00096     //std::cout << corPt << ", scale=" << scale << std::endl;  
00097     if(corPt>highestPt){
00098       nextPt=highestPt;
00099       highestPt=corPt;
00100     }
00101     else if(corPt>nextPt)nextPt=corPt;
00102   }
00103   h_ptCorOnFly.Fill( highestPt );   
00104   h_ptCorOnFly.Fill( nextPt );
00105   //std::cout <<  "corOnFly=" << highestPt  << std::endl;
00106   //std::cout <<  "corOnFly=" << nextPt  << std::endl;
00107 
00108 
00109 }

void CorJetsExample::beginJob ( const edm::EventSetup  )  [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 32 of file CorJetsExample.cc.

References h_ptCal, h_ptCor, h_ptCorOnFly, h_ptGen, and m_file.

00032                                                   {
00033 
00034   // Open the histogram file and book some associated histograms
00035   m_file=new TFile("histo.root","RECREATE"); 
00036   h_ptCal =  TH1F( "h_ptCal",  "p_{T} of leading CaloJets", 50, 0, 1000 );
00037   h_ptGen =  TH1F( "h_ptGen",  "p_{T} of leading GenJets", 50, 0, 1000 );
00038   h_ptCor =  TH1F( "h_ptCor",  "p_{T} of leading CorJets", 50, 0, 1000 );
00039   h_ptCorOnFly =  TH1F( "h_ptCorOnFly",  "p_{T} of leading Jets Corrected on the Fly", 50, 0, 1000 );
00040 }

void CorJetsExample::endJob ( void   )  [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 111 of file CorJetsExample.cc.

References m_file.

00111                             {
00112 
00113   //Write out the histogram file.
00114   m_file->Write(); 
00115 
00116 }


Member Data Documentation

std::string CorJetsExample::CaloJetAlgorithm [private]

Definition at line 23 of file CorJetsExample.h.

Referenced by analyze().

std::string CorJetsExample::CorJetAlgorithm [private]

Definition at line 23 of file CorJetsExample.h.

Referenced by analyze().

std::string CorJetsExample::GenJetAlgorithm [private]

Definition at line 23 of file CorJetsExample.h.

Referenced by analyze().

TH1F CorJetsExample::h_ptCal [private]

Definition at line 24 of file CorJetsExample.h.

Referenced by analyze(), and beginJob().

TH1F CorJetsExample::h_ptCor [private]

Definition at line 24 of file CorJetsExample.h.

Referenced by analyze(), and beginJob().

TH1F CorJetsExample::h_ptCorOnFly [private]

Definition at line 24 of file CorJetsExample.h.

Referenced by analyze(), and beginJob().

TH1F CorJetsExample::h_ptGen [private]

Definition at line 24 of file CorJetsExample.h.

Referenced by analyze(), and beginJob().

std::string CorJetsExample::JetCorrectionService [private]

Definition at line 23 of file CorJetsExample.h.

Referenced by analyze().

TFile* CorJetsExample::m_file [private]

Definition at line 25 of file CorJetsExample.h.

Referenced by beginJob(), and endJob().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:16:49 2009 for CMSSW by  doxygen 1.5.4