CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/RecoTauTag/RecoTau/plugins/RecoTauMVATransform.cc

Go to the documentation of this file.
00001 /*
00002  * ============================================================================
00003  *       Filename:  RecoTauMVATransform.cc
00004  *
00005  *    Description:  Transform TaNC output according to decay mode.
00006  *        Created:  10/22/2010 15:36:12
00007  *         Author:  Evan K. Friis (UC Davis), evan.klose.friis@cern.ch
00008  * ============================================================================
00009  */
00010 
00011 #include <boost/foreach.hpp>
00012 #include <boost/ptr_container/ptr_map.hpp>
00013 #include <memory>
00014 
00015 #include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
00016 #include "RecoTauTag/RecoTau/interface/PFTauDecayModeTools.h"
00017 
00018 #include "DataFormats/TauReco/interface/PFTau.h"
00019 #include "DataFormats/TauReco/interface/PFTauFwd.h"
00020 
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00023 
00024 #include "TGraph.h"
00025 
00026 namespace {
00027 
00028 // Build the transformation function from the PSet format
00029 std::auto_ptr<TGraph> buildTransform(const edm::ParameterSet &pset) {
00030   double min = pset.getParameter<double>("min");
00031   double max = pset.getParameter<double>("max");
00032   const std::vector<double> &values =
00033       pset.getParameter<std::vector<double> >("transform");
00034   double stepSize = (max - min)/(values.size()-1);
00035   std::auto_ptr<TGraph> output(new TGraph(values.size()));
00036   for (size_t step = 0; step < values.size(); ++step) {
00037     double x = min + step*stepSize;
00038     output->SetPoint(step, x, values[step]);
00039   }
00040   return output;
00041 }
00042 
00043 }
00044 
00045 class RecoTauMVATransform : public PFTauDiscriminationProducerBase {
00046   public:
00047     explicit RecoTauMVATransform(const edm::ParameterSet& pset);
00048     ~RecoTauMVATransform() {}
00049 
00050     void beginEvent(const edm::Event&, const edm::EventSetup&);
00051     double discriminate(const reco::PFTauRef&);
00052 
00053   private:
00054     // Map a decay mode to a transformation
00055     typedef boost::ptr_map<reco::PFTau::hadronicDecayMode, TGraph> TransformMap;
00056     TransformMap transforms_;
00057     // Points to raw TaNC output
00058     edm::InputTag input_;
00059     edm::Handle<reco::PFTauDiscriminator> disc_;
00060 };
00061 
00062 
00063 RecoTauMVATransform::RecoTauMVATransform(const edm::ParameterSet& pset)
00064   :PFTauDiscriminationProducerBase(pset) {
00065   input_ = pset.getParameter<edm::InputTag>("toTransform");
00066   typedef std::vector<edm::ParameterSet> VPSet;
00067   const VPSet& transforms = pset.getParameter<VPSet>("transforms");
00068   prediscriminantFailValue_ = -2.0;
00069   BOOST_FOREACH(const edm::ParameterSet &transform, transforms) {
00070     unsigned int nCharged = transform.getParameter<unsigned int>("nCharged");
00071     unsigned int nPiZeros = transform.getParameter<unsigned int>("nPiZeros");
00072     // Get the transform
00073     const edm::ParameterSet &transformImpl =
00074         transform.getParameter<edm::ParameterSet>("transform");
00075     // Get the acutal decay mode
00076     reco::PFTau::hadronicDecayMode decayMode =
00077         reco::tau::translateDecayMode(nCharged, nPiZeros);
00078 
00079     if (!transforms_.count(decayMode)) {
00080       // Add it
00081       transforms_.insert(decayMode, buildTransform(transformImpl));
00082     } else {
00083       edm::LogError("DecayModeNotUnique") << "The tau decay mode with "
00084         "nCharged/nPiZero = " << nCharged << "/" << nPiZeros <<
00085         " dm: " << decayMode <<
00086         " is associated to multiple MVA transforms, "
00087         "the second instantiation is being ignored!!!";
00088     }
00089   }
00090 }
00091 
00092 // Update our discriminator handle at the begninng of the event
00093 void
00094 RecoTauMVATransform::beginEvent(const edm::Event& evt, const edm::EventSetup&) {
00095   evt.getByLabel(input_, disc_);
00096 }
00097 
00098 double RecoTauMVATransform::discriminate(const reco::PFTauRef& tau) {
00099   // Check if we support this decay mode:
00100   TransformMap::const_iterator transformIter =
00101       transforms_.find(tau->decayMode());
00102   // Unsupported DM
00103   if (transformIter == transforms_.end())
00104     return prediscriminantFailValue_;
00105   const TGraph *transform = transformIter->second;
00106   // Get the discriminator output to transform
00107   double value = (*disc_)[tau];
00108   double result = transform->Eval(value);
00109   return result;
00110 }
00111 
00112 #include "FWCore/Framework/interface/MakerMacros.h"
00113 DEFINE_FWK_MODULE(RecoTauMVATransform);