CMS 3D CMS Logo

PFRecoTauDiscriminationAgainstElectronMVA6.cc
Go to the documentation of this file.
1 /* class PFRecoTauDiscriminationAgainstElectronMVA6
2  * created : Nov 2 2015,
3  * revised : ,
4  * Authorss : Fabio Colombo (KIT)
5  */
6 
9 
12 
15 
23 
24 #include <iostream>
25 #include <sstream>
26 #include <fstream>
27 
28 using namespace reco;
29 
31 {
32  public:
35  mva_(),
36  category_output_()
37  {
38  mva_ = std::make_unique<AntiElectronIDMVA6>(cfg);
39 
40  srcGsfElectrons_ = cfg.getParameter<edm::InputTag>("srcGsfElectrons");
41  GsfElectrons_token = consumes<reco::GsfElectronCollection>(srcGsfElectrons_);
42  vetoEcalCracks_ = cfg.getParameter<bool>("vetoEcalCracks");
43 
44  verbosity_ = cfg.getParameter<int>("verbosity");
45 
46  // add category index
47  produces<PFTauDiscriminator>("category");
48  }
49 
50  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
51 
52  double discriminate(const PFTauRef&) const override;
53 
54  void endEvent(edm::Event&) override;
55 
57 
58  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
59 
60 private:
61  bool isInEcalCrack(double) const;
62 
64  std::unique_ptr<AntiElectronIDMVA6> mva_;
65 
70 
71  std::unique_ptr<PFTauDiscriminator> category_output_;
72 
74 
76 };
77 
79 {
80  mva_->beginEvent(evt, es);
81 
82  evt.getByToken(Tau_token, taus_);
83  category_output_.reset(new PFTauDiscriminator(TauRefProd(taus_)));
84 
85  evt.getByToken(GsfElectrons_token, gsfElectrons_);
86 }
87 
89 {
90  double mvaValue = 1.;
91  double category = -1.;
92  bool isGsfElectronMatched = false;
93 
94  float deltaRDummy = 9.9;
95 
96  const float ECALBarrelEndcapEtaBorder = 1.479;
97  float tauEtaAtEcalEntrance = -99.;
98  float sumEtaTimesEnergy = 0.;
99  float sumEnergy = 0.;
100  for ( const auto & pfCandidate : thePFTauRef->signalPFCands() ) {
101  sumEtaTimesEnergy += (pfCandidate->positionAtECALEntrance().eta()*pfCandidate->energy());
102  sumEnergy += pfCandidate->energy();
103  }
104  if ( sumEnergy > 0. ) {
105  tauEtaAtEcalEntrance = sumEtaTimesEnergy/sumEnergy;
106  }
107 
108  float leadChargedPFCandEtaAtEcalEntrance = -99.;
109  float leadChargedPFCandPt = -99.;
110  for ( const auto & pfCandidate : thePFTauRef->signalPFCands() ) {
111  const reco::Track* track = nullptr;
112  if ( pfCandidate->trackRef().isNonnull() ) track = pfCandidate->trackRef().get();
113  else if ( pfCandidate->muonRef().isNonnull() && pfCandidate->muonRef()->innerTrack().isNonnull() ) track = pfCandidate->muonRef()->innerTrack().get();
114  else if ( pfCandidate->muonRef().isNonnull() && pfCandidate->muonRef()->globalTrack().isNonnull() ) track = pfCandidate->muonRef()->globalTrack().get();
115  else if ( pfCandidate->muonRef().isNonnull() && pfCandidate->muonRef()->outerTrack().isNonnull() ) track = pfCandidate->muonRef()->outerTrack().get();
116  else if ( pfCandidate->gsfTrackRef().isNonnull() ) track = pfCandidate->gsfTrackRef().get();
117  if ( track ) {
118  if ( track->pt() > leadChargedPFCandPt ) {
119  leadChargedPFCandEtaAtEcalEntrance = pfCandidate->positionAtECALEntrance().eta();
120  leadChargedPFCandPt = track->pt();
121  }
122  }
123  }
124 
125  if( (*thePFTauRef).leadChargedHadrCand().isNonnull()) {
126 
127  int numSignalGammaCandsInSigCone = 0;
128  const std::vector<reco::CandidatePtr>& signalGammaCands = thePFTauRef->signalGammaCands();
129 
130  for ( const auto & pfGamma : signalGammaCands ) {
131 
132  double dR = deltaR(pfGamma->p4(), thePFTauRef->leadChargedHadrCand()->p4());
133  double signalrad = std::max(0.05, std::min(0.10, 3.0/std::max(1.0, thePFTauRef->pt())));
134 
135  // pfGammas inside the tau signal cone
136  if (dR < signalrad) {
137  numSignalGammaCandsInSigCone += 1;
138  }
139  }
140 
141  // loop over the electrons
142  for ( const auto & theGsfElectron : *gsfElectrons_ ) {
143  if ( theGsfElectron.pt() > 10. ) { // CV: only take electrons above some minimal energy/Pt into account...
144  double deltaREleTau = deltaR(theGsfElectron.p4(), thePFTauRef->p4());
145  deltaRDummy = deltaREleTau;
146  if ( deltaREleTau < 0.3 ) {
147  double mva_match = mva_->MVAValue(*thePFTauRef, theGsfElectron);
148  const reco::PFCandidatePtr& lpfch = thePFTauRef->leadPFChargedHadrCand();
149  bool hasGsfTrack = false;
150  if (lpfch.isNonnull()) {
151  hasGsfTrack = lpfch->gsfTrackRef().isNonnull();
152  }
153  if ( !hasGsfTrack )
154  hasGsfTrack = theGsfElectron.gsfTrack().isNonnull();
155 
157  if ( vetoEcalCracks_ && (isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance)) ) {
158  // add category index
159  category_output_->setValue(tauIndex_, category);
160  // return MVA output value
161  return -99;
162  }
164 
165  if ( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ) { // Barrel
166  if ( numSignalGammaCandsInSigCone == 0 && hasGsfTrack ) {
167  category = 5.;
168  }
169  else if ( numSignalGammaCandsInSigCone >= 1 && hasGsfTrack ) {
170  category = 7.;
171  }
172  } else { // Endcap
173  if ( numSignalGammaCandsInSigCone == 0 && hasGsfTrack ) {
174  category = 13.;
175  }
176  else if ( numSignalGammaCandsInSigCone >= 1 && hasGsfTrack ) {
177  category = 15.;
178  }
179  }
180 
181  mvaValue = std::min(mvaValue, mva_match);
182  isGsfElectronMatched = true;
183  } // deltaR < 0.3
184  } // electron pt > 10
185  } // end of loop over electrons
186 
187  if ( !isGsfElectronMatched ) {
188  mvaValue = mva_->MVAValue(*thePFTauRef);
189  const reco::PFCandidatePtr& lpfch = thePFTauRef->leadPFChargedHadrCand();
190  bool hasGsfTrack = false;
191  if (lpfch.isNonnull()) {
192  hasGsfTrack = lpfch->gsfTrackRef().isNonnull();
193  }
194 
196  if ( vetoEcalCracks_ && (isInEcalCrack(tauEtaAtEcalEntrance) || isInEcalCrack(leadChargedPFCandEtaAtEcalEntrance)) ) {
197  // add category index
198  category_output_->setValue(tauIndex_, category);
199  // return MVA output value
200  return -99;
201  }
203 
204  if ( std::abs(tauEtaAtEcalEntrance) < ECALBarrelEndcapEtaBorder ) { // Barrel
205  if ( numSignalGammaCandsInSigCone == 0 && !hasGsfTrack ) {
206  category = 0.;
207  }
208  else if ( numSignalGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
209  category = 2.;
210  }
211  } else { // Endcap
212  if ( numSignalGammaCandsInSigCone == 0 && !hasGsfTrack ) {
213  category = 8.;
214  }
215  else if ( numSignalGammaCandsInSigCone >= 1 && !hasGsfTrack ) {
216  category = 10.;
217  }
218  }
219  }
220  }
221 
222  if ( verbosity_ ) {
223  edm::LogPrint("PFTauAgainstEleMVA6") << "<PFRecoTauDiscriminationAgainstElectronMVA6::discriminate>:" ;
224  edm::LogPrint("PFTauAgainstEleMVA6") << " tau: Pt = " << thePFTauRef->pt() << ", eta = " << thePFTauRef->eta() << ", phi = " << thePFTauRef->phi();
225  edm::LogPrint("PFTauAgainstEleMVA6") << " deltaREleTau = " << deltaRDummy << ", isGsfElectronMatched = " << isGsfElectronMatched;
226  edm::LogPrint("PFTauAgainstEleMVA6") << " #Prongs = " << thePFTauRef->signalChargedHadrCands().size();
227  edm::LogPrint("PFTauAgainstEleMVA6") << " MVA = " << mvaValue << ", category = " << category;
228  }
229 
230  // add category index
231  category_output_->setValue(tauIndex_, category);
232  // return MVA output value
233  return mvaValue;
234 }
235 
237 {
238  // add all category indices to event
239  evt.put(std::move(category_output_), "category");
240 }
241 
242 bool
244 {
245  double absEta = fabs(eta);
246  return (absEta > 1.460 && absEta < 1.558);
247 }
248 
249 void
251  // pfRecoTauDiscriminationAgainstElectronMVA6
253  desc.add<double>("minMVANoEleMatchWOgWOgsfBL", 0.0);
254  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfTauProducer"));
255  desc.add<double>("minMVANoEleMatchWgWOgsfBL", 0.0);
256  desc.add<std::string>("mvaName_wGwGSF_EC", "gbr_wGwGSF_EC");
257  desc.add<double>("minMVAWgWgsfBL", 0.0);
258  desc.add<std::string>("mvaName_woGwGSF_EC", "gbr_woGwGSF_EC");
259  desc.add<double>("minMVAWOgWgsfEC", 0.0);
260  desc.add<std::string>("mvaName_wGwGSF_BL", "gbr_wGwGSF_BL");
261  desc.add<std::string>("mvaName_woGwGSF_BL", "gbr_woGwGSF_BL");
262  desc.add<bool>("returnMVA", true);
263  desc.add<bool>("loadMVAfromDB", true);
264  {
265  edm::ParameterSetDescription pset_Prediscriminants;
266  pset_Prediscriminants.add<std::string>("BooleanOperator", "and");
267  {
269  psd1.add<double>("cut");
270  psd1.add<edm::InputTag>("Producer");
271  pset_Prediscriminants.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
272  }
273  {
274  // encountered this at
275  // RecoTauTag/Configuration/python/HPSPFTaus_cff.py
277  psd1.add<double>("cut");
278  psd1.add<edm::InputTag>("Producer");
279  pset_Prediscriminants.addOptional<edm::ParameterSetDescription>("decayMode", psd1);
280  }
281  desc.add<edm::ParameterSetDescription>("Prediscriminants", pset_Prediscriminants);
282  }
283  desc.add<std::string>("mvaName_NoEleMatch_woGwoGSF_BL", "gbr_NoEleMatch_woGwoGSF_BL");
284  desc.add<bool>("vetoEcalCracks", true);
285  desc.add<bool>("usePhiAtEcalEntranceExtrapolation", false);
286  desc.add<std::string>("mvaName_NoEleMatch_wGwoGSF_BL", "gbr_NoEleMatch_wGwoGSF_BL");
287  desc.add<double>("minMVANoEleMatchWOgWOgsfEC", 0.0);
288  desc.add<double>("minMVAWOgWgsfBL", 0.0);
289  desc.add<double>("minMVAWgWgsfEC", 0.0);
290  desc.add<int>("verbosity", 0);
291  desc.add<std::string>("mvaName_NoEleMatch_wGwoGSF_EC", "gbr_NoEleMatch_wGwoGSF_EC");
292  desc.add<std::string>("method", "BDTG");
293  desc.add<edm::InputTag>("srcGsfElectrons", edm::InputTag("gedGsfElectrons"));
294  desc.add<std::string>("mvaName_NoEleMatch_woGwoGSF_EC", "gbr_NoEleMatch_woGwoGSF_EC");
295  desc.add<double>("minMVANoEleMatchWgWOgsfEC", 0.0);
296  descriptions.add("pfRecoTauDiscriminationAgainstElectronMVA6", desc);
297 }
298 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
edm::RefProd< TauCollection > TauRefProd
Definition: Tau.h:40
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< reco::GsfElectronCollection > GsfElectrons_token
void beginEvent(const edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:684
double pt() const
track transverse momentum
Definition: TrackBase.h:654
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:480
def move(src, dest)
Definition: eostools.py:511