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  verbosity_ = ( cfg.exists("verbosity") ) ?
37  cfg.getParameter<int>("verbosity") : 0;
38 
39  // add category index
40  produces<PATTauDiscriminator>("category");
41  }
42 
43  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
44 
45  double discriminate(const TauRef&) const override;
46 
47  void endEvent(edm::Event&) override;
48 
50 
51 private:
52  bool isInEcalCrack(double) const;
53 
55  std::unique_ptr<AntiElectronIDMVA6> mva_;
56 
61 
62  std::unique_ptr<PATTauDiscriminator> category_output_;
63 
65 };
66 
68 {
69  mva_->beginEvent(evt, es);
70 
71  evt.getByToken(Tau_token, taus_);
72  category_output_.reset(new PATTauDiscriminator(TauRefProd(taus_)));
73 
74  evt.getByToken(electronToken, Electrons);
75 }
76 
78 {
79  double mvaValue = 1.;
80  double category = -1.;
81  bool isGsfElectronMatched = false;
82  float deltaRDummy = 9.9;
83  const float ECALBarrelEndcapEtaBorder = 1.479;
84  float tauEtaAtEcalEntrance = theTauRef->etaAtEcalEntrance();
85  float leadChargedPFCandEtaAtEcalEntrance = theTauRef->etaAtEcalEntranceLeadChargedCand();
86 
87  if( (*theTauRef).leadChargedHadrCand().isNonnull()) {
88  int numSignalPFGammaCandsInSigCone = 0;
89  const reco::CandidatePtrVector signalGammaCands = theTauRef->signalGammaCands();
90  for ( const auto & gamma : signalGammaCands ) {
91  double dR = deltaR(gamma->p4(), theTauRef->leadChargedHadrCand()->p4());
92  double signalrad = std::max(0.05, std::min(0.10, 3.0/std::max(1.0, theTauRef->pt())));
93  // gammas inside the tau signal cone
94  if (dR < signalrad) {
95  numSignalPFGammaCandsInSigCone += 1;
96  }
97  }
98  // loop over the electrons
99  for ( const auto & theElectron : *Electrons ) {
100  if ( theElectron.pt() > 10. ) { // CV: only take electrons above some minimal energy/Pt into account...
101  double deltaREleTau = deltaR(theElectron.p4(), theTauRef->p4());
102  deltaRDummy = deltaREleTau;
103  if( deltaREleTau < 0.3 ){
104  double mva_match = mva_->MVAValue(*theTauRef, theElectron);
105  bool hasGsfTrack = false;
106  pat::PackedCandidate const* packedLeadTauCand = dynamic_cast<pat::PackedCandidate const*>(theTauRef->leadChargedHadrCand().get());
107  if( abs(packedLeadTauCand->pdgId()) == 11 )
108  hasGsfTrack = true;
109  if ( !hasGsfTrack )
110  hasGsfTrack = theElectron.gsfTrack().isNonnull();
111 
112  // veto taus that go to Ecal crack
113  if ( isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance) ) {
114  // add category index
115  category_output_->setValue(tauIndex_, category);
116  // return MVA output value
117  return -99;
118  }
119  // Veto taus that go to Ecal crack
120  if( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ){ // Barrel
121  if( numSignalPFGammaCandsInSigCone == 0 && hasGsfTrack ){
122  category = 5.;
123  }
124  else if ( numSignalPFGammaCandsInSigCone >= 1 && hasGsfTrack ) {
125  category = 7.;
126  }
127  } else { // Endcap
128  if ( numSignalPFGammaCandsInSigCone == 0 && hasGsfTrack ) {
129  category = 13.;
130  }
131  else if ( numSignalPFGammaCandsInSigCone >= 1 && hasGsfTrack ) {
132  category = 15.;
133  }
134  }
135  mvaValue = std::min(mvaValue, mva_match);
136  isGsfElectronMatched = true;
137  } // deltaR < 0.3
138  } // electron pt > 10
139  } // end of loop over electrons
140 
141  if ( !isGsfElectronMatched ) {
142  mvaValue = mva_->MVAValue(*theTauRef);
143  bool hasGsfTrack = false;
144  pat::PackedCandidate const* packedLeadTauCand = dynamic_cast<pat::PackedCandidate const*>(theTauRef->leadChargedHadrCand().get());
145  if( abs(packedLeadTauCand->pdgId()) == 11 ) hasGsfTrack = true;
146 
147  // veto taus that go to Ecal crack
148  if ( isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance) ) {
149  // add category index
150  category_output_->setValue(tauIndex_, category);
151  // return MVA output value
152  return -99;
153  }
154  // veto taus that go to Ecal crack
155  if( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ){ // Barrel
156  if( numSignalPFGammaCandsInSigCone == 0 && !hasGsfTrack ){
157  category = 0.;
158  }
159  else if ( numSignalPFGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
160  category = 2.;
161  }
162  } else { // Endcap
163  if ( numSignalPFGammaCandsInSigCone == 0 && !hasGsfTrack ) {
164  category = 8.;
165  }
166  else if ( numSignalPFGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
167  category = 10.;
168  }
169  }
170  }
171  }
172  if ( verbosity_ ) {
173  edm::LogPrint("PATTauAgainstEleMVA6") << "<PATTauDiscriminationAgainstElectronMVA6::discriminate>:" ;
174  edm::LogPrint("PATTauAgainstEleMVA6") << " tau: Pt = " << theTauRef->pt() << ", eta = " << theTauRef->eta() << ", phi = " << theTauRef->phi();
175  edm::LogPrint("PATTauAgainstEleMVA6") << " deltaREleTau = " << deltaRDummy << ", isGsfElectronMatched = " << isGsfElectronMatched;
176  edm::LogPrint("PATTauAgainstEleMVA6") << " #Prongs = " << theTauRef->signalChargedHadrCands().size();
177  edm::LogPrint("PATTauAgainstEleMVA6") << " MVA = " << mvaValue << ", category = " << category;
178  }
179  // add category index
180  category_output_->setValue(tauIndex_, category);
181  // return MVA output value
182  return mvaValue;
183 }
184 
186 {
187  // add all category indices to event
188  evt.put(std::move(category_output_), "category");
189 }
190 
191 bool
193 {
194  double absEta = fabs(eta);
195  return (absEta > 1.460 && absEta < 1.558);
196 }
197 
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