CMS 3D CMS Logo

PFRecoTauDecayModeIndexProducer.cc
Go to the documentation of this file.
1 /*
2  * class PFRecoTauDecayModeIndexProducer
3  * created : April 16, 2009
4  * revised : ,
5  * Authors : Evan K. Friis, (UC Davis), Simone Gennai (SNS)
6  *
7  * Associates the decay mode index (see enum in DataFormats/TauReco/interface/PFTauDecayMode.h)
8  * reconstruced in the PFTauDecayMode production to its underlying PFTau, in PFTau discriminator format
9  *
10  */
11 
14 
16 
19 
20 namespace {
21 
22 using namespace reco;
23 
24 class PFRecoTauDecayModeIndexProducer final : public PFTauDiscriminationProducerBase {
25  public:
26  explicit PFRecoTauDecayModeIndexProducer(const edm::ParameterSet& iConfig):PFTauDiscriminationProducerBase(iConfig) {
27  PFTauDecayModeProducer_ = iConfig.getParameter<edm::InputTag>("PFTauDecayModeProducer");
28  }
29  ~PFRecoTauDecayModeIndexProducer() override{}
30  double discriminate(const PFTauRef& thePFTauRef) const override ;
31  void beginEvent(const edm::Event& evt, const edm::EventSetup& evtSetup) override;
32  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
33  private:
34  edm::InputTag PFTauDecayModeProducer_;
35  edm::Handle<PFTauDecayModeAssociation> decayModes_; // holds the PFTauDecayModes for the current event
36 };
37 
38 
39 void PFRecoTauDecayModeIndexProducer::beginEvent(const edm::Event& event, const edm::EventSetup& evtSetup)
40 {
41  // Get the PFTau Decay Modes
42  event.getByLabel(PFTauDecayModeProducer_, decayModes_);
43 }
44 
45 double PFRecoTauDecayModeIndexProducer::discriminate(const PFTauRef& thePFTauRef) const
46 {
47  int theDecayModeIndex = -1;
48 
49  const PFTauDecayMode& theDecayMode = (*decayModes_)[thePFTauRef];
50 
51  // retrieve decay mode
52  theDecayModeIndex = static_cast<int>(theDecayMode.getDecayMode());
53 
54  if (theDecayModeIndex < 0) theDecayModeIndex = -1;
55 
56  return theDecayModeIndex;
57 }
58 
59 }
60 
61 void
63  // pfTauDecayModeIndexProducer
65  desc.add<edm::InputTag>("PFTauDecayModeProducer", edm::InputTag("pfRecoTauDecayModeProducer"));
66  {
68  psd0.add<std::string>("BooleanOperator", "and");
69  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
70  }
71  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
72  descriptions.add("pfTauDecayModeIndexProducer", desc);
73 }
74 
75 DEFINE_FWK_MODULE(PFRecoTauDecayModeIndexProducer);
T getParameter(std::string const &) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
hadronicTauDecayModes getDecayMode() const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
Definition: event.py:1