CMS 3D CMS Logo

PFTauDecayModeCutMultiplexer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PFTauDecayModeCutMultiplexer
4 // Class: PFTauDecayModeCutMultiplexer
5 //
6 /*
7 
8  Description: Applies a different cut to a PFTauDiscriminator, depending on the
9  the reconstructed DecayMode (stored by RecoTauTag/RecoTau/PFTauDecayModeIndexProducer)
10  in PFTauDiscriminator form.
11 
12  Produces a PFTauDiscriminator output with a binary (0 or 1) output.
13 
14  Cuts are specified in the decay mode PSets, which map the cuts
15  to collections of decay mode indices. These decay mode PSets are defined
16  in the same manner as TaNC MVA computers (which also map to specific decay modes)
17 
18 
19 */
20 //
21 // Original Author: Evan K. Friis, UC Davis (friis@physics.ucdavis.edu)
22 // Created: Thurs, April 16, 2009
23 //
24 //
25 
27 
30 
31 using namespace reco;
32 
34  public:
37 
38  struct ComputerAndCut {
40  double userCut;
41  };
42 
43  typedef std::vector<ComputerAndCut> CutList;
44  typedef std::map<int, CutList::iterator> DecayModeToCutMap;
45 
46  double discriminate(const PFTauRef& thePFTau) const override;
47  void beginEvent(const edm::Event& event, const edm::EventSetup& eventSetup) override;
48 
49  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
50  private:
51  // PFTau discriminator continaing the decaymode index of the tau collection
54 
55  // Discriminant to multiplex cut on
58 
59  DecayModeToCutMap computerMap_; //Maps decay mode to MVA implementation
60  CutList computers_;
61 
62  edm::Handle<PFTauDiscriminator> pfTauDecayModeIndices; // holds the decay mode indices for the taus in the current event
63  edm::Handle<PFTauDiscriminator> targetDiscriminant; // holds the discirminant values of the discriminant we will multiplex for the current event
64 };
65 
67 {
68  pfTauDecayModeIndexSrc_ = iConfig.getParameter<edm::InputTag>("PFTauDecayModeSrc");
69  discriminantToMultiplex_ = iConfig.getParameter<edm::InputTag>("PFTauDiscriminantToMultiplex");
70  pfTauDecayModeIndex_token = consumes<PFTauDiscriminator>(pfTauDecayModeIndexSrc_);
71  discriminant_token = consumes<PFTauDiscriminator>(discriminantToMultiplex_);
72  //get the computer/decay mode map
73  std::vector<edm::ParameterSet> decayModeMap = iConfig.getParameter<std::vector<edm::ParameterSet> >("computers");
74  computers_.reserve(decayModeMap.size());
75 
76  // for each decay mode MVA implementation (which may correspond to multiple decay modes, map the decay modes to the correct MVA computer
77  for(std::vector<edm::ParameterSet>::const_iterator iComputer = decayModeMap.begin();
78  iComputer != decayModeMap.end();
79  ++iComputer)
80  {
81  ComputerAndCut toInsert;
82  toInsert.computerName = iComputer->getParameter<std::string>("computerName");
83  toInsert.userCut = iComputer->getParameter<double>("cut");
84  CutList::iterator computerJustAdded = computers_.insert(computers_.end(), toInsert); //add this computer to the end of the list
85 
86  //populate the map
87  std::vector<int> associatedDecayModes = iComputer->getParameter<std::vector<int> >("decayModeIndices");
88  for(std::vector<int>::const_iterator iDecayMode = associatedDecayModes.begin();
89  iDecayMode != associatedDecayModes.end();
90  ++iDecayMode)
91  {
92  //map this integer specifying the decay mode to the MVA comptuer we just added to the list
93  std::pair<DecayModeToCutMap::iterator, bool> insertResult = computerMap_.insert(std::make_pair(*iDecayMode, computerJustAdded));
94 
95  //make sure we aren't double mapping a decay mode
96  if(insertResult.second == false) { //indicates that the current key (decaymode) has already been entered!
97  throw cms::Exception("PFTauDecayModeCutMultiplexer::ctor") << "A tau decay mode: " << *iDecayMode << " has been mapped to two different MVA implementations, "
98  << insertResult.first->second->computerName << " and " << toInsert.computerName
99  << ". Please check the appropriate cfi file." << std::endl;
100  }
101  }
102  }
103 }
104 
105 // ------------ get the relevant decay mode index handles at the beginning of each event ------------
106 void
108 {
109 
112 }
113 
115 {
116  // get decay mode for current tau
117  int decayMode = lrint( (*pfTauDecayModeIndices)[pfTau] ); //convert to int
118 
119  // get value we are trying to multiplex
120  float valueToMultiplex = (*targetDiscriminant)[pfTau];
121 
122  // Get correct cut
123  auto iterToComputer = computerMap_.find(decayMode);
124  if(iterToComputer != computerMap_.end()) //if we don't have a MVA mapped to this decay mode, skip it, it fails.
125  {
126  // use the supplied cut to make a decision
127  if (valueToMultiplex > iterToComputer->second->userCut)
128  return 1.0;
129  else
130  return 0.0;
131  }
132 
133  // no computer associated to this decay mode; it fails
134  return 0.;
135 }
136 
137 void
139  // pfTauDecayModeCutMultiplexer
141 
142  desc.add<edm::InputTag>("PFTauDecayModeSrc");
143  desc.add<edm::InputTag>("PFTauDiscriminantToMultiplex");
144 
145  edm::ParameterSetDescription vpsd_computers;
146  vpsd_computers.add<std::string>("computerName");
147  vpsd_computers.add<double>("cut");
148  vpsd_computers.add<std::vector<int> >("decayModeIndices");
149 
150  // name description defaults items
151  //desc.addVPSet("computers", vpsd_builders, builders_vector);
152  desc.addVPSet("computers", vpsd_computers);
153 
154  fillProducerDescriptions(desc); // inherited from the base
155 
156  descriptions.add("pfTauDecayModeCutMultiplexer", desc);
157 }
158 
T getParameter(std::string const &) const
double discriminate(const PFTauRef &thePFTau) const override
std::vector< ComputerAndCut > CutList
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::map< int, CutList::iterator > DecayModeToCutMap
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillProducerDescriptions(edm::ParameterSetDescription &desc)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< PFTauDiscriminator > discriminant_token
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< PFTauDiscriminator > pfTauDecayModeIndex_token
edm::Handle< PFTauDiscriminator > pfTauDecayModeIndices
fixed size matrix
PFTauDecayModeCutMultiplexer(const edm::ParameterSet &)
void beginEvent(const edm::Event &event, const edm::EventSetup &eventSetup) override
edm::Handle< PFTauDiscriminator > targetDiscriminant
Definition: event.py:1