CMS 3D CMS Logo

PFRecoTauDiscriminationByIsolation.cc
Go to the documentation of this file.
1 #include <functional>
10 
11 #include "TMath.h"
12 #include "TFormula.h"
13 
14 /* class PFRecoTauDiscriminationByIsolation
15  * created : Jul 23 2007,
16  * revised : Thu Aug 13 14:44:40 PDT 2009
17  * contributors : Ludovic Houchu (IPHC, Strasbourg),
18  * Christian Veelken (UC Davis),
19  * Evan K. Friis (UC Davis)
20  * Michalis Bachtis (UW Madison)
21  */
22 
23 using namespace reco;
24 using namespace std;
25 
27 {
28  public:
31  moduleLabel_(pset.getParameter<std::string>("@module_label")),
32  qualityCutsPSet_(pset.getParameter<edm::ParameterSet>("qualityCuts"))
33  {
34  includeTracks_ = pset.getParameter<bool>(
35  "ApplyDiscriminationByTrackerIsolation");
36  includeGammas_ = pset.getParameter<bool>(
37  "ApplyDiscriminationByECALIsolation");
38 
39  calculateWeights_ = pset.exists("ApplyDiscriminationByWeightedECALIsolation") ?
40  pset.getParameter<bool>("ApplyDiscriminationByWeightedECALIsolation") : false;
41 
42  // RIC: multiply neutral isolation by a flat factor.
43  // Useful, for instance, to combine charged and neutral isolations
44  // with different relative weights
45  weightGammas_ = pset.exists("WeightECALIsolation") ?
46  pset.getParameter<double>("WeightECALIsolation") : 1.0;
47 
48  // RIC: allow to relax the isolation completely beyond a given tau pt
49  minPtForNoIso_ = pset.exists("minTauPtForNoIso") ?
50  pset.getParameter<double>("minTauPtForNoIso") : -99.;
51 
52  applyOccupancyCut_ = pset.getParameter<bool>("applyOccupancyCut");
53  maximumOccupancy_ = pset.getParameter<uint32_t>("maximumOccupancy");
54 
55  applySumPtCut_ = pset.getParameter<bool>("applySumPtCut");
56  maximumSumPt_ = pset.getParameter<double>("maximumSumPtCut");
57 
58  applyRelativeSumPtCut_ = pset.getParameter<bool>(
59  "applyRelativeSumPtCut");
60  maximumRelativeSumPt_ = pset.getParameter<double>(
61  "relativeSumPtCut");
62  offsetRelativeSumPt_ = pset.exists("relativeSumPtOffset") ?
63  pset.getParameter<double>("relativeSumPtOffset") : 0.0;
64 
65  storeRawOccupancy_ = pset.exists("storeRawOccupancy") ?
66  pset.getParameter<bool>("storeRawOccupancy") : false;
67  storeRawSumPt_ = pset.exists("storeRawSumPt") ?
68  pset.getParameter<bool>("storeRawSumPt") : false;
69  storeRawPUsumPt_ = pset.exists("storeRawPUsumPt") ?
70  pset.getParameter<bool>("storeRawPUsumPt") : false;
71  storeRawFootprintCorrection_ = pset.exists("storeRawFootprintCorrection") ?
72  pset.getParameter<bool>("storeRawFootprintCorrection") : false;
73  storeRawPhotonSumPt_outsideSignalCone_ = pset.exists("storeRawPhotonSumPt_outsideSignalCone") ?
74  pset.getParameter<bool>("storeRawPhotonSumPt_outsideSignalCone") : false;
75 
76  // Sanity check on requested options. We can't apply cuts and store the
77  // raw output at the same time
78  if ( applySumPtCut_ || applyOccupancyCut_ || applyRelativeSumPtCut_ ) {
79  if ( storeRawSumPt_ || storeRawOccupancy_ || storeRawPUsumPt_ ) {
80  throw cms::Exception("BadIsoConfig")
81  << "A 'store raw' and a 'apply cut' option have been set to true "
82  << "simultaneously. These options are mutually exclusive.";
83  }
84  }
85 
86  // sanity check2 - can't use weighted and unweighted iso at the same time
87  if ( includeGammas_ && calculateWeights_ ) {
88  throw cms::Exception("BasIsoConfig")
89  << "Both 'ApplyDiscriminationByECALIsolation' and 'ApplyDiscriminationByWeightedECALIsolation' "
90  << "have been set to true. These options are mutually exclusive.";
91  }
92 
93  // Can only store one type
94  int numStoreOptions = 0;
95  if ( storeRawSumPt_ ) ++numStoreOptions;
96  if ( storeRawOccupancy_ ) ++numStoreOptions;
97  if ( storeRawPUsumPt_ ) ++numStoreOptions;
98  if ( storeRawFootprintCorrection_ ) ++numStoreOptions;
99  if ( storeRawPhotonSumPt_outsideSignalCone_ ) ++numStoreOptions;
100  if ( numStoreOptions > 1 ) {
101  throw cms::Exception("BadIsoConfig")
102  << "Multiple 'store sum pt' and/or 'store occupancy' options are set."
103  << " These options are mutually exclusive.";
104  }
105 
106  if ( pset.exists("customOuterCone") ) {
107  customIsoCone_ = pset.getParameter<double>("customOuterCone");
108  } else {
109  customIsoCone_ = -1;
110  }
111 
112  applyPhotonPtSumOutsideSignalConeCut_ = ( pset.exists("applyPhotonPtSumOutsideSignalConeCut") ) ?
113  pset.getParameter<bool>("applyPhotonPtSumOutsideSignalConeCut") : false;
114  if ( applyPhotonPtSumOutsideSignalConeCut_ ) {
115  maxAbsPhotonSumPt_outsideSignalCone_ = pset.getParameter<double>("maxAbsPhotonSumPt_outsideSignalCone");
116  maxRelPhotonSumPt_outsideSignalCone_ = pset.getParameter<double>("maxRelPhotonSumPt_outsideSignalCone");
117  }
118 
119  applyFootprintCorrection_ = ( pset.exists("applyFootprintCorrection") ) ?
120  pset.getParameter<bool>("applyFootprintCorrection") : false;
121  if ( applyFootprintCorrection_ || storeRawFootprintCorrection_ ) {
122  edm::VParameterSet cfgFootprintCorrections = pset.getParameter<edm::VParameterSet>("footprintCorrections");
123  for ( edm::VParameterSet::const_iterator cfgFootprintCorrection = cfgFootprintCorrections.begin();
124  cfgFootprintCorrection != cfgFootprintCorrections.end(); ++cfgFootprintCorrection ) {
125  std::string selection = cfgFootprintCorrection->getParameter<std::string>("selection");
126  std::string offset = cfgFootprintCorrection->getParameter<std::string>("offset");
127  std::unique_ptr<FootprintCorrection> footprintCorrection(new FootprintCorrection(selection, offset));
128  footprintCorrections_.push_back(std::move(footprintCorrection));
129  }
130  }
131 
132  // Get the quality cuts specific to the isolation region
133  edm::ParameterSet isolationQCuts = qualityCutsPSet_.getParameterSet(
134  "isolationQualityCuts");
135 
136  qcuts_.reset(new tau::RecoTauQualityCuts(isolationQCuts));
137 
138  vertexAssociator_.reset(
139  new tau::RecoTauVertexAssociator(qualityCutsPSet_,consumesCollector()));
140 
141  applyDeltaBeta_ = pset.exists("applyDeltaBetaCorrection") ?
142  pset.getParameter<bool>("applyDeltaBetaCorrection") : false;
143 
144  if ( applyDeltaBeta_ || calculateWeights_ ) {
145  // Factorize the isolation QCuts into those that are used to
146  // select PU and those that are not.
147  std::pair<edm::ParameterSet, edm::ParameterSet> puFactorizedIsoQCuts =
148  reco::tau::factorizePUQCuts(isolationQCuts);
149 
150  // Determine the pt threshold for the PU tracks
151  // First check if the user specifies explicitly the cut.
152  if ( pset.exists("deltaBetaPUTrackPtCutOverride") ) {
153  puFactorizedIsoQCuts.second.addParameter<double>(
154  "minTrackPt",
155  pset.getParameter<double>("deltaBetaPUTrackPtCutOverride"));
156  } else {
157  // Secondly take it from the minGammaEt
158  puFactorizedIsoQCuts.second.addParameter<double>(
159  "minTrackPt",
160  isolationQCuts.getParameter<double>("minGammaEt"));
161  }
162 
163  pileupQcutsPUTrackSelection_.reset(new tau::RecoTauQualityCuts(
164  puFactorizedIsoQCuts.first));
165 
166  pileupQcutsGeneralQCuts_.reset(new tau::RecoTauQualityCuts(
167  puFactorizedIsoQCuts.second));
168 
169  pfCandSrc_ = pset.getParameter<edm::InputTag>("particleFlowSrc");
170  pfCand_token = consumes<reco::PFCandidateCollection>(pfCandSrc_);
171  vertexSrc_ = pset.getParameter<edm::InputTag>("vertexSrc");
172  vertex_token = consumes<reco::VertexCollection>(vertexSrc_);
173  deltaBetaCollectionCone_ = pset.getParameter<double>(
174  "isoConeSizeForDeltaBeta");
175  std::string deltaBetaFactorFormula =
176  pset.getParameter<string>("deltaBetaFactor");
177  deltaBetaFormula_.reset(
178  new TFormula("DB_corr", deltaBetaFactorFormula.c_str()));
179  }
180 
181  applyRhoCorrection_ = pset.exists("applyRhoCorrection") ?
182  pset.getParameter<bool>("applyRhoCorrection") : false;
183  if ( applyRhoCorrection_ ) {
184  rhoProducer_ = pset.getParameter<edm::InputTag>("rhoProducer");
185  rho_token=consumes<double>(rhoProducer_);
186  rhoConeSize_ = pset.getParameter<double>("rhoConeSize");
187  rhoUEOffsetCorrection_ =
188  pset.getParameter<double>("rhoUEOffsetCorrection");
189  }
190  useAllPFCands_ = pset.exists("UseAllPFCandsForWeights") ?
191  pset.getParameter<bool>("UseAllPFCandsForWeights") : false;
192 
193  verbosity_ = ( pset.exists("verbosity") ) ?
194  pset.getParameter<int>("verbosity") : 0;
195  }
196 
198  {
199  }
200 
201  void beginEvent(const edm::Event& evt, const edm::EventSetup& evtSetup) override;
202  double discriminate(const PFTauRef& pfTau) const override;
203 
204  inline double weightedSum(const std::vector<PFCandidatePtr>& inColl_, double eta, double phi) const {
205  double out = 1.0;
206  for (auto const & inObj_ : inColl_){
207  double sum = (inObj_->pt()*inObj_->pt())/(deltaR2(eta,phi,inObj_->eta(),inObj_->phi()));
208  if(sum > 1.0) out *= sum;
209  }
210  return out;
211  }
212 
213  private:
215 
217  std::auto_ptr<tau::RecoTauQualityCuts> qcuts_;
218 
219  // Inverted QCut which selects tracks with bad DZ/trackWeight
220  std::auto_ptr<tau::RecoTauQualityCuts> pileupQcutsPUTrackSelection_;
221  std::auto_ptr<tau::RecoTauQualityCuts> pileupQcutsGeneralQCuts_;
222 
223  std::auto_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;
224 
237  // RIC:
239 
243 
246  {
248  : selection_(selection),
249  offset_(offset)
250  {}
254  };
255  std::vector<std::unique_ptr<FootprintCorrection> > footprintCorrections_;
256 
257  // Options to store the raw value in the discriminator instead of boolean pass/fail flag
263 
264  /* **********************************************************************
265  **** Pileup Subtraction Parameters ***********************************
266  **********************************************************************/
267 
268  // Delta Beta correction
272  // Keep track of how many vertices are in the event
275  std::vector<reco::PFCandidatePtr> chargedPFCandidatesInEvent_;
276  // Size of cone used to collect PU tracks
278  std::auto_ptr<TFormula> deltaBetaFormula_;
280 
281  // Rho correction
286  double rhoConeSize_;
290 
291  // Flag to enable/disable debug output
293 };
294 
296 {
297  // NB: The use of the PV in this context is necessitated by its use in
298  // applying quality cuts to the different objects in the isolation cone
299  // The vertex associator contains the logic to select the appropriate vertex
300  // We need to pass it the event so it can load the vertices.
301  vertexAssociator_->setEvent(event);
302 
303  // If we are applying the delta beta correction, we need to get the PF
304  // candidates from the event so we can find the PU tracks.
305  if ( applyDeltaBeta_ || calculateWeights_ ) {
306  // Collect all the PF pile up tracks
308  event.getByToken(pfCand_token, pfCandidates);
309  chargedPFCandidatesInEvent_.clear();
310  chargedPFCandidatesInEvent_.reserve(pfCandidates->size());
311  size_t numPFCandidates = pfCandidates->size();
312  for ( size_t i = 0; i < numPFCandidates; ++i ) {
313  reco::PFCandidatePtr pfCandidate(pfCandidates, i);
314  if ( pfCandidate->charge() != 0 ) {
315  chargedPFCandidatesInEvent_.push_back(pfCandidate);
316  }
317  }
318  // Count all the vertices in the event, to parameterize the DB
319  // correction factor
321  event.getByToken(vertex_token, vertices);
322  size_t nVtxThisEvent = vertices->size();
323  deltaBetaFactorThisEvent_ = deltaBetaFormula_->Eval(nVtxThisEvent);
324  }
325 
326  if ( applyRhoCorrection_ ) {
327  edm::Handle<double> rhoHandle_;
328  event.getByToken(rho_token, rhoHandle_);
329  rhoThisEvent_ = (*rhoHandle_ - rhoUEOffsetCorrection_)*
330  (3.14159)*rhoConeSize_*rhoConeSize_;
331  }
332 }
333 
334 double
336 {
337  LogDebug("discriminate") << " tau: Pt = " << pfTau->pt() << ", eta = " << pfTau->eta() << ", phi = " << pfTau->phi();
338  LogDebug("discriminate") << *pfTau ;
339 
340  // collect the objects we are working with (ie tracks, tracks+gammas, etc)
341  std::vector<PFCandidatePtr> isoCharged_;
342  std::vector<PFCandidatePtr> isoNeutral_;
343  std::vector<PFCandidatePtr> isoPU_;
344  PFCandidateCollection isoNeutralWeight_;
345  std::vector<PFCandidatePtr> chPV_;
346  isoCharged_.reserve(pfTau->isolationPFChargedHadrCands().size());
347  isoNeutral_.reserve(pfTau->isolationPFGammaCands().size());
348  isoPU_.reserve(std::min(100UL, chargedPFCandidatesInEvent_.size()));
349  isoNeutralWeight_.reserve(pfTau->isolationPFGammaCands().size());
350 
351  chPV_.reserve(std::min(50UL, chargedPFCandidatesInEvent_.size()));
352 
353  // Get the primary vertex associated to this tau
354  reco::VertexRef pv = vertexAssociator_->associatedVertex(*pfTau);
355  // Let the quality cuts know which the vertex to use when applying selections
356  // on dz, etc.
357  if ( verbosity_ ) {
358  if ( pv.isNonnull() ) {
359  LogTrace("discriminate") << "pv: x = " << pv->position().x() << ", y = " << pv->position().y() << ", z = " << pv->position().z() ;
360  } else {
361  LogTrace("discriminate") << "pv: N/A" ;
362  }
363  if ( pfTau->leadPFChargedHadrCand().isNonnull() ) {
364  LogTrace("discriminate") << "leadPFChargedHadron:"
365  << " Pt = " << pfTau->leadPFChargedHadrCand()->pt() << ","
366  << " eta = " << pfTau->leadPFChargedHadrCand()->eta() << ","
367  << " phi = " << pfTau->leadPFChargedHadrCand()->phi() ;
368  } else {
369  LogTrace("discriminate") << "leadPFChargedHadron: N/A" ;
370  }
371  }
372 
373  // CV: isolation is not well defined in case primary vertex or leading charged hadron do not exist
374  if ( !(pv.isNonnull() && pfTau->leadPFChargedHadrCand().isNonnull()) ) return 0.;
375 
376  qcuts_->setPV(pv);
377  qcuts_->setLeadTrack(pfTau->leadPFChargedHadrCand());
378 
379  if ( applyDeltaBeta_ || calculateWeights_) {
380  pileupQcutsGeneralQCuts_->setPV(pv);
381  pileupQcutsGeneralQCuts_->setLeadTrack(pfTau->leadPFChargedHadrCand());
382  pileupQcutsPUTrackSelection_->setPV(pv);
383  pileupQcutsPUTrackSelection_->setLeadTrack(pfTau->leadPFChargedHadrCand());
384  }
385 
386  // Load the tracks if they are being used.
387  if ( includeTracks_ ) {
388  for( auto const & cand : pfTau->isolationPFChargedHadrCands() ) {
389  if ( qcuts_->filterCandRef(cand) ) {
390  LogTrace("discriminate") << "adding charged iso cand with pt " << cand->pt() ;
391  isoCharged_.push_back(cand);
392  }
393  }
394  }
395  if ( includeGammas_ || calculateWeights_ ) {
396  for( auto const & cand : pfTau->isolationPFGammaCands() ) {
397  if ( qcuts_->filterCandRef(cand) ) {
398  LogTrace("discriminate") << "adding neutral iso cand with pt " << cand->pt() ;
399  isoNeutral_.push_back(cand);
400  }
401  }
402  }
403 
406 
407  // If desired, get PU tracks.
408  if ( applyDeltaBeta_ || calculateWeights_) {
409  // First select by inverted the DZ/track weight cuts. True = invert
410  if ( verbosity_ ) {
411  std::cout << "Initial PFCands: " << chargedPFCandidatesInEvent_.size() << std::endl;
412  }
413 
414  std::vector<PFCandidatePtr> allPU =
415  pileupQcutsPUTrackSelection_->filterCandRefs(
416  chargedPFCandidatesInEvent_, true);
417 
418  std::vector<PFCandidatePtr> allNPU =
419  pileupQcutsPUTrackSelection_->filterCandRefs(
420  chargedPFCandidatesInEvent_);
421  LogTrace("discriminate") << "After track cuts: " << allPU.size() ;
422 
423  // Now apply the rest of the cuts, like pt, and TIP, tracker hits, etc
424  if ( !useAllPFCands_ ) {
425  std::vector<PFCandidatePtr> cleanPU =
426  pileupQcutsGeneralQCuts_->filterCandRefs(allPU);
427 
428  std::vector<PFCandidatePtr> cleanNPU =
429  pileupQcutsGeneralQCuts_->filterCandRefs(allNPU);
430 
431  LogTrace("discriminate") << "After cleaning cuts: " << cleanPU.size() ;
432 
433  // Only select PU tracks inside the isolation cone.
434  DRFilter deltaBetaFilter(pfTau->p4(), 0, deltaBetaCollectionCone_);
435  for ( auto const & cand : cleanPU ) {
436  if ( deltaBetaFilter(cand) ) isoPU_.push_back(cand);
437  }
438 
439  for ( auto const & cand : cleanNPU ) {
440  if ( deltaBetaFilter(cand) ) chPV_.push_back(cand);
441  }
442  LogTrace("discriminate") << "After cone cuts: " << isoPU_.size() << " " << chPV_.size() ;
443  } else {
444  isoPU_ = std::move(allPU);
445  chPV_ = std::move(allNPU);
446  }
447  }
448 
449  if ( calculateWeights_ ) {
450  for ( auto const & isoObject : isoNeutral_ ) {
451  if ( isoObject->charge() != 0 ) {
452  // weight only neutral objects
453  isoNeutralWeight_.push_back(*isoObject);
454  continue;
455  }
456 
457  double eta = isoObject->eta();
458  double phi = isoObject->phi();
459  double sumNPU = 0.5*log(weightedSum(chPV_, eta, phi));
460 
461  double sumPU = 0.5*log(weightedSum(isoPU_, eta, phi));
462  PFCandidate neutral = (*isoObject);
463  if ( (sumNPU + sumPU) > 0 ) neutral.setP4(((sumNPU)/(sumNPU + sumPU))*neutral.p4());
464 
465  isoNeutralWeight_.push_back(neutral);
466  }
467  }
468 
469  // Check if we want a custom iso cone
470  if ( customIsoCone_ >= 0. ) {
471  DRFilter filter(pfTau->p4(), 0, customIsoCone_);
472  DRFilter2 filter2(pfTau->p4(), 0, customIsoCone_);
473  std::vector<PFCandidatePtr> isoCharged_filter;
474  std::vector<PFCandidatePtr> isoNeutral_filter;
475  PFCandidateCollection isoNeutralWeight_filter;
476  // Remove all the objects not in our iso cone
477  for( auto const & isoObject : isoCharged_ ) {
478  if ( filter(isoObject) ) isoCharged_filter.push_back(isoObject);
479  }
480  if(!calculateWeights_){
481  for( auto const & isoObject : isoNeutral_ ) {
482  if ( filter(isoObject) ) isoNeutral_filter.push_back(isoObject);
483  }
484  isoNeutral_ = isoNeutral_filter;
485  }else{
486  for( auto const & isoObject : isoNeutralWeight_){
487  if ( filter2(isoObject) ) isoNeutralWeight_filter.push_back(isoObject);
488  }
489  isoNeutralWeight_ = isoNeutralWeight_filter;
490  }
491  isoCharged_ = isoCharged_filter;
492  }
493 
494  bool failsOccupancyCut = false;
495  bool failsSumPtCut = false;
496  bool failsRelativeSumPtCut = false;
497 
498 //--- nObjects requirement
499  int neutrals = isoNeutral_.size();
500 
501  if ( applyDeltaBeta_ ) {
502  neutrals -= TMath::Nint(deltaBetaFactorThisEvent_*isoPU_.size());
503  }
504  if ( neutrals < 0 ) {
505  neutrals = 0;
506  }
507 
508  size_t nOccupants = isoCharged_.size() + neutrals;
509 
510  failsOccupancyCut = ( nOccupants > maximumOccupancy_ );
511 
512  double footprintCorrection_value = 0.;
513  if ( applyFootprintCorrection_ || storeRawFootprintCorrection_ ) {
514  for ( std::vector<std::unique_ptr<FootprintCorrection> >::const_iterator footprintCorrection = footprintCorrections_.begin();
515  footprintCorrection != footprintCorrections_.end(); ++footprintCorrection ) {
516  if ( (*footprintCorrection)->selection_(*pfTau) ) {
517  footprintCorrection_value = (*footprintCorrection)->offset_(*pfTau);
518  }
519  }
520  }
521 
522  double totalPt = 0.;
523  double puPt = 0.;
524 //--- Sum PT requirement
525  if ( applySumPtCut_ || applyRelativeSumPtCut_ || storeRawSumPt_ || storeRawPUsumPt_ ) {
526  double chargedPt = 0.;
527  double neutralPt = 0.;
528  double weightedNeutralPt = 0.;
529  for ( auto const & isoObject : isoCharged_ ) {
530  chargedPt += isoObject->pt();
531  }
532  if ( !calculateWeights_ ) {
533  for ( auto const & isoObject : isoNeutral_ ) {
534  neutralPt += isoObject->pt();
535  }
536  } else {
537  for ( auto const & isoObject : isoNeutralWeight_ ) {
538  weightedNeutralPt += isoObject.pt();
539  }
540  }
541  for ( auto const & isoObject : isoPU_ ) {
542  puPt += isoObject->pt();
543  }
544  LogTrace("discriminate") << "chargedPt = " << chargedPt ;
545  LogTrace("discriminate") << "neutralPt = " << neutralPt ;
546  LogTrace("discriminate") << "weighted neutral Pt = " << weightedNeutralPt ;
547  LogTrace("discriminate") << "puPt = " << puPt << " (delta-beta corr. = " << (deltaBetaFactorThisEvent_*puPt) << ")" ;
548 
549  if ( calculateWeights_ ) {
550  neutralPt = weightedNeutralPt;
551  }
552 
553  if ( applyDeltaBeta_ ) {
554  neutralPt -= (deltaBetaFactorThisEvent_*puPt);
555  }
556 
557  if ( applyFootprintCorrection_ ) {
558  neutralPt -= footprintCorrection_value;
559  }
560 
561  if ( applyRhoCorrection_ ) {
562  neutralPt -= rhoThisEvent_;
563  }
564 
565  if ( neutralPt < 0. ) {
566  neutralPt = 0.;
567  }
568 
569  totalPt = chargedPt + weightGammas_ * neutralPt;
570  LogTrace("discriminate") << "totalPt = " << totalPt << " (cut = " << maximumSumPt_ << ")" ;
571 
572  failsSumPtCut = (totalPt > maximumSumPt_);
573 
574 //--- Relative Sum PT requirement
575  failsRelativeSumPtCut = (totalPt > ((pfTau->pt() - offsetRelativeSumPt_)*maximumRelativeSumPt_));
576  }
577 
578  bool failsPhotonPtSumOutsideSignalConeCut = false;
579  double photonSumPt_outsideSignalCone = 0.;
580  if ( applyPhotonPtSumOutsideSignalConeCut_ || storeRawPhotonSumPt_outsideSignalCone_ ) {
581  const std::vector<reco::PFCandidatePtr>& signalPFGammas = pfTau->signalPFGammaCands();
582  for ( std::vector<reco::PFCandidatePtr>::const_iterator signalPFGamma = signalPFGammas.begin();
583  signalPFGamma != signalPFGammas.end(); ++signalPFGamma ) {
584  double dR = deltaR(pfTau->eta(), pfTau->phi(), (*signalPFGamma)->eta(), (*signalPFGamma)->phi());
585  if ( dR > pfTau->signalConeSize() ) photonSumPt_outsideSignalCone += (*signalPFGamma)->pt();
586  }
587  if ( photonSumPt_outsideSignalCone > maxAbsPhotonSumPt_outsideSignalCone_ || photonSumPt_outsideSignalCone > (maxRelPhotonSumPt_outsideSignalCone_*pfTau->pt()) ) {
588  failsPhotonPtSumOutsideSignalConeCut = true;
589  }
590  }
591 
592  bool fails = (applyOccupancyCut_ && failsOccupancyCut) ||
593  (applySumPtCut_ && failsSumPtCut) ||
594  (applyRelativeSumPtCut_ && failsRelativeSumPtCut) ||
595  (applyPhotonPtSumOutsideSignalConeCut_ && failsPhotonPtSumOutsideSignalConeCut);
596 
597 
598  if (pfTau->pt() > minPtForNoIso_ && minPtForNoIso_ > 0.){
599  return 1.;
600  LogDebug("discriminate") << "tau pt = " << pfTau->pt() << "\t min cutoff pt = " << minPtForNoIso_ ;
601  }
602 
603  // We did error checking in the constructor, so this is safe.
604  if ( storeRawSumPt_ ) {
605  return totalPt;
606  } else if ( storeRawPUsumPt_ ) {
607  if ( applyDeltaBeta_ ) return puPt;
608  else if ( applyRhoCorrection_ ) return rhoThisEvent_;
609  else return 0.;
610  } else if ( storeRawOccupancy_ ) {
611  return nOccupants;
612  } else if ( storeRawFootprintCorrection_ ) {
613  return footprintCorrection_value;
614  } else if ( storeRawPhotonSumPt_outsideSignalCone_ ) {
615  return photonSumPt_outsideSignalCone;
616  } else {
617  return (fails ? 0. : 1.);
618  }
619 }
620 
#define LogDebug(id)
PFRecoTauDiscriminationByIsolation(const edm::ParameterSet &pset)
T getParameter(std::string const &) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
std::auto_ptr< tau::RecoTauQualityCuts > pileupQcutsGeneralQCuts_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
selection
main part
Definition: corrVsCorr.py:98
void beginEvent(const edm::Event &evt, const edm::EventSetup &evtSetup) override
bool exists(std::string const &parameterName) const
checks if a parameter exists
double weightedSum(const std::vector< PFCandidatePtr > &inColl_, double eta, double phi) const
virtual int charge() const final
electric charge
Definition: LeafCandidate.h:91
std::vector< std::unique_ptr< FootprintCorrection > > footprintCorrections_
edm::EDGetTokenT< reco::PFCandidateCollection > pfCand_token
FootprintCorrection(const std::string &selection, const std::string &offset)
def pv(vc)
Definition: MetAnalyzer.py:6
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
std::vector< reco::PFCandidatePtr > chargedPFCandidatesInEvent_
std::pair< edm::ParameterSet, edm::ParameterSet > factorizePUQCuts(const edm::ParameterSet &inputSet)
T min(T a, T b)
Definition: MathUtil.h:58
#define LogTrace(id)
std::auto_ptr< tau::RecoTauQualityCuts > qcuts_
std::auto_ptr< tau::RecoTauQualityCuts > pileupQcutsPUTrackSelection_
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
double discriminate(const PFTauRef &pfTau) const override
ParameterSet const & getParameterSet(std::string const &) const
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
fixed size matrix
HLT enums.
std::auto_ptr< tau::RecoTauVertexAssociator > vertexAssociator_
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1
edm::EDGetTokenT< reco::VertexCollection > vertex_token