CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoJets/JetAnalyzers/src/FlavorJetCorrectionExample.cc

Go to the documentation of this file.
00001 /* Example analizer to use L5 flavor JetCorrector services
00002    Applies different flavor corrections randomly
00003     F.Ratnikov (UMd)  Nov 16, 2007
00004 */
00005 
00006 #include <string>
00007 
00008 #include "FWCore/Framework/interface/EDAnalyzer.h"
00009 #include "FWCore/Utilities/interface/InputTag.h"
00010 
00011 class FlavorJetCorrectionExample : public edm::EDAnalyzer {
00012  public:
00013   explicit FlavorJetCorrectionExample (const edm::ParameterSet& fParameters);
00014   virtual ~FlavorJetCorrectionExample () {}
00015   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00016  private:
00017   edm::InputTag mInput;
00018   std::string mUDSCorrectorName;
00019   std::string mCCorrectorName;
00020   std::string mBCorrectorName;
00021 };
00022 
00023 
00024 
00025 
00026 #include "FWCore/Framework/interface/Event.h"
00027 #include "DataFormats/Common/interface/Handle.h"
00028 #include "FWCore/Framework/interface/EventSetup.h"
00029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00030 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00031 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
00032 
00033 using namespace std;
00034 using namespace reco;
00035 
00036 FlavorJetCorrectionExample::FlavorJetCorrectionExample (const edm::ParameterSet& fConfig) 
00037   : mInput (fConfig.getParameter <edm::InputTag> ("src")),
00038     mUDSCorrectorName (fConfig.getParameter <std::string> ("UDSQuarksCorrector")),
00039     mCCorrectorName (fConfig.getParameter <std::string> ("CQuarkCorrector")),
00040     mBCorrectorName (fConfig.getParameter <std::string> ("BQuarkCorrector"))
00041 {}
00042 
00043 void FlavorJetCorrectionExample::analyze(const edm::Event& fEvent, const edm::EventSetup& fSetup) {
00044   // get all correctors
00045   const JetCorrector* udsJetCorrector = JetCorrector::getJetCorrector (mUDSCorrectorName, fSetup);
00046   const JetCorrector* cQuarkJetCorrector = JetCorrector::getJetCorrector (mCCorrectorName, fSetup);
00047   const JetCorrector* bQuarkJetCorrector = JetCorrector::getJetCorrector (mBCorrectorName, fSetup);
00048   const JetCorrector* corrector = 0;
00049   
00050   // get input jets (supposed to be MC corrected already)
00051   edm::Handle<CaloJetCollection> jets;                    
00052   fEvent.getByLabel (mInput, jets);
00053   // loop over jets                      
00054   for (unsigned ijet = 0; ijet < jets->size(); ++ijet) {
00055     const CaloJet& jet = (*jets)[ijet];
00056     std::cout << "FlavorJetCorrectionExample::analize-> jet #" << ijet;
00057     if (ijet%3 == 0) { // assume it is light quark
00058       std::cout << ": use USD quark corrections" << std::endl;
00059       corrector = udsJetCorrector;
00060     }
00061     else if (ijet%3 == 1) { // assume it is c quark
00062       std::cout << ": use c quark corrections" << std::endl;
00063       corrector = cQuarkJetCorrector;
00064     }
00065     else { // assume it is b quark
00066       std::cout << ": use b quark corrections" << std::endl;
00067       corrector = bQuarkJetCorrector;
00068     }
00069     // get selected correction for the jet
00070     double correction = corrector->correction (jet);
00071     // dump it
00072     std::cout << "  jet pt/eta/phi: " << jet.pt() << '/' <<  jet.eta() << '/' << jet.phi() 
00073               << " -> correction factor: " << correction 
00074               << ", corrected pt: " << jet.pt()*correction
00075               << std::endl;
00076   }
00077 }
00078 
00079 #include "FWCore/Framework/interface/MakerMacros.h"
00080 DEFINE_FWK_MODULE(FlavorJetCorrectionExample);
00081 
00082