CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/JetMETCorrections/Type1MET/plugins/TauMET.cc

Go to the documentation of this file.
00001 // Package:    TauMET
00002 // Class:      TauMET
00003 // 
00004 // Original Authors:  Alfredo Gurrola, Chi Nhan Nguyen
00005 
00006 #include "JetMETCorrections/Type1MET/plugins/TauMET.h"
00007 
00008 #include "DataFormats/METReco/interface/MET.h"
00009 #include "DataFormats/METReco/interface/METCollection.h"
00010 #include "DataFormats/METReco/interface/CaloMET.h"
00011 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00012 #include "DataFormats/METReco/interface/CorrMETData.h"
00013 
00014 #include "DataFormats/TauReco/interface/PFTau.h"
00015 #include "DataFormats/JetReco/interface/CaloJet.h"
00016 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00017 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
00018 
00019 #include "FWCore/Framework/interface/MakerMacros.h"
00020 
00021 using namespace reco;
00022 using namespace std;
00023 namespace cms 
00024 {
00025 
00026   TauMET::TauMET(const edm::ParameterSet& iConfig) : _algo() {
00027 
00028     _InputTausLabel    = iConfig.getParameter<std::string>("InputTausLabel");
00029     _tauType    = iConfig.getParameter<std::string>("tauType");
00030     _InputCaloJetsLabel    = iConfig.getParameter<std::string>("InputCaloJetsLabel");
00031     _jetPTthreshold      = iConfig.getParameter<double>("jetPTthreshold");
00032     _jetEMfracLimit      = iConfig.getParameter<double>("jetEMfracLimit");
00033     _correctorLabel      = iConfig.getParameter<std::string>("correctorLabel");
00034     _InputMETLabel    = iConfig.getParameter<std::string>("InputMETLabel");
00035     _metType    = iConfig.getParameter<std::string>("metType");
00036     _JetMatchDeltaR      = iConfig.getParameter<double>("JetMatchDeltaR");
00037     _TauMinEt      = iConfig.getParameter<double>("TauMinEt");
00038     _TauEtaMax      = iConfig.getParameter<double>("TauEtaMax");
00039     _UseSeedTrack      = iConfig.getParameter<bool>("UseSeedTrack");
00040     _seedTrackPt      = iConfig.getParameter<double>("seedTrackPt");
00041     _UseTrackIsolation      = iConfig.getParameter<bool>("UseTrackIsolation");
00042     _trackIsolationMinPt      = iConfig.getParameter<double>("trackIsolationMinPt");
00043     _UseECALIsolation      = iConfig.getParameter<bool>("UseECALIsolation");
00044     _gammaIsolationMinPt      = iConfig.getParameter<double>("gammaIsolationMinPt");
00045     _UseProngStructure      = iConfig.getParameter<bool>("UseProngStructure");
00046 
00047     if( _metType == "recoMET" ) {
00048       produces< METCollection   >();
00049     } else {
00050       produces< CaloMETCollection   >();
00051     }
00052 
00053   }
00054   
00055 
00056   TauMET::~TauMET() {
00057     // do anything here that needs to be done at desctruction time
00058     // (e.g. close files, deallocate resources etc.)
00059   }
00060   
00061   // ------------ method called to produce the data  ------------
00062   void TauMET::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00063 
00064     using namespace edm;
00065 
00066     Handle<CaloJetCollection> calojetHandle;
00067     iEvent.getByLabel(_InputCaloJetsLabel, calojetHandle);
00068     const JetCorrector* correctedjets = JetCorrector::getJetCorrector (_correctorLabel, iSetup);
00069 
00070     Handle<PFTauCollection> tauHandle;
00071     iEvent.getByLabel(_InputTausLabel,tauHandle);
00072 
00073     if( _metType == "recoCaloMET" ) {
00074       Handle<CaloMETCollection> metHandle;
00075       iEvent.getByLabel(_InputMETLabel,metHandle);
00076       std::auto_ptr< CaloMETCollection > output( new CaloMETCollection() );
00077       _algo.run(iEvent,iSetup,tauHandle,calojetHandle,_jetPTthreshold,_jetEMfracLimit,*correctedjets,*(metHandle.product()),
00078                 _JetMatchDeltaR,_TauMinEt,
00079                 _TauEtaMax,_UseSeedTrack,_seedTrackPt,_UseTrackIsolation,_trackIsolationMinPt,_UseECALIsolation,
00080                 _gammaIsolationMinPt,_UseProngStructure,&*output);
00081       iEvent.put( output );
00082     } else if( _metType == "recoMET" ) {
00083       Handle<METCollection> metHandle;
00084       iEvent.getByLabel(_InputMETLabel,metHandle);
00085       std::auto_ptr< METCollection > output( new METCollection() );
00086       _algo.run(iEvent,iSetup,tauHandle,calojetHandle,_jetPTthreshold,_jetEMfracLimit,*correctedjets,*(metHandle.product()),
00087                 _JetMatchDeltaR,_TauMinEt,
00088                 _TauEtaMax,_UseSeedTrack,_seedTrackPt,_UseTrackIsolation,_trackIsolationMinPt,_UseECALIsolation,
00089                 _gammaIsolationMinPt,_UseProngStructure,&*output);
00090       iEvent.put( output );
00091     } else {
00092       std::cerr << "Incorrect Met Type!!! " << std::endl;
00093       std::cerr << "Please re-run and set the metType to 'recoCaloMET' or 'recoMET' " << std::endl;
00094       return;
00095     }
00096 
00097   }
00098   
00099   // ------------ method called once each job just before starting event loop  ------------
00100   void TauMET::beginJob() { }
00101   
00102   // ------------ method called once each job just after ending the event loop  ------------
00103   void TauMET::endJob() { }
00104 
00105   DEFINE_FWK_MODULE(TauMET); //define this as a plug-in
00106 
00107 } //end namespace cms
00108 
00109