CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
16 #include <vector>
17 #include <algorithm>
18 #include <numeric>
19 
20 // Boost helpers
21 #include <boost/iterator/transform_iterator.hpp>
22 #include <boost/iterator/indirect_iterator.hpp>
23 #include <boost/mem_fn.hpp>
24 
25 #include <boost/type_traits/is_base_of.hpp>
26 #include <boost/static_assert.hpp>
27 
28 namespace reco { namespace tau {
29 
31  public:
32  bool operator()(const PFCandidatePtr& a, const PFCandidatePtr& b) const {
33  return a->pt() > b->pt();
34  }
35 };
36 
39 template<typename Iterator> std::vector<PFCandidatePtr>
41  int particleId, bool sort=true) {
42  std::vector<PFCandidatePtr> output;
43  for(Iterator iter = begin; iter != end; ++iter) {
44  reco::PFCandidatePtr ptr(*iter);
45  if (ptr->particleId() == particleId)
46  output.push_back(ptr);
47  }
48  if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
49  return output;
50 }
51 
54 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
55  int particleId, bool sort=true);
56 
58 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
59  const std::vector<int>& particleIds,
60  bool sort=true);
61 
63 std::vector<PFCandidatePtr> pfChargedCands(const PFJet& jet, bool sort=true);
64 
66 std::vector<PFCandidatePtr> pfGammas(const PFJet& jet, bool sort=true);
67 
69 std::vector<PFCandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>::const_iterator&, const std::vector<RecoTauPiZero>::const_iterator&);
70 std::vector<PFCandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>&);
71 
73 template<typename RefVectorType, typename BaseView>
74 RefVectorType castView(const edm::Handle<BaseView>& view) {
75  typedef typename RefVectorType::value_type OutputRef;
76  // Double check at compile time that the inheritance is okay. It can still
77  // fail at runtime if you pass it the wrong collection.
78  BOOST_STATIC_ASSERT(
79  (boost::is_base_of<typename BaseView::value_type,
80  typename RefVectorType::member_type>::value));
81  RefVectorType output;
82  size_t nElements = view->size();
83  output.reserve(nElements);
84  // Cast each of our Refs
85  for (size_t i = 0; i < nElements; ++i) {
86  output.push_back(view->refAt(i).template castTo<OutputRef>());
87  }
88  return output;
89 }
90 
91 /*
92  *Given a range over a container of type C, return a new 'end' iterator such
93  *that at max <N> elements are taken. If there are less than N elements in the
94  *array, leave the <end> as it is
95  */
96 template<typename InputIterator> InputIterator takeNElements(
97  const InputIterator& begin, const InputIterator& end, size_t N) {
98  size_t input_size = end - begin;
99  return (N > input_size) ? end : begin + N;
100 }
101 
103 template<typename InputIterator, typename FunctionPtr, typename ReturnType>
104  ReturnType sumPFVector(InputIterator begin, InputIterator end,
105  FunctionPtr func, ReturnType init) {
106  ReturnType output = init;
107  for(InputIterator cand = begin; cand != end; ++cand) {
108  //#define CALL_MEMBER_FN(object,ptrToMember) ((object).*(ptrToMember))
109  output += ((**cand).*(func))();
110  }
111  return output;
112  }
113 
114 template<typename InputIterator> reco::Candidate::LorentzVector sumPFCandP4(
115  InputIterator begin, InputIterator end) {
116  return sumPFVector(begin, end, &PFCandidate::p4,
118 }
119 
121 template<typename InputIterator> double sumPFCandPt(InputIterator begin,
122  InputIterator end) {
123  return sumPFVector(begin, end, &PFCandidate::pt, 0.0);
124  }
125 
127 template<typename InputIterator> int sumPFCandCharge(InputIterator begin,
128  InputIterator end) {
129  return sumPFVector(begin, end, &PFCandidate::charge, 0);
130  }
131 
132 template<typename InputIterator> InputIterator leadPFCand(InputIterator begin,
133  InputIterator end) {
134  double max_pt = 0;
135  InputIterator max_cand = begin;
136  for(InputIterator cand = begin; cand != end; ++cand) {
137  if( (*cand)->pt() > max_pt ) {
138  max_pt = (*cand)->pt();
139  max_cand = cand;
140  }
141  }
142  return max_cand;
143  }
144 }} // end namespace reco::tau
145 #endif
int i
Definition: DBlmapReader.cc:9
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
double sumPFCandPt(InputIterator begin, InputIterator end)
Sum the pT of a collection of PFCandidates.
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
reco::Candidate::LorentzVector sumPFCandP4(InputIterator begin, InputIterator end)
std::vector< PFCandidatePtr > pfGammas(const PFJet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
RefVectorType castView(const edm::Handle< BaseView > &view)
Convert a BaseView (View&lt;T&gt;) to a TRefVector.
int init
Definition: HydjetWrapper.h:62
std::vector< PFCandidatePtr > filterPFCandidates(const Iterator &begin, const Iterator &end, int particleId, bool sort=true)
ReturnType sumPFVector(InputIterator begin, InputIterator end, FunctionPtr func, ReturnType init)
Sum the four vectors in a collection of PFCandidates.
Jets made from PFObjects.
Definition: PFJet.h:21
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
#define end
Definition: vmac.h:37
Container::value_type value_type
InputIterator leadPFCand(InputIterator begin, InputIterator end)
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.
#define N
Definition: blowfish.cc:9
virtual int charge() const GCC11_FINAL
electric charge
double b
Definition: hdecay.h:120
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
#define begin
Definition: vmac.h:30
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
double a
Definition: hdecay.h:121
string func
Definition: statics.py:48
bool operator()(const PFCandidatePtr &a, const PFCandidatePtr &b) const
virtual float pt() const GCC11_FINAL
transverse momentum
int sumPFCandCharge(InputIterator begin, InputIterator end)
Sum the charge of a collection of PFCandidates.