CMS 3D CMS Logo

RecoTauDecayModeCutMultiplexer.cc
Go to the documentation of this file.
3 
6 
8  public:
10 
12  double discriminate(const reco::PFTauRef&) const override;
13  void beginEvent(const edm::Event& event, const edm::EventSetup& eventSetup) override;
14 
15  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
16  private:
17  typedef std::pair<unsigned int, unsigned int> IntPair;
18  typedef std::map<IntPair, double> DecayModeCutMap;
19 
20  DecayModeCutMap decayModeCuts_;
24 };
25 
28  toMultiplex_ = pset.getParameter<edm::InputTag>("toMultiplex");
29  toMultiplex_token = consumes<reco::PFTauDiscriminator>(toMultiplex_);
30  typedef std::vector<edm::ParameterSet> VPSet;
31  const VPSet& decayModes = pset.getParameter<VPSet>("decayModes");
32  // Setup our cut map
33  for(auto const& dm : decayModes) {
34  // Get the mass window for each decay mode
35  decayModeCuts_.insert(std::make_pair(
36  // The decay mode as a key
37  std::make_pair(
38  dm.getParameter<uint32_t>("nCharged"),
39  dm.getParameter<uint32_t>("nPiZeros")),
40  // The selection
41  dm.getParameter<double>("cut")
42  ));
43  }
44 }
45 
46 void
48  const edm::Event& evt, const edm::EventSetup& es) {
50 }
51 
52 double
54  double disc_result = (*handle_)[tau];
55  DecayModeCutMap::const_iterator cutIter =
56  decayModeCuts_.find(std::make_pair(tau->signalChargedHadrCands().size(),
57  tau->signalPiZeroCandidates().size()));
58 
59  // Return null if it doesn't exist
60  if (cutIter == decayModeCuts_.end()) {
62  }
63  // See if the discriminator passes our cuts
64  return disc_result > cutIter->second;
65 }
66 
67 void
69  // recoTauDecayModeCutMultiplexer
71 
72  desc.add<edm::InputTag>("toMultiplex");
73 
74  edm::ParameterSetDescription vpsd_decayModes;
75 
76  vpsd_decayModes.add<uint32_t>("nCharged");
77  vpsd_decayModes.add<uint32_t>("nPiZeros");
78  vpsd_decayModes.add<double>("cut");
79 
80  // name description no defaults items
81  desc.addVPSet("decayModes", vpsd_decayModes);
82 
83  fillProducerDescriptions(desc); // inherited from the base
84 
85  descriptions.add("recoTauDecayModeCutMultiplexer", desc);
86 }
87 
88 
T getParameter(std::string const &) const
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
edm::Handle< reco::PFTauDiscriminator > handle_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::pair< unsigned int, unsigned int > IntPair
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillProducerDescriptions(edm::ParameterSetDescription &desc)
edm::EDGetTokenT< reco::PFTauDiscriminator > toMultiplex_token
void beginEvent(const edm::Event &event, const edm::EventSetup &eventSetup) override
RecoTauDecayModeCutMultiplexer(const edm::ParameterSet &pset)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double discriminate(const reco::PFTauRef &) const override
Definition: event.py:1