CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/RecoTauTag/TauTagTools/plugins/RecoTauObjectEmbedderPlugin.cc

Go to the documentation of this file.
00001 /*
00002  * ===========================================================================
00003  *
00004  *       Filename:  RecoTauObjectEmbedder
00005  *
00006  *    Description:  Embed truth information into (currently unused)
00007  *                  tau variables.  This is a hack to allow some PAT-like
00008  *                  operations on taus without introducing new dependencies.
00009  *
00010  *         Author:  Evan K. Friis (UC Davis)
00011  *
00012  * ===========================================================================
00013  */
00014 
00015 #include <boost/foreach.hpp>
00016 
00017 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h"
00018 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00019 #include "RecoTauTag/RecoTau/interface/PFTauDecayModeTools.h"
00020 #include "DataFormats/Common/interface/Handle.h"
00021 #include "DataFormats/Common/interface/Association.h"
00022 
00023 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00024 #include "DataFormats/JetReco/interface/GenJet.h"
00025 
00026 namespace reco { namespace tau {
00027 
00028 namespace helpers {
00029 unsigned int nCharged(const GenJet& jet) {
00030   unsigned int output = 0;
00031   BOOST_FOREACH(const CandidatePtr &cand, jet.getJetConstituents()) {
00032     if (cand->charge())
00033       ++output;
00034   }
00035   return output;
00036 }
00037 
00038 unsigned int nGammas(const GenJet& jet) {
00039   unsigned int output = 0;
00040   BOOST_FOREACH(const CandidatePtr &cand, jet.getJetConstituents()) {
00041     if (cand->pdgId()==22)
00042       ++output;
00043   }
00044   return output;
00045 }
00046 
00047 unsigned int nCharged(const PFTau& tau) {
00048   return tau.signalPFChargedHadrCands().size();
00049 }
00050 unsigned int nGammas(const PFTau& tau) {
00051   return tau.signalPiZeroCandidates().size();
00052 }
00053 }
00054 
00055 template<typename T>
00056 class RecoTauObjectEmbedder : public RecoTauModifierPlugin {
00057   public:
00058     explicit RecoTauObjectEmbedder(const edm::ParameterSet &pset)
00059         :RecoTauModifierPlugin(pset),
00060         jetMatchSrc_(pset.getParameter<edm::InputTag>("jetTruthMatch")) {}
00061     virtual ~RecoTauObjectEmbedder() {}
00062     virtual void operator()(PFTau&) const;
00063     virtual void beginEvent();
00064   private:
00065     edm::InputTag jetMatchSrc_;
00066     edm::Handle<edm::Association<T> > jetMatch_;
00067 };
00068 
00069 // Update our handle to the matching
00070 template<typename T>
00071 void RecoTauObjectEmbedder<T>::beginEvent() {
00072   evt()->getByLabel(jetMatchSrc_, jetMatch_);
00073 }
00074 
00075 template<typename T>
00076 void RecoTauObjectEmbedder<T>::operator()(PFTau& tau) const {
00077   // Get the matched object that is matched to the same jet as the current tau,
00078   // if it exists
00079   edm::Ref<T> matchedObject = (*jetMatch_)[tau.jetRef()];
00080   if (matchedObject.isNonnull()) {
00081     // Store our matched object information
00082     tau.setalternatLorentzVect(matchedObject->p4());
00083     // Store our generator decay mode
00084     tau.setbremsRecoveryEOverPLead(
00085         reco::tau::translateDecayMode(
00086             helpers::nCharged(*matchedObject),
00087             helpers::nGammas(*matchedObject)/2)
00088         );
00089   } else {
00090     tau.setbremsRecoveryEOverPLead(-10);
00091   }
00092 }
00093 }}  // end namespace reco::tau
00094 
00095 #include "FWCore/Framework/interface/MakerMacros.h"
00096 DEFINE_EDM_PLUGIN(RecoTauModifierPluginFactory,
00097     reco::tau::RecoTauObjectEmbedder<reco::GenJetCollection>,
00098     "RecoTauTruthEmbedder");
00099 
00100 DEFINE_EDM_PLUGIN(RecoTauModifierPluginFactory,
00101     reco::tau::RecoTauObjectEmbedder<reco::PFTauCollection>,
00102     "RecoTauPFTauEmbedder");