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/ptr_container/ptr_map.hpp>
12 #include <memory>
13 
16 
19 
22 
23 #include "TGraph.h"
24 
25 namespace {
26 
27 // Build the transformation function from the PSet format
28 std::unique_ptr<TGraph> buildTransform(const edm::ParameterSet &pset) {
29  double min = pset.getParameter<double>("min");
30  double max = pset.getParameter<double>("max");
31  const std::vector<double> &values =
32  pset.getParameter<std::vector<double> >("transform");
33  double stepSize = (max - min)/(values.size()-1);
34  std::unique_ptr<TGraph> output(new TGraph(values.size()));
35  for (size_t step = 0; step < values.size(); ++step) {
36  double x = min + step*stepSize;
37  output->SetPoint(step, x, values[step]);
38  }
39  return output;
40 }
41 
42 }
43 
45  public:
46  explicit RecoTauMVATransform(const edm::ParameterSet& pset);
47  ~RecoTauMVATransform() override {}
48 
49  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
50  double discriminate(const reco::PFTauRef&) const override;
51 
52  private:
53  // Map a decay mode to a transformation
54  typedef boost::ptr_map<reco::PFTau::hadronicDecayMode, TGraph> TransformMap;
55  TransformMap transforms_;
56  // Points to raw TaNC output
59 };
60 
61 
64  input_ = pset.getParameter<edm::InputTag>("toTransform");
65  typedef std::vector<edm::ParameterSet> VPSet;
66  const VPSet& transforms = pset.getParameter<VPSet>("transforms");
68  for(auto const& transform : transforms) {
69  unsigned int nCharged = transform.getParameter<unsigned int>("nCharged");
70  unsigned int nPiZeros = transform.getParameter<unsigned int>("nPiZeros");
71  // Get the transform
72  const edm::ParameterSet &transformImpl =
73  transform.getParameter<edm::ParameterSet>("transform");
74  // Get the acutal decay mode
76  reco::tau::translateDecayMode(nCharged, nPiZeros);
77 
78  if (!transforms_.count(decayMode)) {
79  // Add it
80  transforms_.insert(decayMode, buildTransform(transformImpl).get());
81  } else {
82  edm::LogError("DecayModeNotUnique") << "The tau decay mode with "
83  "nCharged/nPiZero = " << nCharged << "/" << nPiZeros <<
84  " dm: " << decayMode <<
85  " is associated to multiple MVA transforms, "
86  "the second instantiation is being ignored!!!";
87  }
88  }
89 }
90 
91 // Update our discriminator handle at the begninng of the event
92 void
94  evt.getByLabel(input_, disc_);
95 }
96 
98  // Check if we support this decay mode:
99  TransformMap::const_iterator transformIter =
100  transforms_.find(tau->decayMode());
101  // Unsupported DM
102  if (transformIter == transforms_.end())
104  const TGraph *transform = transformIter->second;
105  // Get the discriminator output to transform
106  double value = (*disc_)[tau];
107  double result = transform->Eval(value);
108  return result;
109 }
110 
T getParameter(std::string const &) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
PFTau::hadronicDecayMode translateDecayMode(unsigned int nCharged, unsigned int nPiZero)
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:535
boost::ptr_map< reco::PFTau::hadronicDecayMode, TGraph > TransformMap
edm::Handle< reco::PFTauDiscriminator > disc_
RecoTauMVATransform(const edm::ParameterSet &pset)
hadronicDecayMode
Definition: PFTau.h:36
step