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 
28 using namespace reco;
29 
31  public:
34 
35  struct ComputerAndCut {
37  double userCut;
38  };
39 
40  typedef std::vector<ComputerAndCut> CutList;
41  typedef std::map<int, CutList::iterator> DecayModeToCutMap;
42 
43  double discriminate(const PFTauRef& thePFTau) const override;
44  void beginEvent(const edm::Event& event, const edm::EventSetup& eventSetup) override;
45 
46  private:
47  // PFTau discriminator continaing the decaymode index of the tau collection
50 
51  // Discriminant to multiplex cut on
54 
55  DecayModeToCutMap computerMap_; //Maps decay mode to MVA implementation
56  CutList computers_;
57 
58  edm::Handle<PFTauDiscriminator> pfTauDecayModeIndices; // holds the decay mode indices for the taus in the current event
59  edm::Handle<PFTauDiscriminator> targetDiscriminant; // holds the discirminant values of the discriminant we will multiplex for the current event
60 };
61 
63 {
64  pfTauDecayModeIndexSrc_ = iConfig.getParameter<edm::InputTag>("PFTauDecayModeSrc");
65  discriminantToMultiplex_ = iConfig.getParameter<edm::InputTag>("PFTauDiscriminantToMultiplex");
66  pfTauDecayModeIndex_token = consumes<PFTauDiscriminator>(pfTauDecayModeIndexSrc_);
67  discriminant_token = consumes<PFTauDiscriminator>(discriminantToMultiplex_);
68  //get the computer/decay mode map
69  std::vector<edm::ParameterSet> decayModeMap = iConfig.getParameter<std::vector<edm::ParameterSet> >("computers");
70  computers_.reserve(decayModeMap.size());
71 
72  // for each decay mode MVA implementation (which may correspond to multiple decay modes, map the decay modes to the correct MVA computer
73  for(std::vector<edm::ParameterSet>::const_iterator iComputer = decayModeMap.begin();
74  iComputer != decayModeMap.end();
75  ++iComputer)
76  {
77  ComputerAndCut toInsert;
78  toInsert.computerName = iComputer->getParameter<std::string>("computerName");
79  toInsert.userCut = iComputer->getParameter<double>("cut");
80  CutList::iterator computerJustAdded = computers_.insert(computers_.end(), toInsert); //add this computer to the end of the list
81 
82  //populate the map
83  std::vector<int> associatedDecayModes = iComputer->getParameter<std::vector<int> >("decayModeIndices");
84  for(std::vector<int>::const_iterator iDecayMode = associatedDecayModes.begin();
85  iDecayMode != associatedDecayModes.end();
86  ++iDecayMode)
87  {
88  //map this integer specifying the decay mode to the MVA comptuer we just added to the list
89  std::pair<DecayModeToCutMap::iterator, bool> insertResult = computerMap_.insert(std::make_pair(*iDecayMode, computerJustAdded));
90 
91  //make sure we aren't double mapping a decay mode
92  if(insertResult.second == false) { //indicates that the current key (decaymode) has already been entered!
93  throw cms::Exception("PFTauDecayModeCutMultiplexer::ctor") << "A tau decay mode: " << *iDecayMode << " has been mapped to two different MVA implementations, "
94  << insertResult.first->second->computerName << " and " << toInsert.computerName
95  << ". Please check the appropriate cfi file." << std::endl;
96  }
97  }
98  }
99 }
100 
101 // ------------ get the relevant decay mode index handles at the beginning of each event ------------
102 void
104 {
105 
108 }
109 
111 {
112  // get decay mode for current tau
113  int decayMode = lrint( (*pfTauDecayModeIndices)[pfTau] ); //convert to int
114 
115  // get value we are trying to multiplex
116  float valueToMultiplex = (*targetDiscriminant)[pfTau];
117 
118  // Get correct cut
119  auto iterToComputer = computerMap_.find(decayMode);
120  if(iterToComputer != computerMap_.end()) //if we don't have a MVA mapped to this decay mode, skip it, it fails.
121  {
122  // use the supplied cut to make a decision
123  if (valueToMultiplex > iterToComputer->second->userCut)
124  return 1.0;
125  else
126  return 0.0;
127  }
128 
129  // no computer associated to this decay mode; it fails
130  return 0.;
131 }
132 
T getParameter(std::string const &) const
double discriminate(const PFTauRef &thePFTau) const override
std::vector< ComputerAndCut > CutList
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::map< int, CutList::iterator > DecayModeToCutMap
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< PFTauDiscriminator > discriminant_token
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