00001 #include <boost/foreach.hpp> 00002 00003 #include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h" 00004 #include "FWCore/Utilities/interface/InputTag.h" 00005 00006 /* class PFRecoTauDiscriminationByInvMass 00007 * created : August 30 2010, 00008 * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki) 00009 * contributors : Evan Friis (UC Davis) 00010 * based on H+ tau ID by Lauri Wendland 00011 */ 00012 00013 class PFRecoTauDiscriminationByInvMass: public PFTauDiscriminationProducerBase { 00014 public: 00015 explicit PFRecoTauDiscriminationByInvMass(const edm::ParameterSet& pset) 00016 :PFTauDiscriminationProducerBase(pset) { 00017 // If select is not set, just return the invariant mass 00018 cut_ = pset.exists("select"); 00019 if (cut_) { 00020 const edm::ParameterSet &select = pset.getParameter<edm::ParameterSet> 00021 ("select"); 00022 // Get default cuts 00023 min_default_ = select.getParameter<double>("min"); 00024 max_default_ = select.getParameter<double>("max"); 00025 // Get decay mode specific cuts 00026 std::vector<std::string> decayModeCutNames = 00027 select.getParameterNamesForType<edm::ParameterSet>(); 00028 BOOST_FOREACH(const std::string& dmName, decayModeCutNames) { 00029 const edm::ParameterSet &dmPSet = 00030 select.getParameter<edm::ParameterSet>(dmName); 00031 unsigned int nCharged = dmPSet.getParameter<unsigned int>("charged"); 00032 unsigned int nPiZero = dmPSet.getParameter<unsigned int>("pizeros"); 00033 double minCut = dmPSet.getParameter<double>("min"); 00034 double maxCut = dmPSet.getParameter<double>("max"); 00035 // Add our dm-specific cut to the map 00036 decayModeCuts_[std::make_pair(nCharged, nPiZero)] = 00037 std::make_pair(minCut, maxCut); 00038 } 00039 } 00040 } 00041 ~PFRecoTauDiscriminationByInvMass(){} 00042 double discriminate(const reco::PFTauRef&); 00043 00044 private: 00045 typedef std::pair<unsigned int, unsigned int> IntPair; 00046 typedef std::pair<double, double> DoublePair; 00047 typedef std::map<IntPair, DoublePair> DecayModeCutMap; 00048 DecayModeCutMap decayModeCuts_; 00049 double min_default_; 00050 double max_default_; 00051 bool cut_; 00052 }; 00053 00054 double 00055 PFRecoTauDiscriminationByInvMass::discriminate(const reco::PFTauRef& tau) { 00056 double mass = tau->mass(); 00057 if (cut_) { 00058 unsigned int charged = tau->signalPFChargedHadrCands().size(); 00059 unsigned int pizeros = tau->signalPiZeroCandidates().size(); 00060 DecayModeCutMap::const_iterator specificCut = decayModeCuts_.find( 00061 std::make_pair(charged, pizeros)); 00062 // Cut does not exist for this decay mode 00063 if (specificCut == decayModeCuts_.end() ) 00064 return (mass > min_default_ && mass < max_default_); 00065 else 00066 return (mass > specificCut->second.first && 00067 mass < specificCut->second.second); 00068 } 00069 // If we dont' cut, just return the mass 00070 return mass; 00071 } 00072 00073 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByInvMass); 00074