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 { namespace tau {
31 
33  public:
34  bool operator()(const CandidatePtr& a, const CandidatePtr& b) const {
35  return a->pt() > b->pt();
36  }
37 };
38 
41 template<typename Iterator> std::vector<CandidatePtr>
43  int pdgId, bool sort=true) {
44  std::vector<CandidatePtr> output;
45  for(Iterator iter = begin; iter != end; ++iter) {
46  reco::CandidatePtr ptr(*iter);
47  if (std::abs(ptr->pdgId()) == pdgId)
48  output.push_back(ptr);
49  }
50  if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
51  return output;
52 }
53 
56 std::vector<CandidatePtr> pfCandidates(const Jet& jet,
57  int particleId, bool sort=true);
58 
60 std::vector<CandidatePtr> pfCandidates(const Jet& jet,
61  const std::vector<int>& particleIds,
62  bool sort=true);
63 
66 std::vector<CandidatePtr> pfCandidatesByPdgId(const Jet& jet,
67  int pdgId, bool sort=true);
68 
70 std::vector<CandidatePtr> pfCandidatesByPdgId(const Jet& jet,
71  const std::vector<int>& pdgIds,
72  bool sort=true);
73 
75 std::vector<CandidatePtr> pfChargedCands(const Jet& jet, bool sort=true);
76 
78 std::vector<CandidatePtr> pfGammas(const Jet& jet, bool sort=true);
79 
81 std::vector<CandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>::const_iterator&, const std::vector<RecoTauPiZero>::const_iterator&);
82 std::vector<CandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>&);
83 
85 template<typename RefVectorType, typename BaseView>
86 RefVectorType castView(const edm::Handle<BaseView>& view) {
87  typedef typename RefVectorType::value_type OutputRef;
88  // Double check at compile time that the inheritance is okay. It can still
89  // fail at runtime if you pass it the wrong collection.
90  static_assert(
91  (boost::is_base_of<typename BaseView::value_type,
92  typename RefVectorType::member_type>::value));
93  RefVectorType output;
94  size_t nElements = view->size();
95  output.reserve(nElements);
96  // Cast each of our Refs
97  for (size_t i = 0; i < nElements; ++i) {
98  output.push_back(view->refAt(i).template castTo<OutputRef>());
99  }
100  return output;
101 }
102 
103 
104 /*
105  *Given a range over a container of type C, return a new 'end' iterator such
106  *that at max <N> elements are taken. If there are less than N elements in the
107  *array, leave the <end> as it is
108  */
109 template<typename InputIterator> InputIterator takeNElements(
110  const InputIterator& begin, const InputIterator& end, size_t N) {
111  size_t input_size = end - begin;
112  return (N > input_size) ? end : begin + N;
113 }
114 
116 template<typename InputIterator, typename FunctionPtr, typename ReturnType>
117  ReturnType sumPFVector(InputIterator begin, InputIterator end,
118  FunctionPtr func, ReturnType init) {
120  for(InputIterator cand = begin; cand != end; ++cand) {
121  //#define CALL_MEMBER_FN(object,ptrToMember) ((object).*(ptrToMember))
122  output += ((**cand).*(func))();
123  }
124  return output;
125  }
126 
127 template<typename InputIterator> reco::Candidate::LorentzVector sumPFCandP4(
128  InputIterator begin, InputIterator end) {
129  return sumPFVector(begin, end, &Candidate::p4,
131 }
132 
134 template<typename InputIterator> double sumPFCandPt(InputIterator begin,
135  InputIterator end) {
136  return sumPFVector(begin, end, &Candidate::pt, 0.0);
137  }
138 
140 template<typename InputIterator> int sumPFCandCharge(InputIterator begin,
141  InputIterator end) {
142  return sumPFVector(begin, end, &Candidate::charge, 0);
143  }
144 
145 template<typename InputIterator> InputIterator leadCand(InputIterator begin,
146  InputIterator end) {
147  double max_pt = 0;
148  InputIterator max_cand = begin;
149  for(InputIterator cand = begin; cand != end; ++cand) {
150  if( (*cand)->pt() > max_pt ) {
151  max_pt = (*cand)->pt();
152  max_cand = cand;
153  }
154  }
155  return max_cand;
156  }
157 
159 
160 }} // end namespace reco::tau
161 #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:67
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: value.py:1
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:120
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:121
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.