CMS 3D CMS Logo

RecoTauCommonUtilities.h
Go to the documentation of this file.
1 #ifndef RecoTauTag_RecoTau_RecoTauCommonUtilities_h
2 #define RecoTauTag_RecoTau_RecoTauCommonUtilities_h
3 
4 /*
5  * RecoTauCommonUtilities - utilities for extracting PFCandidates from
6  * PFJets and summing over collections of PFCandidates.
7  *
8  * Author: Evan K. Friis (UC Davis)
9  *
10  */
11 
12 #include <vector>
13 #include <algorithm>
14 #include <numeric>
15 
22 
23 // Boost helpers
24 #include <boost/iterator/transform_iterator.hpp>
25 #include <boost/iterator/indirect_iterator.hpp>
26 #include <boost/mem_fn.hpp>
27 
28 #include <boost/type_traits/is_base_of.hpp>
29 
30 namespace reco {
31  namespace tau {
32 
34  public:
35  bool operator()(const CandidatePtr& a, const CandidatePtr& b) const { return a->pt() > b->pt(); }
36  };
37 
40  template <typename Iterator>
41  std::vector<CandidatePtr> filterPFCandidates(const Iterator& begin,
42  const Iterator& end,
43  int pdgId,
44  bool sort = true) {
45  std::vector<CandidatePtr> output;
46  for (Iterator iter = begin; iter != end; ++iter) {
47  reco::CandidatePtr ptr(*iter);
48  if (std::abs(ptr->pdgId()) == pdgId)
49  output.push_back(ptr);
50  }
51  if (sort)
52  std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
53  return output;
54  }
55 
58  std::vector<CandidatePtr> pfCandidates(const Jet& jet, int particleId, bool sort = true);
59 
61  std::vector<CandidatePtr> pfCandidates(const Jet& jet, const std::vector<int>& particleIds, bool sort = true);
62 
65  std::vector<CandidatePtr> pfCandidatesByPdgId(const Jet& jet, int pdgId, bool sort = true);
66 
68  std::vector<CandidatePtr> pfCandidatesByPdgId(const Jet& jet, const std::vector<int>& pdgIds, bool sort = true);
69 
71  std::vector<CandidatePtr> pfChargedCands(const Jet& jet, bool sort = true);
72 
74  std::vector<CandidatePtr> pfGammas(const Jet& jet, bool sort = true);
75 
77  std::vector<CandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>::const_iterator&,
78  const std::vector<RecoTauPiZero>::const_iterator&);
79  std::vector<CandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>&);
80 
82  template <typename RefVectorType, typename BaseView>
83  RefVectorType castView(const edm::Handle<BaseView>& view) {
84  typedef typename RefVectorType::value_type OutputRef;
85  // Double check at compile time that the inheritance is okay. It can still
86  // fail at runtime if you pass it the wrong collection.
88  RefVectorType output;
89  size_t nElements = view->size();
90  output.reserve(nElements);
91  // Cast each of our Refs
92  for (size_t i = 0; i < nElements; ++i) {
93  output.push_back(view->refAt(i).template castTo<OutputRef>());
94  }
95  return output;
96  }
97 
98  /*
99  *Given a range over a container of type C, return a new 'end' iterator such
100  *that at max <N> elements are taken. If there are less than N elements in the
101  *array, leave the <end> as it is
102  */
103  template <typename InputIterator>
104  InputIterator takeNElements(const InputIterator& begin, const InputIterator& end, size_t N) {
105  size_t input_size = end - begin;
106  return (N > input_size) ? end : begin + N;
107  }
108 
110  template <typename InputIterator, typename FunctionPtr, typename ReturnType>
111  ReturnType sumPFVector(InputIterator begin, InputIterator end, FunctionPtr func, ReturnType init) {
113  for (InputIterator cand = begin; cand != end; ++cand) {
114  //#define CALL_MEMBER_FN(object,ptrToMember) ((object).*(ptrToMember))
115  output += ((**cand).*(func))();
116  }
117  return output;
118  }
119 
120  template <typename InputIterator>
121  reco::Candidate::LorentzVector sumPFCandP4(InputIterator begin, InputIterator end) {
123  }
124 
126  template <typename InputIterator>
127  double sumPFCandPt(InputIterator begin, InputIterator end) {
128  return sumPFVector(begin, end, &Candidate::pt, 0.0);
129  }
130 
132  template <typename InputIterator>
133  int sumPFCandCharge(InputIterator begin, InputIterator end) {
134  return sumPFVector(begin, end, &Candidate::charge, 0);
135  }
136 
137  template <typename InputIterator>
138  InputIterator leadCand(InputIterator begin, InputIterator end) {
139  double max_pt = 0;
140  InputIterator max_cand = begin;
141  for (InputIterator cand = begin; cand != end; ++cand) {
142  if ((*cand)->pt() > max_pt) {
143  max_pt = (*cand)->pt();
144  max_cand = cand;
145  }
146  }
147  return max_cand;
148  }
149 
151 
152  } // namespace tau
153 } // namespace reco
154 #endif
std::vector< CandidatePtr > pfCandidatesByPdgId(const Jet &jet, int pdgId, bool sort=true)
std::vector< CandidatePtr > pfCandidates(const Jet &jet, int particleId, bool sort=true)
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
math::XYZPointF atECALEntrance(const reco::Candidate *part, double bField)
double sumPFCandPt(InputIterator begin, InputIterator end)
Sum the pT of a collection of PFCandidates.
reco::Candidate::LorentzVector sumPFCandP4(InputIterator begin, InputIterator end)
RefVectorType castView(const edm::Handle< BaseView > &view)
Convert a BaseView (View<T>) to a TRefVector.
int init
Definition: HydjetWrapper.h:64
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
std::vector< CandidatePtr > pfChargedCands(const Jet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
ReturnType sumPFVector(InputIterator begin, InputIterator end, FunctionPtr func, ReturnType init)
Sum the four vectors in a collection of PFCandidates.
InputIterator leadCand(InputIterator begin, InputIterator end)
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.
bool operator()(const CandidatePtr &a, const CandidatePtr &b) const
std::map< DetId, double > ReturnType
Definition: Jet.py:1
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define end
Definition: vmac.h:39
Definition: init.py:1
#define N
Definition: blowfish.cc:9
virtual double pt() const =0
transverse momentum
part
Definition: HCALResponse.h:20
double b
Definition: hdecay.h:118
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
virtual int charge() const =0
electric charge
fixed size matrix
#define begin
Definition: vmac.h:32
double a
Definition: hdecay.h:119
int sumPFCandCharge(InputIterator begin, InputIterator end)
Sum the charge of a collection of PFCandidates.
std::vector< CandidatePtr > filterPFCandidates(const Iterator &begin, const Iterator &end, int pdgId, bool sort=true)
std::vector< CandidatePtr > pfGammas(const Jet &jet, bool sort=true)
Extract all pfGammas from a PFJet.