CMS 3D CMS Logo

PFRecoTauDiscriminationByInvMass.cc
Go to the documentation of this file.
3 
6 
7 /* class PFRecoTauDiscriminationByInvMass
8  * created : August 30 2010,
9  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
10  * contributors : Evan Friis (UC Davis)
11  * based on H+ tau ID by Lauri Wendland
12  */
13 
15  public:
18  // If select is not set, just return the invariant mass
19  const edm::ParameterSet &select = pset.getParameter<edm::ParameterSet>("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  }
39  double discriminate(const reco::PFTauRef&) const override;
40 
41  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
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 };
50 
51 double
53  double mass = tau->mass();
54  unsigned int charged = tau->signalChargedHadrCands().size();
55  unsigned int pizeros = tau->signalPiZeroCandidates().size();
56  DecayModeCutMap::const_iterator specificCut = decayModeCuts_.find(
57  std::make_pair(charged, pizeros));
58  // Cut does not exist for this decay mode
59  if (specificCut == decayModeCuts_.end() )
60  return (mass > min_default_ && mass < max_default_);
61  else
62  return (mass > specificCut->second.first &&
63  mass < specificCut->second.second);
64  // If we dont' cut, just return the mass
65  return mass;
66 }
67 
68 void
70  // pfRecoTauDiscriminationByInvMass
72  {
74  psd0.add<std::string>("BooleanOperator", "and");
75  {
77  psd1.add<double>("cut");
78  psd1.add<edm::InputTag>("Producer");
79  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
80  }
81  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
82  }
83  {
85  psd0.add<double>("max");
86  psd0.add<double>("min");
87  desc.add<edm::ParameterSetDescription>("select", psd0);
88  }
89  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
90  descriptions.add("pfRecoTauDiscriminationByInvMass", desc);
91 }
92 
T getParameter(std::string const &) const
std::map< IntPair, DoublePair > DecayModeCutMap
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
std::vector< std::string > getParameterNamesForType(bool trackiness=true) const
Definition: ParameterSet.h:169
U second(std::pair< T, U > const &p)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::pair< unsigned int, unsigned int > IntPair
PFRecoTauDiscriminationByInvMass(const edm::ParameterSet &pset)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double discriminate(const reco::PFTauRef &) const override