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&);
44 
45  double discriminate(const TauRef&) const;
46 
47  void endEvent(edm::Event&);
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( reco::CandidatePtrVector::const_iterator gamma = signalGammaCands.begin(); gamma != signalGammaCands.end(); ++gamma ){
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( unsigned int ie = 0; ie < Electrons->size(); ++ie ){
100  const pat::Electron& theElectron = Electrons->at(ie);
101  if ( theElectron.pt() > 10. ) { // CV: only take electrons above some minimal energy/Pt into account...
102  double deltaREleTau = deltaR(theElectron.p4(), theTauRef->p4());
103  deltaRDummy = deltaREleTau;
104  if( deltaREleTau < 0.3 ){
105  double mva_match = mva_->MVAValue(*theTauRef, theElectron);
106  bool hasGsfTrack = 0;
107  pat::PackedCandidate const* packedLeadTauCand = dynamic_cast<pat::PackedCandidate const*>(theTauRef->leadChargedHadrCand().get());
108  if( abs(packedLeadTauCand->pdgId()) == 11 )
109  hasGsfTrack = 1;
110  if ( !hasGsfTrack )
111  hasGsfTrack = theElectron.gsfTrack().isNonnull();
112 
113  // veto taus that go to Ecal crack
114  if ( isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance) ) {
115  // add category index
116  category_output_->setValue(tauIndex_, category);
117  // return MVA output value
118  return -99;
119  }
120  // Veto taus that go to Ecal crack
121  if( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ){ // Barrel
122  if( numSignalPFGammaCandsInSigCone == 0 && hasGsfTrack ){
123  category = 5.;
124  }
125  else if ( numSignalPFGammaCandsInSigCone >= 1 && hasGsfTrack ) {
126  category = 7.;
127  }
128  } else { // Endcap
129  if ( numSignalPFGammaCandsInSigCone == 0 && hasGsfTrack ) {
130  category = 13.;
131  }
132  else if ( numSignalPFGammaCandsInSigCone >= 1 && hasGsfTrack ) {
133  category = 15.;
134  }
135  }
136  mvaValue = std::min(mvaValue, mva_match);
137  isGsfElectronMatched = true;
138  } // deltaR < 0.3
139  } // electron pt > 10
140  } // end of loop over electrons
141 
142  if ( !isGsfElectronMatched ) {
143  mvaValue = mva_->MVAValue(*theTauRef);
144  bool hasGsfTrack = 0;
145  pat::PackedCandidate const* packedLeadTauCand = dynamic_cast<pat::PackedCandidate const*>(theTauRef->leadChargedHadrCand().get());
146  if( abs(packedLeadTauCand->pdgId()) == 11 ) hasGsfTrack = 1;
147 
148  // veto taus that go to Ecal crack
149  if ( isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance) ) {
150  // add category index
151  category_output_->setValue(tauIndex_, category);
152  // return MVA output value
153  return -99;
154  }
155  // veto taus that go to Ecal crack
156  if( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ){ // Barrel
157  if( numSignalPFGammaCandsInSigCone == 0 && !hasGsfTrack ){
158  category = 0.;
159  }
160  else if ( numSignalPFGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
161  category = 2.;
162  }
163  } else { // Endcap
164  if ( numSignalPFGammaCandsInSigCone == 0 && !hasGsfTrack ) {
165  category = 8.;
166  }
167  else if ( numSignalPFGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
168  category = 10.;
169  }
170  }
171  }
172  }
173  if ( verbosity_ ) {
174  edm::LogPrint("PATTauAgainstEleMVA6") << "<PATTauDiscriminationAgainstElectronMVA6::discriminate>:" ;
175  edm::LogPrint("PATTauAgainstEleMVA6") << " tau: Pt = " << theTauRef->pt() << ", eta = " << theTauRef->eta() << ", phi = " << theTauRef->phi();
176  edm::LogPrint("PATTauAgainstEleMVA6") << " deltaREleTau = " << deltaRDummy << ", isGsfElectronMatched = " << isGsfElectronMatched;
177  edm::LogPrint("PATTauAgainstEleMVA6") << " #Prongs = " << theTauRef->signalChargedHadrCands().size();
178  edm::LogPrint("PATTauAgainstEleMVA6") << " MVA = " << mvaValue << ", category = " << category;
179  }
180  // add category index
181  category_output_->setValue(tauIndex_, category);
182  // return MVA output value
183  return mvaValue;
184 }
185 
187 {
188  // add all category indices to event
189  evt.put(std::move(category_output_), "category");
190 }
191 
192 bool
194 {
195  double absEta = fabs(eta);
196  return (absEta > 1.460 && absEta < 1.558);
197 }
198 
void beginEvent(const edm::Event &, const edm::EventSetup &)
T getParameter(std::string const &) const
virtual double pt() const final
transverse momentum
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
edm::RefProd< TauCollection > TauRefProd
Definition: Tau.h:39
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:225
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#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
const_iterator begin() const
Definition: PtrVector.h:130
reco::GsfTrackRef gsfTrack() const
override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster ...
edm::EDGetTokenT< pat::ElectronCollection > electronToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const_iterator end() const
Definition: PtrVector.h:135
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:244
T min(T a, T b)
Definition: MathUtil.h:58
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
Analysis-level electron class.
Definition: Electron.h:52
virtual int pdgId() const
PDG identifier.
def move(src, dest)
Definition: eostools.py:510