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
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(const edm::ParameterSet&);
00023 ~PFRecoTauDiscriminationByHPSSelection();
00024 double discriminate(const reco::PFTauRef&);
00025
00026 private:
00027 typedef StringObjectFunction<reco::PFTau> TauFunc;
00028
00029 struct DecayModeCuts
00030 {
00031 DecayModeCuts()
00032 : maxMass_(0)
00033 {}
00034 ~DecayModeCuts() {}
00035 double minMass_;
00036 TauFunc* maxMass_;
00037 double minPi0Mass_;
00038 double maxPi0Mass_;
00039 double assumeStripMass_;
00040 };
00041
00042 typedef std::pair<unsigned int, unsigned int> IntPair;
00043 typedef std::pair<double, double> DoublePair;
00044 typedef std::map<IntPair, DecayModeCuts> DecayModeCutMap;
00045
00046 TauFunc signalConeFun_;
00047 DecayModeCutMap decayModeCuts_;
00048 double matchingCone_;
00049 double minPt_;
00050 };
00051
00052 PFRecoTauDiscriminationByHPSSelection::PFRecoTauDiscriminationByHPSSelection(const edm::ParameterSet& pset)
00053 : PFTauDiscriminationProducerBase(pset),
00054 signalConeFun_(pset.getParameter<std::string>("coneSizeFormula"))
00055 {
00056
00057 matchingCone_ = pset.getParameter<double>("matchingCone");
00058 minPt_ = pset.getParameter<double>("minTauPt");
00059
00060 typedef std::vector<edm::ParameterSet> VPSet;
00061 const VPSet& decayModes = pset.getParameter<VPSet>("decayModes");
00062 BOOST_FOREACH(const edm::ParameterSet &dm, decayModes) {
00063
00064 DecayModeCuts cuts;
00065 cuts.minMass_ = dm.getParameter<double>("minMass");
00066 cuts.maxMass_ = new TauFunc(dm.getParameter<std::string>("maxMass"));
00067 if (dm.exists("minPi0Mass")) {
00068 cuts.minPi0Mass_ = dm.getParameter<double>("minPi0Mass");
00069 cuts.maxPi0Mass_ = dm.getParameter<double>("maxPi0Mass");
00070 } else {
00071 cuts.minPi0Mass_ = -1.e3;
00072 cuts.maxPi0Mass_ = 1.e9;
00073 }
00074 if (dm.exists("assumeStripMass")) {
00075 cuts.assumeStripMass_ = dm.getParameter<double>("assumeStripMass");
00076 } else {
00077 cuts.assumeStripMass_ = -1.0;
00078 }
00079 decayModeCuts_.insert(std::make_pair(
00080
00081 std::make_pair(
00082 dm.getParameter<uint32_t>("nCharged"),
00083 dm.getParameter<uint32_t>("nPiZeros")),
00084 cuts
00085 ));
00086 }
00087 }
00088
00089 PFRecoTauDiscriminationByHPSSelection::~PFRecoTauDiscriminationByHPSSelection()
00090 {
00091 for ( DecayModeCutMap::iterator it = decayModeCuts_.begin();
00092 it != decayModeCuts_.end(); ++it ) {
00093 delete it->second.maxMass_;
00094 }
00095 }
00096
00097 double
00098 PFRecoTauDiscriminationByHPSSelection::discriminate(const reco::PFTauRef& tau)
00099 {
00100
00101
00102
00103 if (tau->pt() < minPt_)
00104 return 0.0;
00105
00106
00107 DecayModeCutMap::const_iterator massWindowIter =
00108 decayModeCuts_.find(std::make_pair(tau->signalPFChargedHadrCands().size(),
00109 tau->signalPiZeroCandidates().size()));
00110
00111 if (massWindowIter == decayModeCuts_.end()) {
00112 return 0.0;
00113 }
00114 const DecayModeCuts& massWindow = massWindowIter->second;
00115
00116 math::XYZTLorentzVector tauP4 = tau->p4();
00117
00118
00119 reco::Candidate::LorentzVector stripsP4;
00120 BOOST_FOREACH(const reco::RecoTauPiZero& cand,
00121 tau->signalPiZeroCandidates()){
00122 math::XYZTLorentzVector candP4 = cand.p4();
00123 stripsP4 += candP4;
00124 }
00125
00126
00127 if (massWindow.assumeStripMass_ >= 0) {
00128 BOOST_FOREACH(const reco::RecoTauPiZero& cand,
00129 tau->signalPiZeroCandidates()){
00130 math::XYZTLorentzVector uncorrected = cand.p4();
00131 math::XYZTLorentzVector corrected =
00132 applyMassConstraint(uncorrected, massWindow.assumeStripMass_);
00133 math::XYZTLorentzVector correction = corrected - uncorrected;
00134 tauP4 += correction;
00135 stripsP4 += correction;
00136 }
00137 }
00138
00139
00140
00141 double maxMass_value = (*massWindow.maxMass_)(*tau);
00142 if (tauP4.M() > maxMass_value || tauP4.M() < massWindow.minMass_) {
00143 return 0.0;
00144 }
00145
00146
00147 if (stripsP4.M() > massWindow.maxPi0Mass_ ||
00148 stripsP4.M() < massWindow.minPi0Mass_) {
00149 return 0.0;
00150 }
00151
00152
00153
00154 if (deltaR(tauP4, tau->jetRef()->p4()) > matchingCone_) {
00155 return 0.0;
00156 }
00157
00158
00159 double cone_size = signalConeFun_(*tau);
00160
00161 BOOST_FOREACH(const reco::PFCandidateRef& cand,
00162 tau->signalPFChargedHadrCands()) {
00163
00164 if (deltaR(cand->p4(), tauP4) > cone_size)
00165 return 0.0;
00166 }
00167
00168 BOOST_FOREACH(const reco::RecoTauPiZero& cand,
00169 tau->signalPiZeroCandidates()) {
00170
00171 if (deltaR(cand.p4(), tauP4) > cone_size)
00172 return 0.0;
00173 }
00174
00175
00176 return 1.0;
00177 }
00178
00179 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByHPSSelection);