CMS 3D CMS Logo

PFRecoTauDiscriminationByInvMass.cc
Go to the documentation of this file.
3 
4 /* class PFRecoTauDiscriminationByInvMass
5  * created : August 30 2010,
6  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
7  * contributors : Evan Friis (UC Davis)
8  * based on H+ tau ID by Lauri Wendland
9  */
10 
12  public:
15  // If select is not set, just return the invariant mass
16  cut_ = pset.exists("select");
17  if (cut_) {
19  ("select");
20  // Get default cuts
21  min_default_ = select.getParameter<double>("min");
22  max_default_ = select.getParameter<double>("max");
23  // Get decay mode specific cuts
24  std::vector<std::string> decayModeCutNames =
26  for(auto const& dmName : decayModeCutNames) {
27  const edm::ParameterSet &dmPSet =
28  select.getParameter<edm::ParameterSet>(dmName);
29  unsigned int nCharged = dmPSet.getParameter<unsigned int>("charged");
30  unsigned int nPiZero = dmPSet.getParameter<unsigned int>("pizeros");
31  double minCut = dmPSet.getParameter<double>("min");
32  double maxCut = dmPSet.getParameter<double>("max");
33  // Add our dm-specific cut to the map
34  decayModeCuts_[std::make_pair(nCharged, nPiZero)] =
35  std::make_pair(minCut, maxCut);
36  }
37  }
38  }
40  double discriminate(const reco::PFTauRef&) const override;
41 
42  private:
43  typedef std::pair<unsigned int, unsigned int> IntPair;
44  typedef std::pair<double, double> DoublePair;
45  typedef std::map<IntPair, DoublePair> DecayModeCutMap;
46  DecayModeCutMap decayModeCuts_;
47  double min_default_;
48  double max_default_;
49  bool cut_;
50 };
51 
52 double
54  double mass = tau->mass();
55  if (cut_) {
56  unsigned int charged = tau->signalPFChargedHadrCands().size();
57  unsigned int pizeros = tau->signalPiZeroCandidates().size();
58  DecayModeCutMap::const_iterator specificCut = decayModeCuts_.find(
59  std::make_pair(charged, pizeros));
60  // Cut does not exist for this decay mode
61  if (specificCut == decayModeCuts_.end() )
62  return (mass > min_default_ && mass < max_default_);
63  else
64  return (mass > specificCut->second.first &&
65  mass < specificCut->second.second);
66  }
67  // If we dont' cut, just return the mass
68  return mass;
69 }
70 
72 
T getParameter(std::string const &) const
std::map< IntPair, DoublePair > DecayModeCutMap
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< std::string > getParameterNamesForType(bool trackiness=true) const
Definition: ParameterSet.h:193
U second(std::pair< T, U > const &p)
std::pair< unsigned int, unsigned int > IntPair
PFRecoTauDiscriminationByInvMass(const edm::ParameterSet &pset)
double discriminate(const reco::PFTauRef &) const override