00001
00002
00003
00004
00005
00006 #include "RecoJets/JetAnalyzers/interface/JetPlotsExample.h"
00007 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00008 #include "DataFormats/JetReco/interface/CaloJet.h"
00009 #include "DataFormats/JetReco/interface/GenJet.h"
00010 #include "FWCore/Framework/interface/Event.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012 #include <TROOT.h>
00013 #include <TSystem.h>
00014 #include <TFile.h>
00015 #include <TCanvas.h>
00016 #include <cmath>
00017 using namespace edm;
00018 using namespace reco;
00019 using namespace std;
00020
00021
00022
00023 JetPlotsExample::JetPlotsExample( const ParameterSet & cfg ) :
00024 CaloJetAlgorithm( cfg.getParameter<string>( "CaloJetAlgorithm" ) ),
00025 GenJetAlgorithm( cfg.getParameter<string>( "GenJetAlgorithm" ) )
00026 {
00027 }
00028
00029 void JetPlotsExample::beginJob( const EventSetup & ) {
00030
00031
00032 m_file=new TFile("histo.root","RECREATE");
00033 h_ptCal = TH1F( "ptCal", "p_{T} of leading CaloJets", 50, 0, 1000 );
00034 h_etaCal = TH1F( "etaCal", "#eta of leading CaloJets", 50, -3, 3 );
00035 h_phiCal = TH1F( "phiCal", "#phi of leading CaloJets", 50, -M_PI, M_PI );
00036 h_ptGen = TH1F( "ptGen", "p_{T} of leading GenJets", 50, 0, 1000 );
00037 h_etaGen = TH1F( "etaGen", "#eta of leading GenJets", 50, -3, 3 );
00038 h_phiGen = TH1F( "phiGen", "#phi of leading GenJets", 50, -M_PI, M_PI );
00039 }
00040
00041 void JetPlotsExample::analyze( const Event& evt, const EventSetup& es ) {
00042
00043
00044 Handle<CaloJetCollection> caloJets;
00045 evt.getByLabel( CaloJetAlgorithm, caloJets );
00046
00047
00048 int jetInd = 0;
00049 for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end() && jetInd<2; ++ cal ) {
00050
00051 h_ptCal.Fill( cal->pt() );
00052 h_etaCal.Fill( cal->eta() );
00053 h_phiCal.Fill( cal->phi() );
00054 jetInd++;
00055 }
00056
00057
00058 Handle<GenJetCollection> genJets;
00059 evt.getByLabel( GenJetAlgorithm, genJets );
00060
00061
00062 jetInd = 0;
00063 for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end() && jetInd<2; ++ gen ) {
00064
00065 h_ptGen.Fill( gen->pt() );
00066 h_etaGen.Fill( gen->eta() );
00067 h_phiGen.Fill( gen->phi() );
00068 jetInd++;
00069 }
00070
00071 }
00072
00073 void JetPlotsExample::endJob() {
00074 if (m_file !=0)
00075 {
00076
00077 m_file->cd();
00078 h_ptCal.Write();
00079 h_etaCal.Write();
00080 h_phiCal.Write();
00081 h_ptGen.Write();
00082 h_etaGen.Write();
00083 h_phiGen.Write();
00084 delete m_file;
00085 m_file=0;
00086 }
00087
00088 }
00089 #include "FWCore/Framework/interface/MakerMacros.h"
00090 DEFINE_FWK_MODULE(JetPlotsExample);