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 
51 private:
52  // PFTau discriminator continaing the decaymode index of the tau collection
55 
56  // Discriminant to multiplex cut on
59 
60  DecayModeToCutMap computerMap_; //Maps decay mode to MVA implementation
62 
64  pfTauDecayModeIndices; // holds the decay mode indices for the taus in the current event
66  targetDiscriminant; // holds the discirminant values of the discriminant we will multiplex for the current event
67 };
68 
71  pfTauDecayModeIndexSrc_ = iConfig.getParameter<edm::InputTag>("PFTauDecayModeSrc");
72  discriminantToMultiplex_ = iConfig.getParameter<edm::InputTag>("PFTauDiscriminantToMultiplex");
73  pfTauDecayModeIndex_token = consumes<PFTauDiscriminator>(pfTauDecayModeIndexSrc_);
74  discriminant_token = consumes<PFTauDiscriminator>(discriminantToMultiplex_);
75  //get the computer/decay mode map
76  std::vector<edm::ParameterSet> decayModeMap = iConfig.getParameter<std::vector<edm::ParameterSet> >("computers");
77  computers_.reserve(decayModeMap.size());
78 
79  // for each decay mode MVA implementation (which may correspond to multiple decay modes, map the decay modes to the correct MVA computer
80  for (std::vector<edm::ParameterSet>::const_iterator iComputer = decayModeMap.begin(); iComputer != decayModeMap.end();
81  ++iComputer) {
82  ComputerAndCut toInsert;
83  toInsert.computerName = iComputer->getParameter<std::string>("computerName");
84  toInsert.userCut = iComputer->getParameter<double>("cut");
85  CutList::iterator computerJustAdded =
86  computers_.insert(computers_.end(), toInsert); //add this computer to the end of the list
87 
88  //populate the map
89  std::vector<int> associatedDecayModes = iComputer->getParameter<std::vector<int> >("decayModeIndices");
90  for (std::vector<int>::const_iterator iDecayMode = associatedDecayModes.begin();
91  iDecayMode != associatedDecayModes.end();
92  ++iDecayMode) {
93  //map this integer specifying the decay mode to the MVA comptuer we just added to the list
94  std::pair<DecayModeToCutMap::iterator, bool> insertResult =
95  computerMap_.insert(std::make_pair(*iDecayMode, computerJustAdded));
96 
97  //make sure we aren't double mapping a decay mode
98  if (insertResult.second == false) { //indicates that the current key (decaymode) has already been entered!
99  throw cms::Exception("PFTauDecayModeCutMultiplexer::ctor")
100  << "A tau decay mode: " << *iDecayMode << " has been mapped to two different MVA implementations, "
101  << insertResult.first->second->computerName << " and " << toInsert.computerName
102  << ". Please check the appropriate cfi file." << std::endl;
103  }
104  }
105  }
106 }
107 
108 // ------------ get the relevant decay mode index handles at the beginning of each event ------------
112 }
113 
115  // get decay mode for current tau
116  int decayMode = lrint((*pfTauDecayModeIndices)[pfTau]); //convert to int
117 
118  // get value we are trying to multiplex
119  float valueToMultiplex = (*targetDiscriminant)[pfTau];
120 
121  // Get correct cut
122  auto iterToComputer = computerMap_.find(decayMode);
123  if (iterToComputer != computerMap_.end()) //if we don't have a MVA mapped to this decay mode, skip it, it fails.
124  {
125  // use the supplied cut to make a decision
126  if (valueToMultiplex > iterToComputer->second->userCut)
127  return 1.0;
128  else
129  return 0.0;
130  }
131 
132  // no computer associated to this decay mode; it fails
133  return 0.;
134 }
135 
137  // pfTauDecayModeCutMultiplexer
139 
140  desc.add<edm::InputTag>("PFTauDecayModeSrc");
141  desc.add<edm::InputTag>("PFTauDiscriminantToMultiplex");
142 
143  edm::ParameterSetDescription vpsd_computers;
144  vpsd_computers.add<std::string>("computerName");
145  vpsd_computers.add<double>("cut");
146  vpsd_computers.add<std::vector<int> >("decayModeIndices");
147 
148  // name description defaults items
149  //desc.addVPSet("computers", vpsd_builders, builders_vector);
150  desc.addVPSet("computers", vpsd_computers);
151 
152  fillProducerDescriptions(desc); // inherited from the base
153 
154  descriptions.add("pfTauDecayModeCutMultiplexer", desc);
155 }
156 
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< ComputerAndCut > CutList
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::map< int, CutList::iterator > DecayModeToCutMap
double discriminate(const PFTauRef &thePFTau) const override
int iEvent
Definition: GenABIO.cc:224
static void fillProducerDescriptions(edm::ParameterSetDescription &desc)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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