CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTauTag/RecoTau/plugins/PFRecoTauDiscriminationByHPSSelection.cc

Go to the documentation of this file.
00001 #include <boost/foreach.hpp>
00002 
00003 #include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
00004 #include "FWCore/Utilities/interface/InputTag.h"
00005 
00006 #include "CommonTools/Utils/interface/StringObjectFunction.h"
00007 #include "DataFormats/Math/interface/deltaR.h"
00008 
00009 namespace {
00010 // Apply a hypothesis on the mass of the strips.
00011 math::XYZTLorentzVector applyMassConstraint(
00012     const math::XYZTLorentzVector& vec,double mass) {
00013   double factor = sqrt(vec.energy()*vec.energy()-mass*mass)/vec.P();
00014   return math::XYZTLorentzVector(
00015       vec.px()*factor,vec.py()*factor,vec.pz()*factor,vec.energy());
00016 }
00017 }
00018 
00019 class PFRecoTauDiscriminationByHPSSelection
00020   : public PFTauDiscriminationProducerBase {
00021   public:
00022     explicit PFRecoTauDiscriminationByHPSSelection(
00023         const edm::ParameterSet& pset);
00024 
00025     ~PFRecoTauDiscriminationByHPSSelection() {}
00026     double discriminate(const reco::PFTauRef&);
00027 
00028   private:
00029     struct DecayModeCuts {
00030       double minMass_;
00031       double maxMass_;
00032       double minPi0Mass_;
00033       double maxPi0Mass_;
00034       double assumeStripMass_;
00035     };
00036 
00037     typedef StringObjectFunction<reco::PFTau> TauFunc;
00038     typedef std::pair<unsigned int, unsigned int> IntPair;
00039     typedef std::pair<double, double> DoublePair;
00040     typedef std::map<IntPair, DecayModeCuts> DecayModeCutMap;
00041 
00042     TauFunc signalConeFun_;
00043     DecayModeCutMap decayModeCuts_;
00044     double matchingCone_;
00045     double minPt_;
00046 };
00047 
00048 PFRecoTauDiscriminationByHPSSelection::PFRecoTauDiscriminationByHPSSelection(
00049     const edm::ParameterSet& pset):PFTauDiscriminationProducerBase(pset),
00050     signalConeFun_(pset.getParameter<std::string>("coneSizeFormula")) {
00051   // Get the matchign cut
00052   matchingCone_ = pset.getParameter<double>("matchingCone");
00053   minPt_ = pset.getParameter<double>("minTauPt");
00054   // Get the mass cuts for each decay mode
00055   typedef std::vector<edm::ParameterSet> VPSet;
00056   const VPSet& decayModes = pset.getParameter<VPSet>("decayModes");
00057   BOOST_FOREACH(const edm::ParameterSet &dm, decayModes) {
00058     // The mass window(s)
00059     DecayModeCuts cuts;
00060     cuts.minMass_ = dm.getParameter<double>("minMass");
00061     cuts.maxMass_ = dm.getParameter<double>("maxMass");
00062     if (dm.exists("minPi0Mass")) {
00063       cuts.minPi0Mass_ = dm.getParameter<double>("minPi0Mass");
00064       cuts.maxPi0Mass_ = dm.getParameter<double>("maxPi0Mass");
00065     } else {
00066       cuts.minPi0Mass_ = -1.0;
00067       cuts.maxPi0Mass_ = 1e9;
00068     }
00069     if (dm.exists("assumeStripMass")) {
00070       cuts.assumeStripMass_ = dm.getParameter<double>("assumeStripMass");
00071     } else {
00072       cuts.assumeStripMass_ = -1.0;
00073     }
00074     decayModeCuts_.insert(std::make_pair(
00075             // The decay mode as a key
00076             std::make_pair(
00077                 dm.getParameter<uint32_t>("nCharged"),
00078                 dm.getParameter<uint32_t>("nPiZeros")),
00079             cuts
00080           ));
00081   }
00082 }
00083 
00084 double
00085 PFRecoTauDiscriminationByHPSSelection::discriminate(const reco::PFTauRef& tau) {
00086   // Check if we pass the min pt
00087   if (tau->pt() < minPt_)
00088     return 0.0;
00089 
00090   // See if we select this decay mode
00091   DecayModeCutMap::const_iterator massWindowIter =
00092       decayModeCuts_.find(std::make_pair(tau->signalPFChargedHadrCands().size(),
00093                                          tau->signalPiZeroCandidates().size()));
00094   // Check if decay mode is supported
00095   if (massWindowIter == decayModeCuts_.end()) {
00096     return 0.0;
00097   }
00098   const DecayModeCuts &massWindow = massWindowIter->second;
00099 
00100   math::XYZTLorentzVector tauP4 = tau->p4();
00101   // Find the total pizero p4
00102   reco::Candidate::LorentzVector stripsP4;
00103   BOOST_FOREACH(const reco::RecoTauPiZero& cand, 
00104       tau->signalPiZeroCandidates()){
00105     math::XYZTLorentzVector candP4 = cand.p4();
00106     stripsP4 += candP4;
00107   }
00108 
00109   // Apply strip mass assumption corrections
00110   if (massWindow.assumeStripMass_ >= 0) {
00111     BOOST_FOREACH(const reco::RecoTauPiZero& cand, 
00112         tau->signalPiZeroCandidates()){
00113       math::XYZTLorentzVector uncorrected = cand.p4();
00114       math::XYZTLorentzVector corrected = 
00115         applyMassConstraint(uncorrected, massWindow.assumeStripMass_);
00116       math::XYZTLorentzVector correction = corrected - uncorrected;
00117       tauP4 += correction;
00118       stripsP4 += correction;
00119     }
00120   }
00121 
00122   // Check if tau fails mass cut
00123   if (tauP4.M() > massWindow.maxMass_ || tauP4.M() < massWindow.minMass_) {
00124     return 0.0;
00125   }
00126 
00127   // Check if it fails the pi 0 IM cut
00128   if (stripsP4.M() > massWindow.maxPi0Mass_ ||
00129       stripsP4.M() < massWindow.minPi0Mass_) {
00130     return 0.0;
00131   }
00132 
00133   // Check if tau passes matching cone cut
00134   if (deltaR(tauP4, tau->jetRef()->p4()) > matchingCone_) {
00135     return 0.0;
00136   }
00137 
00138   // Check if tau passes cone cut
00139   double cone_size = signalConeFun_(*tau);
00140   // Check if any charged objects fail the signal cone cut
00141   BOOST_FOREACH(const reco::PFCandidateRef& cand,
00142                 tau->signalPFChargedHadrCands()) {
00143     if (deltaR(cand->p4(), tauP4) > cone_size)
00144       return 0.0;
00145   }
00146   // Now check the pizeros
00147   BOOST_FOREACH(const reco::RecoTauPiZero& cand,
00148                 tau->signalPiZeroCandidates()) {
00149     if (deltaR(cand.p4(), tauP4) > cone_size)
00150       return 0.0;
00151   }
00152 
00153   // Otherwise, we pass!
00154   return 1.0;
00155 }
00156 
00157 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByHPSSelection);