CMS 3D CMS Logo

PATTauDiscriminationAgainstElectronMVA6.cc
Go to the documentation of this file.
1 /* class PATTauDiscriminationAgainstElectronMVA6
2  * created : Apr 14 2016,
3  * revised : ,
4  * Authorss : Anne-Catherine Le Bihan (IPHC)
5  */
6 
9 
12 
17 
18 #include <iostream>
19 #include <sstream>
20 #include <fstream>
21 
22 using namespace pat;
23 
25 {
26  public:
29  mva_(),
30  category_output_()
31  {
32  mva_ = std::make_unique<AntiElectronIDMVA6>(cfg);
33 
34  srcElectrons = cfg.getParameter<edm::InputTag>("srcElectrons");
35  electronToken = consumes<pat::ElectronCollection>(srcElectrons);
36  vetoEcalCracks_ = cfg.getParameter<bool>("vetoEcalCracks");
37  verbosity_ = ( cfg.exists("verbosity") ) ?
38  cfg.getParameter<int>("verbosity") : 0;
39 
40  // add category index
41  produces<PATTauDiscriminator>("category");
42  }
43 
44  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
45 
46  double discriminate(const TauRef&) const override;
47 
48  void endEvent(edm::Event&) override;
49 
51 
52 private:
53  bool isInEcalCrack(double) const;
54 
56  std::unique_ptr<AntiElectronIDMVA6> mva_;
57 
62 
63  std::unique_ptr<PATTauDiscriminator> category_output_;
64 
66 
68 };
69 
71 {
72  mva_->beginEvent(evt, es);
73 
74  evt.getByToken(Tau_token, taus_);
75  category_output_.reset(new PATTauDiscriminator(TauRefProd(taus_)));
76 
77  evt.getByToken(electronToken, Electrons);
78 }
79 
81 {
82  double mvaValue = 1.;
83  double category = -1.;
84  bool isGsfElectronMatched = false;
85  float deltaRDummy = 9.9;
86  const float ECALBarrelEndcapEtaBorder = 1.479;
87  float tauEtaAtEcalEntrance = theTauRef->etaAtEcalEntrance();
88  float leadChargedPFCandEtaAtEcalEntrance = theTauRef->etaAtEcalEntranceLeadChargedCand();
89 
90  if( (*theTauRef).leadChargedHadrCand().isNonnull()) {
91  int numSignalPFGammaCandsInSigCone = 0;
92  const reco::CandidatePtrVector signalGammaCands = theTauRef->signalGammaCands();
93  for ( const auto & gamma : signalGammaCands ) {
94  double dR = deltaR(gamma->p4(), theTauRef->leadChargedHadrCand()->p4());
95  double signalrad = std::max(0.05, std::min(0.10, 3.0/std::max(1.0, theTauRef->pt())));
96  // gammas inside the tau signal cone
97  if (dR < signalrad) {
98  numSignalPFGammaCandsInSigCone += 1;
99  }
100  }
101  // loop over the electrons
102  for ( const auto & theElectron : *Electrons ) {
103  if ( theElectron.pt() > 10. ) { // CV: only take electrons above some minimal energy/Pt into account...
104  double deltaREleTau = deltaR(theElectron.p4(), theTauRef->p4());
105  deltaRDummy = deltaREleTau;
106  if( deltaREleTau < 0.3 ){
107  double mva_match = mva_->MVAValue(*theTauRef, theElectron);
108  bool hasGsfTrack = false;
109  pat::PackedCandidate const* packedLeadTauCand = dynamic_cast<pat::PackedCandidate const*>(theTauRef->leadChargedHadrCand().get());
110  if( abs(packedLeadTauCand->pdgId()) == 11 )
111  hasGsfTrack = true;
112  if ( !hasGsfTrack )
113  hasGsfTrack = theElectron.gsfTrack().isNonnull();
114 
115  // veto taus that go to Ecal crack
116  if ( vetoEcalCracks_ && (isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance)) ) {
117  // add category index
118  category_output_->setValue(tauIndex_, category);
119  // return MVA output value
120  return -99;
121  }
122  // Veto taus that go to Ecal crack
123  if( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ){ // Barrel
124  if( numSignalPFGammaCandsInSigCone == 0 && hasGsfTrack ){
125  category = 5.;
126  }
127  else if ( numSignalPFGammaCandsInSigCone >= 1 && hasGsfTrack ) {
128  category = 7.;
129  }
130  } else { // Endcap
131  if ( numSignalPFGammaCandsInSigCone == 0 && hasGsfTrack ) {
132  category = 13.;
133  }
134  else if ( numSignalPFGammaCandsInSigCone >= 1 && hasGsfTrack ) {
135  category = 15.;
136  }
137  }
138  mvaValue = std::min(mvaValue, mva_match);
139  isGsfElectronMatched = true;
140  } // deltaR < 0.3
141  } // electron pt > 10
142  } // end of loop over electrons
143 
144  if ( !isGsfElectronMatched ) {
145  mvaValue = mva_->MVAValue(*theTauRef);
146  bool hasGsfTrack = false;
147  pat::PackedCandidate const* packedLeadTauCand = dynamic_cast<pat::PackedCandidate const*>(theTauRef->leadChargedHadrCand().get());
148  if( abs(packedLeadTauCand->pdgId()) == 11 ) hasGsfTrack = true;
149 
150  // veto taus that go to Ecal crack
151  if ( vetoEcalCracks_ && (isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance)) ) {
152  // add category index
153  category_output_->setValue(tauIndex_, category);
154  // return MVA output value
155  return -99;
156  }
157  // veto taus that go to Ecal crack
158  if( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ){ // Barrel
159  if( numSignalPFGammaCandsInSigCone == 0 && !hasGsfTrack ){
160  category = 0.;
161  }
162  else if ( numSignalPFGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
163  category = 2.;
164  }
165  } else { // Endcap
166  if ( numSignalPFGammaCandsInSigCone == 0 && !hasGsfTrack ) {
167  category = 8.;
168  }
169  else if ( numSignalPFGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
170  category = 10.;
171  }
172  }
173  }
174  }
175  if ( verbosity_ ) {
176  edm::LogPrint("PATTauAgainstEleMVA6") << "<PATTauDiscriminationAgainstElectronMVA6::discriminate>:" ;
177  edm::LogPrint("PATTauAgainstEleMVA6") << " tau: Pt = " << theTauRef->pt() << ", eta = " << theTauRef->eta() << ", phi = " << theTauRef->phi();
178  edm::LogPrint("PATTauAgainstEleMVA6") << " deltaREleTau = " << deltaRDummy << ", isGsfElectronMatched = " << isGsfElectronMatched;
179  edm::LogPrint("PATTauAgainstEleMVA6") << " #Prongs = " << theTauRef->signalChargedHadrCands().size();
180  edm::LogPrint("PATTauAgainstEleMVA6") << " MVA = " << mvaValue << ", category = " << category;
181  }
182  // add category index
183  category_output_->setValue(tauIndex_, category);
184  // return MVA output value
185  return mvaValue;
186 }
187 
189 {
190  // add all category indices to event
191  evt.put(std::move(category_output_), "category");
192 }
193 
194 bool
196 {
197  double absEta = fabs(eta);
198  return (absEta > 1.460 && absEta < 1.558);
199 }
200 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
edm::RefProd< TauCollection > TauRefProd
Definition: Tau.h:39
int pdgId() const override
PDG identifier.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: HeavyIon.h:7
edm::EDGetTokenT< pat::ElectronCollection > electronToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:245
T min(T a, T b)
Definition: MathUtil.h:58
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void beginEvent(const edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:510