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  * $Id $
11  */
12 
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 
32 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
33  int particleId, bool sort=true);
34 
36 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
37  const std::vector<int>& particleIds,
38  bool sort=true);
39 
41 std::vector<PFCandidatePtr> pfChargedCands(const PFJet& jet, bool sort=true);
42 
44 std::vector<PFCandidatePtr> pfGammas(const PFJet& jet, bool sort=true);
45 
47 std::vector<PFCandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>&);
48 
50 template<typename RefVectorType, typename BaseView>
51 RefVectorType castView(const edm::Handle<BaseView>& view) {
52  typedef typename RefVectorType::value_type OutputRef;
53  // Double check at compile time that the inheritance is okay. It can still
54  // fail at runtime if you pass it the wrong collection.
55  BOOST_STATIC_ASSERT(
56  (boost::is_base_of<typename BaseView::value_type,
57  typename RefVectorType::member_type>::value));
58  RefVectorType output;
59  size_t nElements = view->size();
60  output.reserve(nElements);
61  // Cast each of our Refs
62  for (size_t i = 0; i < nElements; ++i) {
63  output.push_back(view->refAt(i).template castTo<OutputRef>());
64  }
65  return output;
66 }
67 
68 /*
69  *Given a range over a container of type C, return a new 'end' iterator such
70  *that at max <N> elements are taken. If there are less than N elements in the
71  *array, leave the <end> as it is
72  */
73 template<typename InputIterator> InputIterator takeNElements(
74  const InputIterator& begin, const InputIterator& end, size_t N) {
75  size_t input_size = end - begin;
76  return (N > input_size) ? end : begin + N;
77 }
78 
80 template<typename InputIterator, typename FunctionPtr, typename ReturnType>
81  ReturnType sumPFVector(InputIterator begin, InputIterator end,
82  FunctionPtr func, ReturnType init) {
83  ReturnType output = init;
84  for(InputIterator cand = begin; cand != end; ++cand) {
85  //#define CALL_MEMBER_FN(object,ptrToMember) ((object).*(ptrToMember))
86  output += ((**cand).*(func))();
87  }
88  return output;
89  }
90 
91 template<typename InputIterator> reco::Candidate::LorentzVector sumPFCandP4(
92  InputIterator begin, InputIterator end) {
93  return sumPFVector(begin, end, &PFCandidate::p4,
95 }
96 
98 template<typename InputIterator> double sumPFCandPt(InputIterator begin,
99  InputIterator end) {
100  return sumPFVector(begin, end, &PFCandidate::pt, 0.0);
101  }
102 
104 template<typename InputIterator> int sumPFCandCharge(InputIterator begin,
105  InputIterator end) {
106  return sumPFVector(begin, end, &PFCandidate::charge, 0);
107  }
108 
109 template<typename InputIterator> InputIterator leadPFCand(InputIterator begin,
110  InputIterator end) {
111  double max_pt = 0;
112  InputIterator max_cand = begin;
113  for(InputIterator cand = begin; cand != end; ++cand) {
114  if( (*cand)->pt() > max_pt ) {
115  max_pt = (*cand)->pt();
116  max_cand = cand;
117  }
118  }
119  return max_cand;
120  }
121 }} // end namespace reco::tau
122 #endif
int i
Definition: DBlmapReader.cc:9
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
std::vector< PFCandidatePtr > flattenPiZeros(const std::vector< RecoTauPiZero > &)
Flatten a list of pi zeros into a list of there constituent PFCandidates.
double sumPFCandPt(InputIterator begin, InputIterator end)
Sum the PT of a collection of PFCandidates.
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:63
Container::value_type value_type
ReturnType sumPFVector(InputIterator begin, InputIterator end, FunctionPtr func, ReturnType init)
Sum the four vectors in a collection of PFCandidates.
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
virtual int charge() const
electric charge
#define end
Definition: vmac.h:38
InputIterator leadPFCand(InputIterator begin, InputIterator end)
virtual double pt() const
transverse momentum
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:39
#define begin
Definition: vmac.h:31
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
int sumPFCandCharge(InputIterator begin, InputIterator end)
Sum the PT of a collection of PFCandidates.