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 {
31  namespace tau {
32 
33  typedef std::vector<RecoTauPiZero> PiZeroList;
34 
36  public:
39  // Build a tau from a jet
41  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons,
42  const std::vector<RecoTauPiZero>& piZeros,
43  const std::vector<CandidatePtr>& regionalExtras) const override;
44 
45  private:
46  std::unique_ptr<RecoTauQualityCuts> qcuts_;
47 
49 
51 
52  /* String function to extract values from PFJets */
54 
55  // Cone defintions
64 
68 
69  void setTauQuantities(reco::PFTau& aTau,
73  double minRelPhotonSumPt_outsideSignalCone = 1.e+9) const;
74  };
75 
76  // ctor - initialize all of our variables
79  qcuts_(std::make_unique<RecoTauQualityCuts>(
80  pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts"))),
81  usePFLeptonsAsChargedHadrons_(pset.getParameter<bool>("usePFLeptons")),
82  leadObjecPtThreshold_(pset.getParameter<double>("leadObjectPt")),
83  matchingCone_(pset.getParameter<std::string>("matchingCone")),
84  signalConeChargedHadrons_(pset.getParameter<std::string>("signalConeChargedHadrons")),
85  isoConeChargedHadrons_(pset.getParameter<std::string>("isoConeChargedHadrons")),
86  signalConePiZeros_(pset.getParameter<std::string>("signalConePiZeros")),
87  signalConeSizeToStore_( //MB: stored inside PFTau and used to compte photon-pt-out-of-signal-cone and DM hypothsis
88  pset.getParameter<std::string>("signalConePiZeros")),
89  isoConePiZeros_(pset.getParameter<std::string>("isoConePiZeros")),
90  signalConeNeutralHadrons_(pset.getParameter<std::string>("signalConeNeutralHadrons")),
91  isoConeNeutralHadrons_(pset.getParameter<std::string>("isoConeNeutralHadrons")),
92  maxSignalConeChargedHadrons_(pset.getParameter<int>("maxSignalConeChargedHadrons")),
93  minAbsPhotonSumPt_insideSignalCone_(pset.getParameter<double>("minAbsPhotonSumPt_insideSignalCone")),
94  minRelPhotonSumPt_insideSignalCone_(pset.getParameter<double>("minRelPhotonSumPt_insideSignalCone"))
95 
96  {}
97 
98  namespace xclean {
99  // define template specialization for cross-cleaning
100  template <>
102  const cone::CandPtrDRFilterIter& signalTracksBegin, const cone::CandPtrDRFilterIter& signalTracksEnd) {
103  // Get the list of objects we need to clean
104  for (cone::CandPtrDRFilterIter signalTrack = signalTracksBegin; signalTrack != signalTracksEnd; ++signalTrack) {
105  toRemove_.insert(reco::CandidatePtr(*signalTrack));
106  }
107  }
108 
109  template <>
110  inline void CrossCleanPtrs<PiZeroList::const_iterator>::initialize(const PiZeroList::const_iterator& piZerosBegin,
111  const PiZeroList::const_iterator& piZerosEnd) {
112  for (auto const& ptr : flattenPiZeros(piZerosBegin, piZerosEnd)) {
113  toRemove_.insert(CandidatePtr(ptr));
114  }
115  }
116  } // namespace xclean
117 
123  //MB: Set tau quantities which are computed by RecoTauConstructor::get()
124  // method with PFRecoTauChargedHadrons not used here
125 
126  // Set charge
127  int charge = 0;
128  int leadCharge = aTau.leadChargedHadrCand().isNonnull() ? aTau.leadChargedHadrCand()->charge() : 0;
129  const std::vector<reco::CandidatePtr>& pfChs = aTau.signalChargedHadrCands();
130  for (const reco::CandidatePtr& pfCh : pfChs) {
131  charge += pfCh->charge();
132  }
133  charge = charge == 0 ? leadCharge : charge;
134  aTau.setCharge(charge);
135 
136  // Set PDG id
137  aTau.setPdgId(aTau.charge() < 0 ? 15 : -15);
138 
139  // Set Decay Mode
141  unsigned int nPiZeros = 0;
142  const std::vector<RecoTauPiZero>& piZeros = aTau.signalPiZeroCandidates();
143  for (const RecoTauPiZero& piZero : piZeros) {
144  double photonSumPt_insideSignalCone = 0.;
145  double photonSumPt_outsideSignalCone = 0.;
146  int numPhotons = piZero.numberOfDaughters();
147  for (int idxPhoton = 0; idxPhoton < numPhotons; ++idxPhoton) {
148  const reco::Candidate* photon = piZero.daughter(idxPhoton);
149  double dR = deltaR(photon->p4(), aTau.p4());
150  if (dR < aTau.signalConeSize()) {
151  photonSumPt_insideSignalCone += photon->pt();
152  } else {
153  photonSumPt_outsideSignalCone += photon->pt();
154  }
155  }
156  if (photonSumPt_insideSignalCone > minAbsPhotonSumPt_insideSignalCone ||
157  photonSumPt_insideSignalCone > (minRelPhotonSumPt_insideSignalCone * aTau.pt()) ||
158  photonSumPt_outsideSignalCone > minAbsPhotonSumPt_outsideSignalCone ||
159  photonSumPt_outsideSignalCone > (minRelPhotonSumPt_outsideSignalCone * aTau.pt()))
160  ++nPiZeros;
161  }
162  // Find the maximum number of PiZeros our parameterization can hold
163  const unsigned int maxPiZeros = PFTau::kOneProngNPiZero;
164  // Determine our track index
165  unsigned int nCharged = pfChs.size();
166  if (nCharged > 0) {
167  unsigned int trackIndex = (nCharged - 1) * (maxPiZeros + 1);
168  // Check if we handle the given number of tracks
169  if (trackIndex >= PFTau::kRareDecayMode)
171  else
172  dm = static_cast<PFTau::hadronicDecayMode>(trackIndex + std::min(nPiZeros, maxPiZeros));
173  } else {
174  dm = PFTau::kNull;
175  }
176  aTau.setDecayMode(dm);
177  return;
178  }
179 
181  const reco::JetBaseRef& jet,
182  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons,
183  const std::vector<RecoTauPiZero>& piZeros,
184  const std::vector<CandidatePtr>& regionalExtras) const {
185  //std::cout << "<RecoTauBuilderConePlugin::operator()>:" << std::endl;
186  //std::cout << " jet: Pt = " << jet->pt() << ", eta = " << jet->eta() << ", phi = " << jet->phi() << std::endl;
187 
188  // Get access to our cone tools
189  using namespace cone;
190  // Define output. We only produce one tau per jet for the cone algo.
192 
193  // Our tau builder - the true indicates to automatically copy gamma candidates
194  // from the pizeros.
196  getPFCands(),
197  true,
201  // Setup our quality cuts to use the current vertex (supplied by base class)
202  qcuts_->setPV(primaryVertex(jet));
203 
204  typedef std::vector<CandidatePtr> CandPtrs;
205 
206  // Get the PF Charged hadrons + quality cuts
207  CandPtrs pfchs;
209  pfchs = qcuts_->filterCandRefs(pfCandidatesByPdgId(*jet, 211));
210  } else {
211  // Check if we want to include electrons in muons in "charged hadron"
212  // collection. This is the preferred behavior, as the PF lepton selections
213  // are very loose.
214  pfchs = qcuts_->filterCandRefs(pfChargedCands(*jet));
215  }
216 
217  // CV: sort collection of PF Charged hadrons by descending Pt
218  std::sort(pfchs.begin(), pfchs.end(), SortPFCandsDescendingPt());
219 
220  // Get the PF gammas
221  CandPtrs pfGammas = qcuts_->filterCandRefs(pfCandidatesByPdgId(*jet, 22));
222  // Neutral hadrons
223  CandPtrs pfnhs = qcuts_->filterCandRefs(pfCandidatesByPdgId(*jet, 130));
224 
225  // All the extra junk
226  CandPtrs regionalJunk = qcuts_->filterCandRefs(regionalExtras);
227 
228  /***********************************************
229  ****** Lead Candidate Finding **********
230  ***********************************************/
231 
232  // Define our matching cone and filters
233  double matchingCone = matchingCone_(*jet);
234  CandPtrDRFilter matchingConeFilter(jet->p4(), 0, matchingCone);
235 
236  // Find the maximum PFCharged hadron in the matching cone. The call to
237  // PFCandidates always a sorted list, so we can just take the first if it
238  // if it exists.
239  CandidatePtr leadPFCH;
240  CandPtrs::iterator leadPFCH_iter = std::find_if(pfchs.begin(), pfchs.end(), matchingConeFilter);
241 
242  if (leadPFCH_iter != pfchs.end()) {
243  leadPFCH = *leadPFCH_iter;
244  // Set leading candidate
245  tau.setleadChargedHadrCand(leadPFCH);
246  } else {
247  // If there is no leading charged candidate at all, return nothing - the
248  // producer class that owns the plugin will build a null tau if desired.
249  return output;
250  }
251 
252  // Find the leading neutral candidate
253  CandidatePtr leadPFGamma;
254  CandPtrs::iterator leadPFGamma_iter = std::find_if(pfGammas.begin(), pfGammas.end(), matchingConeFilter);
255 
256  if (leadPFGamma_iter != pfGammas.end()) {
257  leadPFGamma = *leadPFGamma_iter;
258  // Set leading neutral candidate
259  tau.setleadNeutralCand(leadPFGamma);
260  }
261 
262  CandidatePtr leadPFCand;
263  // Always use the leadPFCH if it is above our threshold
264  if (leadPFCH.isNonnull() && leadPFCH->pt() > leadObjecPtThreshold_) {
265  leadPFCand = leadPFCH;
266  } else if (leadPFGamma.isNonnull() && leadPFGamma->pt() > leadObjecPtThreshold_) {
267  // Otherwise use the lead Gamma if it is above threshold
268  leadPFCand = leadPFGamma;
269  } else {
270  // If both are too low PT, just take the charged one
271  leadPFCand = leadPFCH;
272  }
273 
274  tau.setleadCand(leadPFCand);
275 
276  // Our cone axis is defined about the lead charged hadron
277  reco::Candidate::LorentzVector coneAxis = leadPFCH->p4();
278 
279  /***********************************************
280  ****** Cone Construction **********
281  ***********************************************/
282 
283  // Define the signal and isolation cone sizes for this jet and build filters
284  // to select elements in the given DeltaR regions
285 
286  CandPtrDRFilter signalConePFCHFilter(coneAxis, -0.1, signalConeChargedHadrons_(*jet));
287  CandPtrDRFilter signalConePFNHFilter(coneAxis, -0.1, signalConeNeutralHadrons_(*jet));
288  PiZeroDRFilter signalConePiZeroFilter(coneAxis, -0.1, signalConePiZeros_(*jet));
289 
290  CandPtrDRFilter isoConePFCHFilter(coneAxis, signalConeChargedHadrons_(*jet), isoConeChargedHadrons_(*jet));
291  CandPtrDRFilter isoConePFGammaFilter(coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
292  CandPtrDRFilter isoConePFNHFilter(coneAxis, signalConeNeutralHadrons_(*jet), isoConeNeutralHadrons_(*jet));
293  PiZeroDRFilter isoConePiZeroFilter(coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
294 
295  // Additionally make predicates to select the different PF object types
296  // of the regional junk objects to add to the iso cone.
298 
299  xclean::FilterCandByAbsPdgId pfchCandSelector(211);
300  xclean::FilterCandByAbsPdgId pfgammaCandSelector(22);
301  xclean::FilterCandByAbsPdgId pfnhCandSelector(130);
302 
303  // Predicate to select the regional junk in the iso cone by PF id
304  RegionalJunkConeAndIdFilter pfChargedJunk(pfchCandSelector, // select charged stuff from junk
305  isoConePFCHFilter // only take those in iso cone
306  );
307 
308  RegionalJunkConeAndIdFilter pfGammaJunk(pfgammaCandSelector, // select gammas from junk
309  isoConePFGammaFilter // only take those in iso cone
310  );
311 
312  RegionalJunkConeAndIdFilter pfNeutralJunk(pfnhCandSelector, // select neutral stuff from junk
313  isoConePFNHFilter // select stuff in iso cone
314  );
315 
316  // Build filter iterators select the signal charged stuff.
317  CandPtrDRFilterIter signalPFCHCands_begin(signalConePFCHFilter, pfchs.begin(), pfchs.end());
318  CandPtrDRFilterIter signalPFCHCands_end(signalConePFCHFilter, pfchs.end(), pfchs.end());
319  CandPtrs signalPFCHs;
320  int numSignalPFCHs = 0;
321  CandPtrs isolationPFCHs;
322  for (CandPtrDRFilterIter iter = signalPFCHCands_begin; iter != signalPFCHCands_end; ++iter) {
323  if (numSignalPFCHs < maxSignalConeChargedHadrons_ || maxSignalConeChargedHadrons_ == -1) {
324  //std::cout << "adding signalPFCH #" << numSignalPFCHs << ": Pt = " << (*iter)->pt() << ", eta = " << (*iter)->eta() << ", phi = " << (*iter)->phi() << std::endl;
325  signalPFCHs.push_back(*iter);
326  ++numSignalPFCHs;
327  } else {
328  //std::cout << "maxSignalConeChargedHadrons reached"
329  isolationPFCHs.push_back(*iter);
330  }
331  }
332  CandPtrs::const_iterator signalPFCHs_begin = signalPFCHs.begin();
333  CandPtrs::const_iterator signalPFCHs_end = signalPFCHs.end();
334 
335  // Cross clean pi zero content using signal cone charged hadron constituents.
336  xclean::CrossCleanPiZeros<CandPtrDRFilterIter> piZeroXCleaner(signalPFCHCands_begin, signalPFCHCands_end);
337  std::vector<reco::RecoTauPiZero> cleanPiZeros = piZeroXCleaner(piZeros);
338 
339  // For the rest of the constituents, we need to filter any constituents that
340  // are already contained in the pizeros (i.e. electrons)
341  xclean::CrossCleanPtrs<PiZeroList::const_iterator> pfCandXCleaner(cleanPiZeros.begin(), cleanPiZeros.end());
342 
343  auto isolationPFCHCands_begin(boost::make_filter_iterator(
344  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner), pfchs.begin(), pfchs.end()));
345  auto isolationPFCHCands_end(boost::make_filter_iterator(
346  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner), pfchs.end(), pfchs.end()));
347  for (auto iter = isolationPFCHCands_begin; iter != isolationPFCHCands_end; ++iter) {
348  isolationPFCHs.push_back(*iter);
349  }
350  CandPtrs::const_iterator isolationPFCHs_begin = isolationPFCHs.begin();
351  CandPtrs::const_iterator isolationPFCHs_end = isolationPFCHs.end();
352 
353  // Build signal charged hadrons
354  tau.addPFCands(
355  RecoTauConstructor::kSignal, RecoTauConstructor::kChargedHadron, signalPFCHs_begin, signalPFCHs_end);
356 
357  tau.addPFCands(RecoTauConstructor::kSignal,
359  boost::make_filter_iterator(
360  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner), pfnhs.begin(), pfnhs.end()),
361  boost::make_filter_iterator(
362  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner), pfnhs.end(), pfnhs.end()));
363 
364  // Build signal PiZeros
365  tau.addPiZeros(RecoTauConstructor::kSignal,
366  PiZeroDRFilterIter(signalConePiZeroFilter, cleanPiZeros.begin(), cleanPiZeros.end()),
367  PiZeroDRFilterIter(signalConePiZeroFilter, cleanPiZeros.end(), cleanPiZeros.end()));
368 
369  // Build isolation charged hadrons
370  tau.addPFCands(
371  RecoTauConstructor::kIsolation, RecoTauConstructor::kChargedHadron, isolationPFCHs_begin, isolationPFCHs_end);
372 
373  // Add all the stuff in the isolation cone that wasn't in the jet constituents
376  boost::make_filter_iterator(pfChargedJunk, regionalJunk.begin(), regionalJunk.end()),
377  boost::make_filter_iterator(pfChargedJunk, regionalJunk.end(), regionalJunk.end()));
378 
379  // Build isolation neutral hadrons
382  boost::make_filter_iterator(
383  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner), pfnhs.begin(), pfnhs.end()),
384  boost::make_filter_iterator(
385  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner), pfnhs.end(), pfnhs.end()));
386 
387  // Add regional stuff not in jet
390  boost::make_filter_iterator(pfNeutralJunk, regionalJunk.begin(), regionalJunk.end()),
391  boost::make_filter_iterator(pfNeutralJunk, regionalJunk.end(), regionalJunk.end()));
392 
393  // Build isolation PiZeros
395  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.begin(), cleanPiZeros.end()),
396  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.end(), cleanPiZeros.end()));
397 
398  // Add regional stuff not in jet
401  boost::make_filter_iterator(pfGammaJunk, regionalJunk.begin(), regionalJunk.end()),
402  boost::make_filter_iterator(pfGammaJunk, regionalJunk.end(), regionalJunk.end()));
403 
404  // Put our built tau in the output - 'false' indicates don't build the
405  // leading candidates, we already did that explicitly above.
406 
407  std::unique_ptr<reco::PFTau> tauPtr = tau.get(false);
408 
409  // Set event vertex position for tau
410  reco::VertexRef primaryVertexRef = primaryVertex(*tauPtr);
411  if (primaryVertexRef.isNonnull())
412  tauPtr->setVertex(primaryVertexRef->position());
413 
414  // Set missing tau quantities
416 
417  output.push_back(std::move(tauPtr));
418  return output;
419  }
420  } // namespace tau
421 } // namespace reco
422 
std::vector< CandidatePtr > pfCandidatesByPdgId(const Jet &jet, int pdgId, bool sort=true)
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
std::unique_ptr< RecoTauQualityCuts > qcuts_
double signalConeSize() const
Size of signal cone.
Definition: PFTau.h:174
std::vector< std::unique_ptr< reco::PFTau > > output_type
double pt() const final
transverse momentum
return_type operator()(const reco::JetBaseRef &jet, const std::vector< reco::PFRecoTauChargedHadron > &chargedHadrons, const std::vector< RecoTauPiZero > &piZeros, const std::vector< CandidatePtr > &regionalExtras) const override
const std::vector< RecoTauPiZero > & signalPiZeroCandidates() const
Retrieve the association of signal region gamma candidates into candidate PiZeros.
Definition: PFTau.cc:240
const CandidatePtr & leadChargedHadrCand() const
Definition: PFTau.cc:63
const edm::Handle< edm::View< reco::Candidate > > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float charge(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:73
RecoTauBuilderConePlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
std::vector< CandidatePtr > pfChargedCands(const Jet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
const LorentzVector & p4() const final
four-momentum Lorentz vector
StringObjectFunction< reco::PFTau > signalConeSizeToStore_
std::vector< CandidatePtr > 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.
void setCharge(Charge q) final
set electric charge
reco::VertexRef primaryVertex(const reco::JetBaseRef &jet) const
Get primary vertex associated to this jet.
std::vector< reco::CandidatePtr > CandPtrs
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
void setDecayMode(const hadronicDecayMode &)
Definition: PFTau.cc:327
boost::filter_iterator< CandPtrDRFilter, std::vector< CandidatePtr >::const_iterator > CandPtrDRFilterIter
Definition: ConeTools.h:47
boost::filter_iterator< PiZeroDRFilter, std::vector< RecoTauPiZero >::const_iterator > PiZeroDRFilterIter
Definition: ConeTools.h:54
void initialize(const PtrIter &chargedHadronsBegin, const PtrIter &chargedHadronsEnd)
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
DeltaRPtrFilter< CandidatePtr > CandPtrDRFilter
Definition: ConeTools.h:46
DeltaRFilter< RecoTauPiZero > PiZeroDRFilter
Definition: ConeTools.h:53
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
ParameterSet const & getParameterSet(ParameterSetID const &id)
const std::vector< reco::CandidatePtr > & signalChargedHadrCands() const
Charged hadrons in signal region.
Definition: PFTau.cc:76
fixed size matrix
hadronicDecayMode
Definition: PFTau.h:38
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: output.py:1
StringObjectFunction< reco::Jet > JetFunc
void setPdgId(int pdgId) final
void initialize(const PtrIter &particlesBegin, const PtrIter &particlesEnd)
std::vector< RecoTauPiZero > PiZeroList
def move(src, dest)
Definition: eostools.py:511
Transform a pizero to remove given candidates.
int charge() const final
electric charge
std::vector< CandidatePtr > pfGammas(const Jet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
PredicateAND< P1, P2 > makePredicateAND(const P1 &p1, const P2 &p2)