CMS 3D CMS Logo

RecoTauBuilderConePlugin.cc
Go to the documentation of this file.
1 /*
2  * RecoTauBuilderConePlugin
3  *
4  * Build a PFTau using cones defined in DeltaR.
5  *
6  * Original Authors: Ludovic Houchu, Simone Gennai
7  * Modifications: Evan K. Friis
8  *
9  */
10 
11 #include <vector>
12 #include <algorithm>
13 
18 
21 
27 
29 
30 namespace reco { namespace tau {
31 
32 typedef std::vector<RecoTauPiZero> PiZeroList;
33 
35  public:
38  // Build a tau from a jet
40  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons,
41  const std::vector<RecoTauPiZero>& piZeros,
42  const std::vector<PFCandidatePtr>& regionalExtras) const override;
43  private:
45 
47 
49 
50  /* String function to extract values from PFJets */
52 
53  // Cone defintions
54  JetFunc matchingCone_;
59  JetFunc isoConePiZeros_;
62 
66 
67  void setTauQuantities(reco::PFTau& aTau,
71  double minRelPhotonSumPt_outsideSignalCone = 1.e+9) const;
72 
73 };
74 
75 // ctor - initialize all of our variables
79  "qualityCuts").getParameterSet("signalQualityCuts")),
80  usePFLeptonsAsChargedHadrons_(pset.getParameter<bool>("usePFLeptons")),
81  leadObjecPtThreshold_(pset.getParameter<double>("leadObjectPt")),
82  matchingCone_(pset.getParameter<std::string>("matchingCone")),
83  signalConeChargedHadrons_(pset.getParameter<std::string>(
84  "signalConeChargedHadrons")),
86  pset.getParameter<std::string>("isoConeChargedHadrons")),
88  pset.getParameter<std::string>("signalConePiZeros")),
89  signalConeSizeToStore_(//MB: stored inside PFTau and used to compte photon-pt-out-of-signal-cone and DM hypothsis
90  pset.getParameter<std::string>("signalConePiZeros")),
92  pset.getParameter<std::string>("isoConePiZeros")),
94  pset.getParameter<std::string>("signalConeNeutralHadrons")),
96  pset.getParameter<std::string>("isoConeNeutralHadrons")),
98  pset.getParameter<int>("maxSignalConeChargedHadrons")),
100  pset.getParameter<double>("minAbsPhotonSumPt_insideSignalCone")),
102  pset.getParameter<double>("minRelPhotonSumPt_insideSignalCone"))
103 
104 {}
105 
106 namespace xclean
107 {
108  // define template specialization for cross-cleaning
109  template<>
111  {
112  // Get the list of objects we need to clean
113  for ( cone::PFCandPtrDRFilterIter signalTrack = signalTracksBegin; signalTrack != signalTracksEnd; ++signalTrack ) {
114  toRemove_.insert(reco::CandidatePtr(*signalTrack));
115  }
116  }
117 
118  template<>
119  inline void CrossCleanPtrs<PiZeroList::const_iterator>::initialize(const PiZeroList::const_iterator& piZerosBegin, const PiZeroList::const_iterator& piZerosEnd)
120  {
121  BOOST_FOREACH( const PFCandidatePtr &ptr, flattenPiZeros(piZerosBegin, piZerosEnd) ) {
122  toRemove_.insert(CandidatePtr(ptr));
123  }
124  }
125 }
126 
132  //MB: Set tau quantities which are computed by RecoTauConstructor::get()
133  // method with PFRecoTauChargedHadrons not used here
134 
135  // Set charge
136  int charge = 0;
137  int leadCharge = aTau.leadPFChargedHadrCand().isNonnull() ? aTau.leadPFChargedHadrCand()->charge() : 0;
138  const std::vector<reco::PFCandidatePtr>& pfChs = aTau.signalPFChargedHadrCands();
139  for(const reco::PFCandidatePtr& pfCh : pfChs){
140  charge += pfCh->charge();
141  }
142  charge = charge==0 ? leadCharge : charge;
143  aTau.setCharge(charge);
144 
145  // Set PDG id
146  aTau.setPdgId(aTau.charge() < 0 ? 15 : -15);
147 
148  // Set Decay Mode
150  unsigned int nPiZeros = 0;
151  const std::vector<RecoTauPiZero>& piZeros = aTau.signalPiZeroCandidates();
152  for(const RecoTauPiZero& piZero : piZeros) {
153  double photonSumPt_insideSignalCone = 0.;
154  double photonSumPt_outsideSignalCone = 0.;
155  int numPhotons = piZero.numberOfDaughters();
156  for(int idxPhoton = 0; idxPhoton < numPhotons; ++idxPhoton) {
157  const reco::Candidate* photon = piZero.daughter(idxPhoton);
158  double dR = deltaR(photon->p4(), aTau.p4());
159  if(dR < aTau.signalConeSize() ){
160  photonSumPt_insideSignalCone += photon->pt();
161  } else {
162  photonSumPt_outsideSignalCone += photon->pt();
163  }
164  }
165  if( photonSumPt_insideSignalCone > minAbsPhotonSumPt_insideSignalCone || photonSumPt_insideSignalCone > (minRelPhotonSumPt_insideSignalCone*aTau.pt()) ||
166  photonSumPt_outsideSignalCone > minAbsPhotonSumPt_outsideSignalCone || photonSumPt_outsideSignalCone > (minRelPhotonSumPt_outsideSignalCone*aTau.pt()) ) ++nPiZeros;
167  }
168  // Find the maximum number of PiZeros our parameterization can hold
169  const unsigned int maxPiZeros = PFTau::kOneProngNPiZero;
170  // Determine our track index
171  unsigned int nCharged = pfChs.size();
172  if(nCharged>0){
173  unsigned int trackIndex = (nCharged - 1)*(maxPiZeros + 1);
174  // Check if we handle the given number of tracks
175  if(trackIndex >= PFTau::kRareDecayMode)
176  dm = PFTau::kRareDecayMode;
177  else
178  dm = static_cast<PFTau::hadronicDecayMode>(trackIndex + std::min(nPiZeros,maxPiZeros) );
179  }
180  else{
181  dm = PFTau::kNull;
182  }
183  aTau.setDecayMode(dm);
184  return;
185 }
186 
188  const reco::PFJetRef& jet,
189  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons,
190  const std::vector<RecoTauPiZero>& piZeros,
191  const std::vector<PFCandidatePtr>& regionalExtras) const {
192  //std::cout << "<RecoTauBuilderConePlugin::operator()>:" << std::endl;
193  //std::cout << " jet: Pt = " << jet->pt() << ", eta = " << jet->eta() << ", phi = " << jet->phi() << std::endl;
194 
195  // Get access to our cone tools
196  using namespace cone;
197  // Define output. We only produce one tau per jet for the cone algo.
199 
200  // Our tau builder - the true indicates to automatically copy gamma candidates
201  // from the pizeros.
202  RecoTauConstructor tau(jet, getPFCands(), true,
206  // Setup our quality cuts to use the current vertex (supplied by base class)
207  qcuts_.setPV(primaryVertex(jet));
208 
209  typedef std::vector<PFCandidatePtr> PFCandPtrs;
210 
211  // Get the PF Charged hadrons + quality cuts
212  PFCandPtrs pfchs;
215  } else {
216  // Check if we want to include electrons in muons in "charged hadron"
217  // collection. This is the preferred behavior, as the PF lepton selections
218  // are very loose.
219  pfchs = qcuts_.filterCandRefs(pfChargedCands(*jet));
220  }
221 
222  // CV: sort collection of PF Charged hadrons by descending Pt
223  std::sort(pfchs.begin(), pfchs.end(), SortPFCandsDescendingPt());
224 
225  // Get the PF gammas
226  PFCandPtrs pfGammas = qcuts_.filterCandRefs(
228  // Neutral hadrons
229  PFCandPtrs pfnhs = qcuts_.filterCandRefs(
231 
232  // All the extra junk
233  PFCandPtrs regionalJunk = qcuts_.filterCandRefs(regionalExtras);
234 
235  /***********************************************
236  ****** Lead Candidate Finding **********
237  ***********************************************/
238 
239  // Define our matching cone and filters
240  double matchingCone = matchingCone_(*jet);
241  PFCandPtrDRFilter matchingConeFilter(jet->p4(), 0, matchingCone);
242 
243  // Find the maximum PFCharged hadron in the matching cone. The call to
244  // PFCandidates always a sorted list, so we can just take the first if it
245  // if it exists.
246  PFCandidatePtr leadPFCH;
247  PFCandPtrs::iterator leadPFCH_iter =
248  std::find_if(pfchs.begin(), pfchs.end(), matchingConeFilter);
249 
250  if (leadPFCH_iter != pfchs.end()) {
251  leadPFCH = *leadPFCH_iter;
252  // Set leading candidate
253  tau.setleadPFChargedHadrCand(leadPFCH);
254  } else {
255  // If there is no leading charged candidate at all, return nothing - the
256  // producer class that owns the plugin will build a null tau if desired.
257  return output.release();
258  }
259 
260  // Find the leading neutral candidate
261  PFCandidatePtr leadPFGamma;
262  PFCandPtrs::iterator leadPFGamma_iter =
263  std::find_if(pfGammas.begin(), pfGammas.end(), matchingConeFilter);
264 
265  if (leadPFGamma_iter != pfGammas.end()) {
266  leadPFGamma = *leadPFGamma_iter;
267  // Set leading neutral candidate
268  tau.setleadPFNeutralCand(leadPFGamma);
269  }
270 
272  // Always use the leadPFCH if it is above our threshold
273  if (leadPFCH.isNonnull() && leadPFCH->pt() > leadObjecPtThreshold_) {
274  leadPFCand = leadPFCH;
275  } else if (leadPFGamma.isNonnull() &&
276  leadPFGamma->pt() > leadObjecPtThreshold_) {
277  // Otherwise use the lead Gamma if it is above threshold
278  leadPFCand = leadPFGamma;
279  } else {
280  // If both are too low PT, just take the charged one
281  leadPFCand = leadPFCH;
282  }
283 
284  tau.setleadPFCand(leadPFCand);
285 
286  // Our cone axis is defined about the lead charged hadron
287  reco::Candidate::LorentzVector coneAxis = leadPFCH->p4();
288 
289  /***********************************************
290  ****** Cone Construction **********
291  ***********************************************/
292 
293  // Define the signal and isolation cone sizes for this jet and build filters
294  // to select elements in the given DeltaR regions
295 
296  PFCandPtrDRFilter signalConePFCHFilter(
297  coneAxis, -0.1, signalConeChargedHadrons_(*jet));
298  PFCandPtrDRFilter signalConePFNHFilter(
299  coneAxis, -0.1, signalConeNeutralHadrons_(*jet));
300  PiZeroDRFilter signalConePiZeroFilter(
301  coneAxis, -0.1, signalConePiZeros_(*jet));
302 
303  PFCandPtrDRFilter isoConePFCHFilter(
304  coneAxis, signalConeChargedHadrons_(*jet), isoConeChargedHadrons_(*jet));
305  PFCandPtrDRFilter isoConePFGammaFilter(
306  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
307  PFCandPtrDRFilter isoConePFNHFilter(
308  coneAxis, signalConeNeutralHadrons_(*jet), isoConeNeutralHadrons_(*jet));
309  PiZeroDRFilter isoConePiZeroFilter(
310  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
311 
312  // Additionally make predicates to select the different PF object types
313  // of the regional junk objects to add to the iso cone.
315 
319 
320  // Predicate to select the regional junk in the iso cone by PF id
321  RegionalJunkConeAndIdFilter pfChargedJunk(
322  pfchCandSelector, // select charged stuff from junk
323  isoConePFCHFilter // only take those in iso cone
324  );
325 
326  RegionalJunkConeAndIdFilter pfGammaJunk(
327  pfgammaCandSelector, // select gammas from junk
328  isoConePFGammaFilter // only take those in iso cone
329  );
330 
331  RegionalJunkConeAndIdFilter pfNeutralJunk(
332  pfnhCandSelector, // select neutral stuff from junk
333  isoConePFNHFilter // select stuff in iso cone
334  );
335 
336  // Build filter iterators select the signal charged stuff.
337  PFCandPtrDRFilterIter signalPFCHCands_begin(
338  signalConePFCHFilter, pfchs.begin(), pfchs.end());
339  PFCandPtrDRFilterIter signalPFCHCands_end(
340  signalConePFCHFilter, pfchs.end(), pfchs.end());
341  PFCandPtrs signalPFCHs;
342  int numSignalPFCHs = 0;
343  PFCandPtrs isolationPFCHs;
344  int numIsolationPFCHs = 0;
345  for ( PFCandPtrDRFilterIter iter = signalPFCHCands_begin; iter != signalPFCHCands_end; ++iter ) {
346  if ( numSignalPFCHs < maxSignalConeChargedHadrons_ || maxSignalConeChargedHadrons_ == -1 ) {
347  //std::cout << "adding signalPFCH #" << numSignalPFCHs << ": Pt = " << (*iter)->pt() << ", eta = " << (*iter)->eta() << ", phi = " << (*iter)->phi() << std::endl;
348  signalPFCHs.push_back(*iter);
349  ++numSignalPFCHs;
350  } else {
351  //std::cout << "maxSignalConeChargedHadrons reached"
352  // << " --> adding isolationPFCH #" << numIsolationPFCHs << ": Pt = " << (*iter)->pt() << ", eta = " << (*iter)->eta() << ", phi = " << (*iter)->phi() << std::endl;
353  isolationPFCHs.push_back(*iter);
354  ++numIsolationPFCHs;
355  }
356  }
357  PFCandPtrs::const_iterator signalPFCHs_begin = signalPFCHs.begin();
358  PFCandPtrs::const_iterator signalPFCHs_end = signalPFCHs.end();
359 
360  // Cross clean pi zero content using signal cone charged hadron constituents.
362  signalPFCHCands_begin, signalPFCHCands_end);
363  std::vector<reco::RecoTauPiZero> cleanPiZeros = piZeroXCleaner(piZeros);
364 
365  // For the rest of the constituents, we need to filter any constituents that
366  // are already contained in the pizeros (i.e. electrons)
367  xclean::CrossCleanPtrs<PiZeroList::const_iterator> pfCandXCleaner(cleanPiZeros.begin(), cleanPiZeros.end());
368 
369  auto isolationPFCHCands_begin(
370  boost::make_filter_iterator(
371  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
372  pfchs.begin(), pfchs.end()));
373  auto isolationPFCHCands_end(
374  boost::make_filter_iterator(
375  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
376  pfchs.end(), pfchs.end()));
377  for ( auto iter = isolationPFCHCands_begin; iter != isolationPFCHCands_end; ++iter ) {
378  //std::cout << "adding isolationPFCH #" << numIsolationPFCHs << ": Pt = " << (*iter)->pt() << ", eta = " << (*iter)->eta() << ", phi = " << (*iter)->phi() << std::endl;
379  isolationPFCHs.push_back(*iter);
380  ++numIsolationPFCHs;
381  }
382  PFCandPtrs::const_iterator isolationPFCHs_begin = isolationPFCHs.begin();
383  PFCandPtrs::const_iterator isolationPFCHs_end = isolationPFCHs.end();
384 
385  // Build signal charged hadrons
388  signalPFCHs_begin, signalPFCHs_end);
389 
392  boost::make_filter_iterator(
393  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
394  pfnhs.begin(), pfnhs.end()),
395  boost::make_filter_iterator(
396  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
397  pfnhs.end(), pfnhs.end()));
398 
399  // Build signal PiZeros
401  PiZeroDRFilterIter(signalConePiZeroFilter,
402  cleanPiZeros.begin(), cleanPiZeros.end()),
403  PiZeroDRFilterIter(signalConePiZeroFilter,
404  cleanPiZeros.end(), cleanPiZeros.end()));
405 
406  // Build isolation charged hadrons
409  isolationPFCHs_begin, isolationPFCHs_end);
410 
411  // Add all the stuff in the isolation cone that wasn't in the jet constituents
414  boost::make_filter_iterator(
415  pfChargedJunk, regionalJunk.begin(), regionalJunk.end()),
416  boost::make_filter_iterator(
417  pfChargedJunk, regionalJunk.end(), regionalJunk.end())
418  );
419 
420  // Build isolation neutral hadrons
423  boost::make_filter_iterator(
424  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
425  pfnhs.begin(), pfnhs.end()),
426  boost::make_filter_iterator(
427  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
428  pfnhs.end(), pfnhs.end()));
429 
430  // Add regional stuff not in jet
433  boost::make_filter_iterator(
434  pfNeutralJunk, regionalJunk.begin(), regionalJunk.end()),
435  boost::make_filter_iterator(
436  pfNeutralJunk, regionalJunk.end(), regionalJunk.end())
437  );
438 
439  // Build isolation PiZeros
441  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.begin(),
442  cleanPiZeros.end()),
443  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.end(),
444  cleanPiZeros.end()));
445 
446  // Add regional stuff not in jet
449  boost::make_filter_iterator(
450  pfGammaJunk, regionalJunk.begin(), regionalJunk.end()),
451  boost::make_filter_iterator(
452  pfGammaJunk, regionalJunk.end(), regionalJunk.end())
453  );
454 
455  // Put our built tau in the output - 'false' indicates don't build the
456  // leading candidates, we already did that explicitly above.
457 
458  std::auto_ptr<reco::PFTau> tauPtr = tau.get(false);
459 
460  // Set event vertex position for tau
461  reco::VertexRef primaryVertexRef = primaryVertex(*tauPtr);
462  if ( primaryVertexRef.isNonnull() )
463  tauPtr->setVertex(primaryVertexRef->position());
464 
465  // Set missing tau quantities
466  setTauQuantities(*tauPtr,
469  );
470 
471  output.push_back(tauPtr);
472  return output.release();
473 }
474 }} // end namespace reco::tauk
475 
479  "RecoTauBuilderConePlugin");
static AlgebraicMatrix initialize()
const edm::Handle< PFCandidateCollection > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
std::auto_ptr< reco::PFTau > get(bool setupLeadingCandidates=true)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
const PFCandidatePtr & leadPFChargedHadrCand() const
Definition: PFTau.cc:67
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
std::vector< PFCandidatePtr > pfGammas(const PFJet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
ParameterSet const & getParameterSet(ParameterSetID const &id)
std::vector< reco::PFCandidatePtr > PFCandPtrs
boost::filter_iterator< PiZeroDRFilter, std::vector< RecoTauPiZero >::const_iterator > PiZeroDRFilterIter
Definition: ConeTools.h:58
boost::ptr_vector< reco::PFTau > output_type
double pt() const final
transverse momentum
std::auto_ptr< output_type > return_type
int charge() const final
electric charge
Definition: LeafCandidate.h:91
RecoTauBuilderConePlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
StringObjectFunction< reco::PFTau > signalConeSizeToStore_
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
void addPiZeros(Region region, const InputIterator &begin, const InputIterator &end)
Add a list of pizeros to the input collection.
void setCharge(Charge q) final
set electric charge
Definition: LeafCandidate.h:93
return_type operator()(const reco::PFJetRef &jet, const std::vector< reco::PFRecoTauChargedHadron > &chargedHadrons, const std::vector< RecoTauPiZero > &piZeros, const std::vector< PFCandidatePtr > &regionalExtras) const override
reco::VertexRef primaryVertex(const reco::PFJetRef &jet) const
Get primary vertex associated to this jet.
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
const std::vector< RecoTauPiZero > & signalPiZeroCandidates() const
Retrieve the association of signal region gamma candidates into candidate PiZeros.
Definition: PFTau.cc:127
void setDecayMode(const hadronicDecayMode &)
Definition: PFTau.cc:214
void setTauQuantities(reco::PFTau &aTau, double minAbsPhotonSumPt_insideSignalCone=2.5, double minRelPhotonSumPt_insideSignalCone=0., double minAbsPhotonSumPt_outsideSignalCone=1.e+9, double minRelPhotonSumPt_outsideSignalCone=1.e+9) const
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
T min(T a, T b)
Definition: MathUtil.h:58
void setleadPFCand(const T &cand)
Set leading PF candidate.
void addPFCands(Region region, ParticleType type, const InputIterator &begin, const InputIterator &end)
InputIterator leadPFCand(InputIterator begin, InputIterator end)
StringObjectFunction< reco::PFJet > JetFunc
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
std::vector< PFCandidatePtr > flattenPiZeros(const std::vector< RecoTauPiZero >::const_iterator &, const std::vector< RecoTauPiZero >::const_iterator &)
Flatten a list of pi zeros into a list of there constituent PFCandidates.
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
virtual double pt() const =0
transverse momentum
DeltaRFilter< RecoTauPiZero > PiZeroDRFilter
Definition: ConeTools.h:56
double signalConeSize() const
Size of signal cone.
Definition: PFTau.h:159
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void setleadPFChargedHadrCand(const T &cand)
Set leading PFChargedHadron candidate.
fixed size matrix
hadronicDecayMode
Definition: PFTau.h:36
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
void setleadPFNeutralCand(const T &cand)
Set leading PFGamma candidate.
boost::filter_iterator< PFCandPtrDRFilter, std::vector< PFCandidatePtr >::const_iterator > PFCandPtrDRFilterIter
Definition: ConeTools.h:50
#define DEFINE_EDM_PLUGIN(factory, type, name)
void setPdgId(int pdgId) final
DeltaRPtrFilter< PFCandidatePtr > PFCandPtrDRFilter
Definition: ConeTools.h:48
std::vector< RecoTauPiZero > PiZeroList
def move(src, dest)
Definition: eostools.py:510
const std::vector< reco::PFCandidatePtr > & signalPFChargedHadrCands() const
Charged hadrons in signal region.
Definition: PFTau.cc:80
Transform a pizero to remove given candidates.
PredicateAND< P1, P2 > makePredicateAND(const P1 &p1, const P2 &p2)