CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTauTag/RecoTau/plugins/RecoTauDiscriminantInvariantWidth.cc

Go to the documentation of this file.
00001 /*
00002  * RecoTauDiscriminantInvariantWidth
00003  *
00004  * Compute a (hopefully) p_T independent quantity related to the
00005  * opening angle.
00006  *
00007  * Author: Evan K. Friis (UC Davis)
00008  *
00009  */
00010 
00011 #include <boost/foreach.hpp>
00012 #include "RecoTauTag/RecoTau/interface/RecoTauDiscriminantPlugins.h"
00013 #include "RecoTauTag/RecoTau/interface/PFTauDecayModeTools.h"
00014 #include "DataFormats/Math/interface/deltaR.h"
00015 #include <TFormula.h>
00016 #include "RecoTauTag/RecoTau/interface/RecoTauDiscriminantFunctions.h"
00017 #include "CommonTools/Utils/interface/StringObjectFunction.h"
00018 
00019 namespace reco { namespace tau {
00020 
00021 class RecoTauDiscriminantInvariantWidth : public RecoTauDiscriminantPlugin {
00022   public:
00023     explicit RecoTauDiscriminantInvariantWidth(
00024         const edm::ParameterSet& pset);
00025     std::vector<double> operator()(const reco::PFTauRef& tau) const;
00026   private:
00027     typedef StringObjectFunction<PFTau> TauFunc;
00028     typedef boost::shared_ptr<TauFunc> TauFuncPtr;
00029     typedef std::pair<TauFuncPtr, TauFuncPtr> MeanAndWidthFuncs;
00030 
00031     std::map<reco::PFTau::hadronicDecayMode, MeanAndWidthFuncs> transforms_;
00032     MeanAndWidthFuncs defaultTransform_;
00033 };
00034 
00035 RecoTauDiscriminantInvariantWidth::RecoTauDiscriminantInvariantWidth(
00036     const edm::ParameterSet& pset):RecoTauDiscriminantPlugin(pset) {
00037   typedef std::vector<edm::ParameterSet> VPSet;
00038   // Add each of the transformations
00039   BOOST_FOREACH(const edm::ParameterSet& dm,
00040       pset.getParameter<VPSet>("decayModes")) {
00041     uint32_t nCharged = dm.getParameter<uint32_t>("nCharged");
00042     uint32_t nPiZeros = dm.getParameter<uint32_t>("nPiZeros");
00043     MeanAndWidthFuncs functions;
00044     functions.first.reset(new TauFunc(dm.getParameter<std::string>("mean")));
00045     functions.second.reset(new TauFunc(dm.getParameter<std::string>("rms")));
00046     transforms_[translateDecayMode(nCharged, nPiZeros)] = functions;
00047   }
00048   defaultTransform_.first.reset(
00049       new TauFunc(pset.getParameter<std::string>("defaultMean")));
00050   defaultTransform_.second.reset(
00051       new TauFunc(pset.getParameter<std::string>("defaultRMS")));
00052 }
00053 
00054 std::vector<double> RecoTauDiscriminantInvariantWidth::operator()(
00055     const reco::PFTauRef& tau) const {
00056   double weightedDeltaR = disc::OpeningDeltaR(*tau);
00057 
00058   std::map<reco::PFTau::hadronicDecayMode, MeanAndWidthFuncs>::const_iterator
00059     transform = transforms_.find(tau->decayMode());
00060 
00061   const TauFunc* meanFunc = defaultTransform_.first.get();
00062   const TauFunc* rmsFunc = defaultTransform_.second.get();
00063 
00064   if (transform != transforms_.end()) {
00065     meanFunc = transform->second.first.get();
00066     rmsFunc = transform->second.second.get();
00067   }
00068 
00069   double mean = (*meanFunc)(*tau);
00070   double rms = (*rmsFunc)(*tau);
00071 
00072   double result = (rms > 0) ? (weightedDeltaR - mean)/rms : -1.;
00073 
00074   return std::vector<double>(1, result);
00075 }
00076 
00077 }} // end namespace reco::tau
00078 
00079 #include "FWCore/Framework/interface/MakerMacros.h"
00080 DEFINE_EDM_PLUGIN(RecoTauDiscriminantPluginFactory,
00081     reco::tau::RecoTauDiscriminantInvariantWidth,
00082     "RecoTauDiscriminantInvariantWidth");