CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  applyOccupancyCut_ = pset.getParameter<bool>("applyOccupancyCut");
43  maximumOccupancy_ = pset.getParameter<uint32_t>("maximumOccupancy");
44 
45  applySumPtCut_ = pset.getParameter<bool>("applySumPtCut");
46  maximumSumPt_ = pset.getParameter<double>("maximumSumPtCut");
47 
48  applyRelativeSumPtCut_ = pset.getParameter<bool>(
49  "applyRelativeSumPtCut");
50  maximumRelativeSumPt_ = pset.getParameter<double>(
51  "relativeSumPtCut");
52  offsetRelativeSumPt_ = pset.exists("relativeSumPtOffset") ?
53  pset.getParameter<double>("relativeSumPtOffset") : 0.0;
54 
55  storeRawOccupancy_ = pset.exists("storeRawOccupancy") ?
56  pset.getParameter<bool>("storeRawOccupancy") : false;
57  storeRawSumPt_ = pset.exists("storeRawSumPt") ?
58  pset.getParameter<bool>("storeRawSumPt") : false;
59  storeRawPUsumPt_ = pset.exists("storeRawPUsumPt") ?
60  pset.getParameter<bool>("storeRawPUsumPt") : false;
61  storeRawFootprintCorrection_ = pset.exists("storeRawFootprintCorrection") ?
62  pset.getParameter<bool>("storeRawFootprintCorrection") : false;
63  storeRawPhotonSumPt_outsideSignalCone_ = pset.exists("storeRawPhotonSumPt_outsideSignalCone") ?
64  pset.getParameter<bool>("storeRawPhotonSumPt_outsideSignalCone") : false;
65 
66  // Sanity check on requested options. We can't apply cuts and store the
67  // raw output at the same time
68  if ( applySumPtCut_ || applyOccupancyCut_ || applyRelativeSumPtCut_ ) {
69  if ( storeRawSumPt_ || storeRawOccupancy_ || storeRawPUsumPt_ ) {
70  throw cms::Exception("BadIsoConfig")
71  << "A 'store raw' and a 'apply cut' option have been set to true "
72  << "simultaneously. These options are mutually exclusive.";
73  }
74  }
75 
76  // sanity check2 - can't use weighted and unweighted iso at the same time
77  if ( includeGammas_ && calculateWeights_ ) {
78  throw cms::Exception("BasIsoConfig")
79  << "Both 'ApplyDiscriminationByECALIsolation' and 'ApplyDiscriminationByWeightedECALIsolation' "
80  << "have been set to true. These options are mutually exclusive.";
81  }
82 
83  // Can only store one type
84  int numStoreOptions = 0;
85  if ( storeRawSumPt_ ) ++numStoreOptions;
86  if ( storeRawOccupancy_ ) ++numStoreOptions;
87  if ( storeRawPUsumPt_ ) ++numStoreOptions;
88  if ( storeRawFootprintCorrection_ ) ++numStoreOptions;
89  if ( storeRawPhotonSumPt_outsideSignalCone_ ) ++numStoreOptions;
90  if ( numStoreOptions > 1 ) {
91  throw cms::Exception("BadIsoConfig")
92  << "Multiple 'store sum pt' and/or 'store occupancy' options are set."
93  << " These options are mutually exclusive.";
94  }
95 
96  if ( pset.exists("customOuterCone") ) {
97  customIsoCone_ = pset.getParameter<double>("customOuterCone");
98  } else {
99  customIsoCone_ = -1;
100  }
101 
102  applyPhotonPtSumOutsideSignalConeCut_ = ( pset.exists("applyPhotonPtSumOutsideSignalConeCut") ) ?
103  pset.getParameter<bool>("applyPhotonPtSumOutsideSignalConeCut") : false;
104  if ( applyPhotonPtSumOutsideSignalConeCut_ ) {
105  maxAbsPhotonSumPt_outsideSignalCone_ = pset.getParameter<double>("maxAbsPhotonSumPt_outsideSignalCone");
106  maxRelPhotonSumPt_outsideSignalCone_ = pset.getParameter<double>("maxRelPhotonSumPt_outsideSignalCone");
107  }
108 
109  applyFootprintCorrection_ = ( pset.exists("applyFootprintCorrection") ) ?
110  pset.getParameter<bool>("applyFootprintCorrection") : false;
111  if ( applyFootprintCorrection_ || storeRawFootprintCorrection_ ) {
112  edm::VParameterSet cfgFootprintCorrections = pset.getParameter<edm::VParameterSet>("footprintCorrections");
113  for ( edm::VParameterSet::const_iterator cfgFootprintCorrection = cfgFootprintCorrections.begin();
114  cfgFootprintCorrection != cfgFootprintCorrections.end(); ++cfgFootprintCorrection ) {
115  std::string selection = cfgFootprintCorrection->getParameter<std::string>("selection");
116  std::string offset = cfgFootprintCorrection->getParameter<std::string>("offset");
117  std::unique_ptr<FootprintCorrection> footprintCorrection(new FootprintCorrection(selection, offset));
118  footprintCorrections_.push_back(std::move(footprintCorrection));
119  }
120  }
121 
122  // Get the quality cuts specific to the isolation region
123  edm::ParameterSet isolationQCuts = qualityCutsPSet_.getParameterSet(
124  "isolationQualityCuts");
125 
126  qcuts_.reset(new tau::RecoTauQualityCuts(isolationQCuts));
127 
128  vertexAssociator_.reset(
129  new tau::RecoTauVertexAssociator(qualityCutsPSet_,consumesCollector()));
130 
131  applyDeltaBeta_ = pset.exists("applyDeltaBetaCorrection") ?
132  pset.getParameter<bool>("applyDeltaBetaCorrection") : false;
133 
134  if ( applyDeltaBeta_ || calculateWeights_ ) {
135  // Factorize the isolation QCuts into those that are used to
136  // select PU and those that are not.
137  std::pair<edm::ParameterSet, edm::ParameterSet> puFactorizedIsoQCuts =
138  reco::tau::factorizePUQCuts(isolationQCuts);
139 
140  // Determine the pt threshold for the PU tracks
141  // First check if the user specifies explicitly the cut.
142  if ( pset.exists("deltaBetaPUTrackPtCutOverride") ) {
143  puFactorizedIsoQCuts.second.addParameter<double>(
144  "minTrackPt",
145  pset.getParameter<double>("deltaBetaPUTrackPtCutOverride"));
146  } else {
147  // Secondly take it from the minGammaEt
148  puFactorizedIsoQCuts.second.addParameter<double>(
149  "minTrackPt",
150  isolationQCuts.getParameter<double>("minGammaEt"));
151  }
152 
153  pileupQcutsPUTrackSelection_.reset(new tau::RecoTauQualityCuts(
154  puFactorizedIsoQCuts.first));
155 
156  pileupQcutsGeneralQCuts_.reset(new tau::RecoTauQualityCuts(
157  puFactorizedIsoQCuts.second));
158 
159  pfCandSrc_ = pset.getParameter<edm::InputTag>("particleFlowSrc");
160  pfCand_token = consumes<reco::PFCandidateCollection>(pfCandSrc_);
161  vertexSrc_ = pset.getParameter<edm::InputTag>("vertexSrc");
162  vertex_token = consumes<reco::VertexCollection>(vertexSrc_);
163  deltaBetaCollectionCone_ = pset.getParameter<double>(
164  "isoConeSizeForDeltaBeta");
165  std::string deltaBetaFactorFormula =
166  pset.getParameter<string>("deltaBetaFactor");
167  deltaBetaFormula_.reset(
168  new TFormula("DB_corr", deltaBetaFactorFormula.c_str()));
169  }
170 
171  applyRhoCorrection_ = pset.exists("applyRhoCorrection") ?
172  pset.getParameter<bool>("applyRhoCorrection") : false;
173  if ( applyRhoCorrection_ ) {
174  rhoProducer_ = pset.getParameter<edm::InputTag>("rhoProducer");
175  rho_token=consumes<double>(rhoProducer_);
176  rhoConeSize_ = pset.getParameter<double>("rhoConeSize");
177  rhoUEOffsetCorrection_ =
178  pset.getParameter<double>("rhoUEOffsetCorrection");
179  }
180  useAllPFCands_ = pset.exists("UseAllPFCandsForWeights") ?
181  pset.getParameter<bool>("UseAllPFCandsForWeights") : false;
182 
183  verbosity_ = ( pset.exists("verbosity") ) ?
184  pset.getParameter<int>("verbosity") : 0;
185  }
186 
188  {
189  }
190 
191  void beginEvent(const edm::Event& evt, const edm::EventSetup& evtSetup) override;
192  double discriminate(const PFTauRef& pfTau) const override;
193 
194  inline double weightedSum(std::vector<PFCandidatePtr> inColl_, double eta, double phi) const {
195  double out = 1.0;
196  for (auto const & inObj_ : inColl_){
197  double sum = (inObj_->pt()*inObj_->pt())/(deltaR2(eta,phi,inObj_->eta(),inObj_->phi()));
198  if(sum > 1.0) out *= sum;
199  }
200  return out;
201  }
202 
203  private:
205 
207  std::auto_ptr<tau::RecoTauQualityCuts> qcuts_;
208 
209  // Inverted QCut which selects tracks with bad DZ/trackWeight
210  std::auto_ptr<tau::RecoTauQualityCuts> pileupQcutsPUTrackSelection_;
211  std::auto_ptr<tau::RecoTauQualityCuts> pileupQcutsGeneralQCuts_;
212 
213  std::auto_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;
214 
226 
230 
233  {
235  : selection_(selection),
236  offset_(offset)
237  {}
241  };
242  std::vector<std::unique_ptr<FootprintCorrection> > footprintCorrections_;
243 
244  // Options to store the raw value in the discriminator instead of boolean pass/fail flag
250 
251  /* **********************************************************************
252  **** Pileup Subtraction Parameters ***********************************
253  **********************************************************************/
254 
255  // Delta Beta correction
259  // Keep track of how many vertices are in the event
262  std::vector<reco::PFCandidatePtr> chargedPFCandidatesInEvent_;
263  // Size of cone used to collect PU tracks
265  std::auto_ptr<TFormula> deltaBetaFormula_;
267 
268  // Rho correction
273  double rhoConeSize_;
277 
278  // Flag to enable/disable debug output
280 };
281 
283 {
284  // NB: The use of the PV in this context is necessitated by its use in
285  // applying quality cuts to the different objects in the isolation cone
286  // The vertex associator contains the logic to select the appropriate vertex
287  // We need to pass it the event so it can load the vertices.
288  vertexAssociator_->setEvent(event);
289 
290  // If we are applying the delta beta correction, we need to get the PF
291  // candidates from the event so we can find the PU tracks.
292  if ( applyDeltaBeta_ || calculateWeights_ ) {
293  // Collect all the PF pile up tracks
295  event.getByToken(pfCand_token, pfCandidates);
296  chargedPFCandidatesInEvent_.clear();
297  chargedPFCandidatesInEvent_.reserve(pfCandidates->size());
298  size_t numPFCandidates = pfCandidates->size();
299  for ( size_t i = 0; i < numPFCandidates; ++i ) {
300  reco::PFCandidatePtr pfCandidate(pfCandidates, i);
301  if ( pfCandidate->charge() != 0 ) {
302  chargedPFCandidatesInEvent_.push_back(pfCandidate);
303  }
304  }
305  // Count all the vertices in the event, to parameterize the DB
306  // correction factor
308  event.getByToken(vertex_token, vertices);
309  size_t nVtxThisEvent = vertices->size();
310  deltaBetaFactorThisEvent_ = deltaBetaFormula_->Eval(nVtxThisEvent);
311  }
312 
313  if ( applyRhoCorrection_ ) {
314  edm::Handle<double> rhoHandle_;
315  event.getByToken(rho_token, rhoHandle_);
316  rhoThisEvent_ = (*rhoHandle_ - rhoUEOffsetCorrection_)*
317  (3.14159)*rhoConeSize_*rhoConeSize_;
318  }
319 }
320 
321 double
323 {
324  LogDebug("discriminate") << " tau: Pt = " << pfTau->pt() << ", eta = " << pfTau->eta() << ", phi = " << pfTau->phi();
325  LogDebug("discriminate") << *pfTau ;
326 
327  // collect the objects we are working with (ie tracks, tracks+gammas, etc)
328  std::vector<PFCandidatePtr> isoCharged_;
329  std::vector<PFCandidatePtr> isoNeutral_;
330  std::vector<PFCandidatePtr> isoPU_;
331  PFCandidateCollection isoNeutralWeight_;
332  std::vector<PFCandidatePtr> chPV_;
333  isoCharged_.reserve(pfTau->isolationPFChargedHadrCands().size());
334  isoNeutral_.reserve(pfTau->isolationPFGammaCands().size());
335  isoPU_.reserve(chargedPFCandidatesInEvent_.size());
336  isoNeutralWeight_.reserve(pfTau->isolationPFGammaCands().size());
337 
338  chPV_.reserve(chargedPFCandidatesInEvent_.size());
339 
340  // Get the primary vertex associated to this tau
341  reco::VertexRef pv = vertexAssociator_->associatedVertex(*pfTau);
342  // Let the quality cuts know which the vertex to use when applying selections
343  // on dz, etc.
344  if ( verbosity_ ) {
345  if ( pv.isNonnull() ) {
346  LogTrace("discriminate") << "pv: x = " << pv->position().x() << ", y = " << pv->position().y() << ", z = " << pv->position().z() ;
347  } else {
348  LogTrace("discriminate") << "pv: N/A" ;
349  }
350  if ( pfTau->leadPFChargedHadrCand().isNonnull() ) {
351  LogTrace("discriminate") << "leadPFChargedHadron:"
352  << " Pt = " << pfTau->leadPFChargedHadrCand()->pt() << ","
353  << " eta = " << pfTau->leadPFChargedHadrCand()->eta() << ","
354  << " phi = " << pfTau->leadPFChargedHadrCand()->phi() ;
355  } else {
356  LogTrace("discriminate") << "leadPFChargedHadron: N/A" ;
357  }
358  }
359 
360  // CV: isolation is not well defined in case primary vertex or leading charged hadron do not exist
361  if ( !(pv.isNonnull() && pfTau->leadPFChargedHadrCand().isNonnull()) ) return 0.;
362 
363  qcuts_->setPV(pv);
364  qcuts_->setLeadTrack(pfTau->leadPFChargedHadrCand());
365 
366  if ( applyDeltaBeta_ || calculateWeights_) {
367  pileupQcutsGeneralQCuts_->setPV(pv);
368  pileupQcutsGeneralQCuts_->setLeadTrack(pfTau->leadPFChargedHadrCand());
369  pileupQcutsPUTrackSelection_->setPV(pv);
370  pileupQcutsPUTrackSelection_->setLeadTrack(pfTau->leadPFChargedHadrCand());
371  }
372 
373  // Load the tracks if they are being used.
374  if ( includeTracks_ ) {
375  for( auto const & cand : pfTau->isolationPFChargedHadrCands() ) {
376  if ( qcuts_->filterCandRef(cand) ) {
377  LogTrace("discriminate") << "adding charged iso cand with pt " << cand->pt() ;
378  isoCharged_.push_back(cand);
379  }
380  }
381  }
382  if ( includeGammas_ || calculateWeights_ ) {
383  for( auto const & cand : pfTau->isolationPFGammaCands() ) {
384  if ( qcuts_->filterCandRef(cand) ) {
385  LogTrace("discriminate") << "adding neutral iso cand with pt " << cand->pt() ;
386  isoNeutral_.push_back(cand);
387  }
388  }
389  }
390 
393 
394  // If desired, get PU tracks.
395  if ( applyDeltaBeta_ || calculateWeights_) {
396  // First select by inverted the DZ/track weight cuts. True = invert
397  if ( verbosity_ ) {
398  std::cout << "Initial PFCands: " << chargedPFCandidatesInEvent_.size() << std::endl;
399  }
400 
401  std::vector<PFCandidatePtr> allPU =
402  pileupQcutsPUTrackSelection_->filterCandRefs(
403  chargedPFCandidatesInEvent_, true);
404 
405  std::vector<PFCandidatePtr> allNPU =
406  pileupQcutsPUTrackSelection_->filterCandRefs(
407  chargedPFCandidatesInEvent_);
408  LogTrace("discriminate") << "After track cuts: " << allPU.size() ;
409 
410  // Now apply the rest of the cuts, like pt, and TIP, tracker hits, etc
411  if ( !useAllPFCands_ ) {
412  std::vector<PFCandidatePtr> cleanPU =
413  pileupQcutsGeneralQCuts_->filterCandRefs(allPU);
414 
415  std::vector<PFCandidatePtr> cleanNPU =
416  pileupQcutsGeneralQCuts_->filterCandRefs(allNPU);
417 
418  LogTrace("discriminate") << "After cleaning cuts: " << cleanPU.size() ;
419 
420  // Only select PU tracks inside the isolation cone.
421  DRFilter deltaBetaFilter(pfTau->p4(), 0, deltaBetaCollectionCone_);
422  for ( auto const & cand : cleanPU ) {
423  if ( deltaBetaFilter(cand) ) isoPU_.push_back(cand);
424  }
425 
426  for ( auto const & cand : cleanNPU ) {
427  if ( deltaBetaFilter(cand) ) chPV_.push_back(cand);
428  }
429  LogTrace("discriminate") << "After cone cuts: " << isoPU_.size() << " " << chPV_.size() ;
430  } else {
431  isoPU_ = allPU;
432  chPV_ = allNPU;
433  }
434  }
435 
436  if ( calculateWeights_ ) {
437  for ( auto const & isoObject : isoNeutral_ ) {
438  if ( isoObject->charge() != 0 ) {
439  // weight only neutral objects
440  isoNeutralWeight_.push_back(*isoObject);
441  continue;
442  }
443 
444  double eta = isoObject->eta();
445  double phi = isoObject->phi();
446  double sumNPU = 0.5*log(weightedSum(chPV_, eta, phi));
447 
448  double sumPU = 0.5*log(weightedSum(isoPU_, eta, phi));
449  PFCandidate neutral = (*isoObject);
450  if ( (sumNPU + sumPU) > 0 ) neutral.setP4(((sumNPU)/(sumNPU + sumPU))*neutral.p4());
451 
452  isoNeutralWeight_.push_back(neutral);
453  }
454  }
455 
456  // Check if we want a custom iso cone
457  if ( customIsoCone_ >= 0. ) {
458  DRFilter filter(pfTau->p4(), 0, customIsoCone_);
459  DRFilter2 filter2(pfTau->p4(), 0, customIsoCone_);
460  std::vector<PFCandidatePtr> isoCharged_filter;
461  std::vector<PFCandidatePtr> isoNeutral_filter;
462  PFCandidateCollection isoNeutralWeight_filter;
463  // Remove all the objects not in our iso cone
464  for( auto const & isoObject : isoCharged_ ) {
465  if ( filter(isoObject) ) isoCharged_filter.push_back(isoObject);
466  }
467  if(!calculateWeights_){
468  for( auto const & isoObject : isoNeutral_ ) {
469  if ( filter(isoObject) ) isoNeutral_filter.push_back(isoObject);
470  }
471  isoNeutral_ = isoNeutral_filter;
472  }else{
473  for( auto const & isoObject : isoNeutralWeight_){
474  if ( filter2(isoObject) ) isoNeutralWeight_filter.push_back(isoObject);
475  }
476  isoNeutralWeight_ = isoNeutralWeight_filter;
477  }
478  isoCharged_ = isoCharged_filter;
479  }
480 
481  bool failsOccupancyCut = false;
482  bool failsSumPtCut = false;
483  bool failsRelativeSumPtCut = false;
484 
485 //--- nObjects requirement
486  int neutrals = isoNeutral_.size();
487 
488  if ( applyDeltaBeta_ ) {
489  neutrals -= TMath::Nint(deltaBetaFactorThisEvent_*isoPU_.size());
490  }
491  if ( neutrals < 0 ) {
492  neutrals = 0;
493  }
494 
495  size_t nOccupants = isoCharged_.size() + neutrals;
496 
497  failsOccupancyCut = ( nOccupants > maximumOccupancy_ );
498 
499  double footprintCorrection_value = 0.;
500  if ( applyFootprintCorrection_ || storeRawFootprintCorrection_ ) {
501  for ( std::vector<std::unique_ptr<FootprintCorrection> >::const_iterator footprintCorrection = footprintCorrections_.begin();
502  footprintCorrection != footprintCorrections_.end(); ++footprintCorrection ) {
503  if ( (*footprintCorrection)->selection_(*pfTau) ) {
504  footprintCorrection_value = (*footprintCorrection)->offset_(*pfTau);
505  }
506  }
507  }
508 
509  double totalPt = 0.;
510  double puPt = 0.;
511 //--- Sum PT requirement
512  if ( applySumPtCut_ || applyRelativeSumPtCut_ || storeRawSumPt_ || storeRawPUsumPt_ ) {
513  double chargedPt = 0.;
514  double neutralPt = 0.;
515  double weightedNeutralPt = 0.;
516  for ( auto const & isoObject : isoCharged_ ) {
517  chargedPt += isoObject->pt();
518  }
519  if ( !calculateWeights_ ) {
520  for ( auto const & isoObject : isoNeutral_ ) {
521  neutralPt += isoObject->pt();
522  }
523  } else {
524  for ( auto const & isoObject : isoNeutralWeight_ ) {
525  weightedNeutralPt += isoObject.pt();
526  }
527  }
528  for ( auto const & isoObject : isoPU_ ) {
529  puPt += isoObject->pt();
530  }
531  LogTrace("discriminate") << "chargedPt = " << chargedPt ;
532  LogTrace("discriminate") << "neutralPt = " << neutralPt ;
533  LogTrace("discriminate") << "weighted neutral Pt = " << weightedNeutralPt ;
534  LogTrace("discriminate") << "puPt = " << puPt << " (delta-beta corr. = " << (deltaBetaFactorThisEvent_*puPt) << ")" ;
535 
536  if ( calculateWeights_ ) {
537  neutralPt = weightedNeutralPt;
538  }
539 
540  if ( applyDeltaBeta_ ) {
541  neutralPt -= (deltaBetaFactorThisEvent_*puPt);
542  }
543 
544  if ( applyFootprintCorrection_ ) {
545  neutralPt -= footprintCorrection_value;
546  }
547 
548  if ( applyRhoCorrection_ ) {
549  neutralPt -= rhoThisEvent_;
550  }
551 
552  if ( neutralPt < 0. ) {
553  neutralPt = 0.;
554  }
555 
556  totalPt = chargedPt + neutralPt;
557  LogTrace("discriminate") << "totalPt = " << totalPt << " (cut = " << maximumSumPt_ << ")" ;
558 
559  failsSumPtCut = (totalPt > maximumSumPt_);
560 
561 //--- Relative Sum PT requirement
562  failsRelativeSumPtCut = (totalPt > ((pfTau->pt() - offsetRelativeSumPt_)*maximumRelativeSumPt_));
563  }
564 
565  bool failsPhotonPtSumOutsideSignalConeCut = false;
566  double photonSumPt_outsideSignalCone = 0.;
567  if ( applyPhotonPtSumOutsideSignalConeCut_ || storeRawPhotonSumPt_outsideSignalCone_ ) {
568  const std::vector<reco::PFCandidatePtr>& signalPFGammas = pfTau->signalPFGammaCands();
569  for ( std::vector<reco::PFCandidatePtr>::const_iterator signalPFGamma = signalPFGammas.begin();
570  signalPFGamma != signalPFGammas.end(); ++signalPFGamma ) {
571  double dR = deltaR(pfTau->eta(), pfTau->phi(), (*signalPFGamma)->eta(), (*signalPFGamma)->phi());
572  if ( dR > pfTau->signalConeSize() ) photonSumPt_outsideSignalCone += (*signalPFGamma)->pt();
573  }
574  if ( photonSumPt_outsideSignalCone > maxAbsPhotonSumPt_outsideSignalCone_ || photonSumPt_outsideSignalCone > (maxRelPhotonSumPt_outsideSignalCone_*pfTau->pt()) ) {
575  failsPhotonPtSumOutsideSignalConeCut = true;
576  }
577  }
578 
579  bool fails = (applyOccupancyCut_ && failsOccupancyCut) ||
580  (applySumPtCut_ && failsSumPtCut) ||
581  (applyRelativeSumPtCut_ && failsRelativeSumPtCut) ||
582  (applyPhotonPtSumOutsideSignalConeCut_ && failsPhotonPtSumOutsideSignalConeCut);
583 
584  // We did error checking in the constructor, so this is safe.
585  if ( storeRawSumPt_ ) {
586  return totalPt;
587  } else if ( storeRawPUsumPt_ ) {
588  if ( applyDeltaBeta_ ) return puPt;
589  else if ( applyRhoCorrection_ ) return rhoThisEvent_;
590  else return 0.;
591  } else if ( storeRawOccupancy_ ) {
592  return nOccupants;
593  } else if ( storeRawFootprintCorrection_ ) {
594  return footprintCorrection_value;
595  } else if ( storeRawPhotonSumPt_outsideSignalCone_ ) {
596  return photonSumPt_outsideSignalCone;
597  } else {
598  return (fails ? 0. : 1.);
599  }
600 }
601 
#define LogDebug(id)
PFRecoTauDiscriminationByIsolation(const edm::ParameterSet &pset)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
double weightedSum(std::vector< PFCandidatePtr > inColl_, double eta, double phi) 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
virtual void setP4(const LorentzVector &p4)
set 4-momentum
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
std::vector< std::unique_ptr< FootprintCorrection > > footprintCorrections_
edm::EDGetTokenT< reco::PFCandidateCollection > pfCand_token
FootprintCorrection(const std::string &selection, const std::string &offset)
def move
Definition: eostools.py:510
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
std::vector< reco::PFCandidatePtr > chargedPFCandidatesInEvent_
std::pair< edm::ParameterSet, edm::ParameterSet > factorizePUQCuts(const edm::ParameterSet &inputSet)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
double discriminate(const PFTauRef &pfTau) const override
#define LogTrace(id)
std::auto_ptr< tau::RecoTauQualityCuts > qcuts_
std::auto_ptr< tau::RecoTauQualityCuts > pileupQcutsPUTrackSelection_
tuple out
Definition: dbtoconf.py:99
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
ParameterSet const & getParameterSet(std::string const &) const
Geom::Phi< T > phi() const
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
std::auto_ptr< tau::RecoTauVertexAssociator > vertexAssociator_
tuple cout
Definition: gather_cfg.py:121
moduleLabel_(iConfig.getParameter< string >("@module_label"))
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
edm::EDGetTokenT< reco::VertexCollection > vertex_token
tuple log
Definition: cmsBatch.py:341