CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauDiscriminantInvariantWidth.cc
Go to the documentation of this file.
1 /*
2  * RecoTauDiscriminantInvariantWidth
3  *
4  * Compute a (hopefully) p_T independent quantity related to the
5  * opening angle.
6  *
7  * Author: Evan K. Friis (UC Davis)
8  *
9  */
10 
11 #include <boost/foreach.hpp>
15 #include <TFormula.h>
18 
19 namespace reco { namespace tau {
20 
22  public:
24  const edm::ParameterSet& pset);
25  std::vector<double> operator()(const reco::PFTauRef& tau) const;
26  private:
28  typedef boost::shared_ptr<TauFunc> TauFuncPtr;
29  typedef std::pair<TauFuncPtr, TauFuncPtr> MeanAndWidthFuncs;
30 
31  std::map<reco::PFTau::hadronicDecayMode, MeanAndWidthFuncs> transforms_;
33 };
34 
37  typedef std::vector<edm::ParameterSet> VPSet;
38  // Add each of the transformations
39  BOOST_FOREACH(const edm::ParameterSet& dm,
40  pset.getParameter<VPSet>("decayModes")) {
41  uint32_t nCharged = dm.getParameter<uint32_t>("nCharged");
42  uint32_t nPiZeros = dm.getParameter<uint32_t>("nPiZeros");
43  MeanAndWidthFuncs functions;
44  functions.first.reset(new TauFunc(dm.getParameter<std::string>("mean")));
45  functions.second.reset(new TauFunc(dm.getParameter<std::string>("rms")));
46  transforms_[translateDecayMode(nCharged, nPiZeros)] = functions;
47  }
48  defaultTransform_.first.reset(
49  new TauFunc(pset.getParameter<std::string>("defaultMean")));
50  defaultTransform_.second.reset(
51  new TauFunc(pset.getParameter<std::string>("defaultRMS")));
52 }
53 
55  const reco::PFTauRef& tau) const {
56  double weightedDeltaR = disc::OpeningDeltaR(*tau);
57 
58  std::map<reco::PFTau::hadronicDecayMode, MeanAndWidthFuncs>::const_iterator
59  transform = transforms_.find(tau->decayMode());
60 
61  const TauFunc* meanFunc = defaultTransform_.first.get();
62  const TauFunc* rmsFunc = defaultTransform_.second.get();
63 
64  if (transform != transforms_.end()) {
65  meanFunc = transform->second.first.get();
66  rmsFunc = transform->second.second.get();
67  }
68 
69  double mean = (*meanFunc)(*tau);
70  double rms = (*rmsFunc)(*tau);
71 
72  double result = (rms > 0) ? (weightedDeltaR - mean)/rms : -1.;
73 
74  return std::vector<double>(1, result);
75 }
76 
77 }} // end namespace reco::tau
78 
82  "RecoTauDiscriminantInvariantWidth");
T getParameter(std::string const &) const
PFTau::hadronicDecayMode translateDecayMode(unsigned int nCharged, unsigned int nPiZero)
std::vector< double > operator()(const reco::PFTauRef &tau) const
std::map< reco::PFTau::hadronicDecayMode, MeanAndWidthFuncs > transforms_
tuple result
Definition: query.py:137
unsigned int nCharged(const GenJet &jet)
#define DEFINE_EDM_PLUGIN(factory, type, name)