CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h

Go to the documentation of this file.
00001 #ifndef RecoTauTag_RecoTau_RecoTauCommonUtilities_h
00002 #define RecoTauTag_RecoTau_RecoTauCommonUtilities_h
00003 
00004 /*
00005  * RecoTauCommonUtilities - utilities for extracting PFCandidates from
00006  * PFJets and summing over collections of PFCandidates.
00007  *
00008  * Author: Evan K. Friis (UC Davis)
00009  *
00010  * $Id $
00011  */
00012 
00013 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00014 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00015 #include "DataFormats/TauReco/interface/PFTau.h"
00016 #include <vector>
00017 #include <algorithm>
00018 #include <numeric>
00019 
00020 // Boost helpers
00021 #include <boost/iterator/transform_iterator.hpp>
00022 #include <boost/iterator/indirect_iterator.hpp>
00023 #include <boost/mem_fn.hpp>
00024 
00025 #include <boost/type_traits/is_base_of.hpp>
00026 #include <boost/static_assert.hpp>
00027 
00028 namespace reco { namespace tau {
00029 
00030 class SortPFCandsDescendingPt {
00031   public:
00032     bool operator()(const PFCandidatePtr& a, const PFCandidatePtr& b) const {
00033       return a->pt() > b->pt();
00034     }
00035 };
00036 
00039 template<typename Iterator> std::vector<PFCandidatePtr>
00040 filterPFCandidates(const Iterator& begin, const Iterator& end,
00041     int particleId, bool sort=true) {
00042   std::vector<PFCandidatePtr> output;
00043   for(Iterator iter = begin; iter != end; ++iter) {
00044     reco::PFCandidatePtr ptr(*iter);
00045     if (ptr->particleId() == particleId)
00046       output.push_back(ptr);
00047   }
00048   if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
00049   return output;
00050 }
00051 
00054 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
00055                                          int particleId, bool sort=true);
00056 
00058 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
00059                                          const std::vector<int>& particleIds,
00060                                          bool sort=true);
00061 
00063 std::vector<PFCandidatePtr> pfChargedCands(const PFJet& jet, bool sort=true);
00064 
00066 std::vector<PFCandidatePtr> pfGammas(const PFJet& jet, bool sort=true);
00067 
00069 std::vector<PFCandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>&);
00070 
00072 template<typename RefVectorType, typename BaseView>
00073 RefVectorType castView(const edm::Handle<BaseView>& view) {
00074   typedef typename RefVectorType::value_type OutputRef;
00075   // Double check at compile time that the inheritance is okay.  It can still
00076   // fail at runtime if you pass it the wrong collection.
00077   BOOST_STATIC_ASSERT(
00078       (boost::is_base_of<typename BaseView::value_type,
00079                          typename RefVectorType::member_type>::value));
00080   RefVectorType output;
00081   size_t nElements = view->size();
00082   output.reserve(nElements);
00083   // Cast each of our Refs
00084   for (size_t i = 0; i < nElements; ++i) {
00085     output.push_back(view->refAt(i).template castTo<OutputRef>());
00086   }
00087   return output;
00088 }
00089 
00090 /*
00091  *Given a range over a container of type C, return a new 'end' iterator such
00092  *that at max <N> elements are taken.  If there are less than N elements in the
00093  *array, leave the <end>  as it is
00094  */
00095 template<typename InputIterator> InputIterator takeNElements(
00096     const InputIterator& begin, const InputIterator& end, size_t N) {
00097   size_t input_size = end - begin;
00098   return (N > input_size) ? end : begin + N;
00099 }
00100 
00102 template<typename InputIterator, typename FunctionPtr, typename ReturnType>
00103   ReturnType sumPFVector(InputIterator begin, InputIterator end,
00104       FunctionPtr func, ReturnType init) {
00105     ReturnType output = init;
00106     for(InputIterator cand = begin; cand != end; ++cand) {
00107       //#define CALL_MEMBER_FN(object,ptrToMember)  ((object).*(ptrToMember))
00108       output += ((**cand).*(func))();
00109     }
00110     return output;
00111   }
00112 
00113 template<typename InputIterator> reco::Candidate::LorentzVector sumPFCandP4(
00114     InputIterator begin, InputIterator end) {
00115   return sumPFVector(begin, end, &PFCandidate::p4,
00116       reco::Candidate::LorentzVector());
00117 }
00118 
00120 template<typename InputIterator> double sumPFCandPt(InputIterator begin,
00121     InputIterator end) {
00122     return sumPFVector(begin, end, &PFCandidate::pt, 0.0);
00123   }
00124 
00126 template<typename InputIterator> int sumPFCandCharge(InputIterator begin,
00127     InputIterator end) {
00128     return sumPFVector(begin, end, &PFCandidate::charge, 0);
00129   }
00130 
00131 template<typename InputIterator> InputIterator leadPFCand(InputIterator begin,
00132     InputIterator end) {
00133     double max_pt = 0;
00134     InputIterator max_cand = begin;
00135     for(InputIterator cand = begin; cand != end; ++cand) {
00136       if( (*cand)->pt() > max_pt ) {
00137         max_pt = (*cand)->pt();
00138         max_cand = cand;
00139       }
00140     }
00141     return max_cand;
00142   }
00143 }}  // end namespace reco::tau
00144 #endif