CMS 3D CMS Logo

RecoTauMVATransform.cc
Go to the documentation of this file.
1 /*
2  * ============================================================================
3  * Filename: RecoTauMVATransform.cc
4  *
5  * Description: Transform TaNC output according to decay mode.
6  * Created: 10/22/2010 15:36:12
7  * Author: Evan K. Friis (UC Davis), evan.klose.friis@cern.ch
8  * ============================================================================
9  */
10 
11 #include <boost/foreach.hpp>
12 #include <boost/ptr_container/ptr_map.hpp>
13 #include <memory>
14 
17 
20 
23 
24 #include "TGraph.h"
25 
26 namespace {
27 
28 // Build the transformation function from the PSet format
29 std::auto_ptr<TGraph> buildTransform(const edm::ParameterSet &pset) {
30  double min = pset.getParameter<double>("min");
31  double max = pset.getParameter<double>("max");
32  const std::vector<double> &values =
33  pset.getParameter<std::vector<double> >("transform");
34  double stepSize = (max - min)/(values.size()-1);
35  std::auto_ptr<TGraph> output(new TGraph(values.size()));
36  for (size_t step = 0; step < values.size(); ++step) {
37  double x = min + step*stepSize;
38  output->SetPoint(step, x, values[step]);
39  }
40  return output;
41 }
42 
43 }
44 
46  public:
47  explicit RecoTauMVATransform(const edm::ParameterSet& pset);
49 
50  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
51  double discriminate(const reco::PFTauRef&) const override;
52 
53  private:
54  // Map a decay mode to a transformation
55  typedef boost::ptr_map<reco::PFTau::hadronicDecayMode, TGraph> TransformMap;
56  TransformMap transforms_;
57  // Points to raw TaNC output
60 };
61 
62 
65  input_ = pset.getParameter<edm::InputTag>("toTransform");
66  typedef std::vector<edm::ParameterSet> VPSet;
67  const VPSet& transforms = pset.getParameter<VPSet>("transforms");
69  BOOST_FOREACH(const edm::ParameterSet &transform, transforms) {
70  unsigned int nCharged = transform.getParameter<unsigned int>("nCharged");
71  unsigned int nPiZeros = transform.getParameter<unsigned int>("nPiZeros");
72  // Get the transform
73  const edm::ParameterSet &transformImpl =
74  transform.getParameter<edm::ParameterSet>("transform");
75  // Get the acutal decay mode
77  reco::tau::translateDecayMode(nCharged, nPiZeros);
78 
79  if (!transforms_.count(decayMode)) {
80  // Add it
81  transforms_.insert(decayMode, buildTransform(transformImpl));
82  } else {
83  edm::LogError("DecayModeNotUnique") << "The tau decay mode with "
84  "nCharged/nPiZero = " << nCharged << "/" << nPiZeros <<
85  " dm: " << decayMode <<
86  " is associated to multiple MVA transforms, "
87  "the second instantiation is being ignored!!!";
88  }
89  }
90 }
91 
92 // Update our discriminator handle at the begninng of the event
93 void
95  evt.getByLabel(input_, disc_);
96 }
97 
99  // Check if we support this decay mode:
100  TransformMap::const_iterator transformIter =
101  transforms_.find(tau->decayMode());
102  // Unsupported DM
103  if (transformIter == transforms_.end())
105  const TGraph *transform = transformIter->second;
106  // Get the discriminator output to transform
107  double value = (*disc_)[tau];
108  double result = transform->Eval(value);
109  return result;
110 }
111 
T getParameter(std::string const &) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
PFTau::hadronicDecayMode translateDecayMode(unsigned int nCharged, unsigned int nPiZero)
T x() const
Cartesian x coordinate.
double discriminate(const reco::PFTauRef &) const override
void beginEvent(const edm::Event &, const edm::EventSetup &) override
Definition: value.py:1
T min(T a, T b)
Definition: MathUtil.h:58
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:416
boost::ptr_map< reco::PFTau::hadronicDecayMode, TGraph > TransformMap
edm::Handle< reco::PFTauDiscriminator > disc_
RecoTauMVATransform(const edm::ParameterSet &pset)
hadronicDecayMode
Definition: PFTau.h:36
step