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 
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>&);
70 
72 template<typename RefVectorType, typename BaseView>
73 RefVectorType castView(const edm::Handle<BaseView>& view) {
74  typedef typename RefVectorType::value_type OutputRef;
75  // Double check at compile time that the inheritance is okay. It can still
76  // fail at runtime if you pass it the wrong collection.
77  BOOST_STATIC_ASSERT(
78  (boost::is_base_of<typename BaseView::value_type,
79  typename RefVectorType::member_type>::value));
80  RefVectorType output;
81  size_t nElements = view->size();
82  output.reserve(nElements);
83  // Cast each of our Refs
84  for (size_t i = 0; i < nElements; ++i) {
85  output.push_back(view->refAt(i).template castTo<OutputRef>());
86  }
87  return output;
88 }
89 
90 /*
91  *Given a range over a container of type C, return a new 'end' iterator such
92  *that at max <N> elements are taken. If there are less than N elements in the
93  *array, leave the <end> as it is
94  */
95 template<typename InputIterator> InputIterator takeNElements(
96  const InputIterator& begin, const InputIterator& end, size_t N) {
97  size_t input_size = end - begin;
98  return (N > input_size) ? end : begin + N;
99 }
100 
102 template<typename InputIterator, typename FunctionPtr, typename ReturnType>
103  ReturnType sumPFVector(InputIterator begin, InputIterator end,
104  FunctionPtr func, ReturnType init) {
105  ReturnType output = init;
106  for(InputIterator cand = begin; cand != end; ++cand) {
107  //#define CALL_MEMBER_FN(object,ptrToMember) ((object).*(ptrToMember))
108  output += ((**cand).*(func))();
109  }
110  return output;
111  }
112 
113 template<typename InputIterator> reco::Candidate::LorentzVector sumPFCandP4(
114  InputIterator begin, InputIterator end) {
115  return sumPFVector(begin, end, &PFCandidate::p4,
117 }
118 
120 template<typename InputIterator> double sumPFCandPt(InputIterator begin,
121  InputIterator end) {
122  return sumPFVector(begin, end, &PFCandidate::pt, 0.0);
123  }
124 
126 template<typename InputIterator> int sumPFCandCharge(InputIterator begin,
127  InputIterator end) {
128  return sumPFVector(begin, end, &PFCandidate::charge, 0);
129  }
130 
131 template<typename InputIterator> InputIterator leadPFCand(InputIterator begin,
132  InputIterator end) {
133  double max_pt = 0;
134  InputIterator max_cand = begin;
135  for(InputIterator cand = begin; cand != end; ++cand) {
136  if( (*cand)->pt() > max_pt ) {
137  max_pt = (*cand)->pt();
138  max_cand = cand;
139  }
140  }
141  return max_cand;
142  }
143 }} // end namespace reco::tau
144 #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
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:22
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
virtual int charge() const
electric charge
#define end
Definition: vmac.h:38
Container::value_type value_type
InputIterator leadPFCand(InputIterator begin, InputIterator end)
virtual double pt() const
transverse momentum
double b
Definition: hdecay.h:120
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.
double a
Definition: hdecay.h:121
bool operator()(const PFCandidatePtr &a, const PFCandidatePtr &b) const
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
int sumPFCandCharge(InputIterator begin, InputIterator end)
Sum the PT of a collection of PFCandidates.