CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauCommonUtilities.cc
Go to the documentation of this file.
2 
7 
8 #include <algorithm>
9 
10 typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
11 typedef PFCandPtrs::iterator PFCandIter;
12 
13 namespace reco { namespace tau {
14 
15 std::vector<PFCandidatePtr>
16 flattenPiZeros(const std::vector<RecoTauPiZero>::const_iterator& piZerosBegin, const std::vector<RecoTauPiZero>::const_iterator& piZerosEnd) {
17  std::vector<PFCandidatePtr> output;
18 
19  for(std::vector<RecoTauPiZero>::const_iterator piZero = piZerosBegin;
20  piZero != piZerosEnd; ++piZero) {
21  for(size_t iDaughter = 0; iDaughter < piZero->numberOfDaughters();
22  ++iDaughter) {
23  output.push_back(PFCandidatePtr(piZero->daughterPtr(iDaughter)));
24  }
25  }
26  return output;
27 }
28 
29 std::vector<PFCandidatePtr>
30 flattenPiZeros(const std::vector<RecoTauPiZero>& piZeros) {
31  return flattenPiZeros(piZeros.begin(), piZeros.end());
32 }
33 
34 std::vector<reco::PFCandidatePtr> pfCandidates(const reco::PFJet& jet,
35  int particleId, bool sort) {
36  PFCandPtrs pfCands = jet.getPFConstituents();
37  PFCandPtrs selectedPFCands = filterPFCandidates(
38  pfCands.begin(), pfCands.end(), particleId, sort);
39  return selectedPFCands;
40 }
41 
42 std::vector<reco::PFCandidatePtr> pfCandidates(const reco::PFJet& jet,
43  const std::vector<int>& particleIds, bool sort) {
44  PFCandPtrs&& pfCands = jet.getPFConstituents();
46  // Get each desired candidate type, unsorted for now
47  for(std::vector<int>::const_iterator particleId = particleIds.begin();
48  particleId != particleIds.end(); ++particleId) {
49  PFCandPtrs&& selectedPFCands = filterPFCandidates(pfCands.begin(), pfCands.end(), *particleId, false);
50  output.insert(output.end(), selectedPFCands.begin(), selectedPFCands.end());
51  }
52  if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
53  return output;
54 }
55 
56 std::vector<reco::PFCandidatePtr> pfGammas(const reco::PFJet& jet, bool sort) {
57  return pfCandidates(jet, reco::PFCandidate::gamma, sort);
58 }
59 
60 std::vector<reco::PFCandidatePtr> pfChargedCands(const reco::PFJet& jet,
61  bool sort) {
62  PFCandPtrs&& pfCands = jet.getPFConstituents();
64  PFCandPtrs&& mus = filterPFCandidates(pfCands.begin(), pfCands.end(), reco::PFCandidate::mu, false);
65  PFCandPtrs&& es = filterPFCandidates(pfCands.begin(), pfCands.end(), reco::PFCandidate::e, false);
66  PFCandPtrs&& chs = filterPFCandidates(pfCands.begin(), pfCands.end(), reco::PFCandidate::h, false);
67  output.reserve(mus.size() + es.size() + chs.size());
68  output.insert(output.end(), mus.begin(), mus.end());
69  output.insert(output.end(), es.begin(), es.end());
70  output.insert(output.end(), chs.begin(), chs.end());
71  if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
72  return output;
73 }
74 
75 } }
std::vector< PFCandidatePtr > pfGammas(const PFJet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
std::vector< reco::PFCandidatePtr > PFCandPtrs
std::vector< PFCandidatePtr > filterPFCandidates(const Iterator &begin, const Iterator &end, int particleId, bool sort=true)
Jets made from PFObjects.
Definition: PFJet.h:21
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
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.
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:52
PFCandPtrs::iterator PFCandIter
edm::Ptr< PFCandidate > PFCandidatePtr
persistent Ptr to a PFCandidate