CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/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 "DataFormats/VertexReco/interface/VertexFwd.h"
00017 #include <vector>
00018 #include <algorithm>
00019 #include <numeric>
00020 
00021 // Boost helpers
00022 #include <boost/iterator/transform_iterator.hpp>
00023 #include <boost/iterator/indirect_iterator.hpp>
00024 #include <boost/mem_fn.hpp>
00025 
00026 #include <boost/type_traits/is_base_of.hpp>
00027 #include <boost/static_assert.hpp>
00028 
00029 namespace reco { namespace tau {
00030 
00031 class SortPFCandsDescendingPt {
00032   public:
00033     bool operator()(const PFCandidatePtr& a, const PFCandidatePtr& b) const {
00034       return a->pt() > b->pt();
00035     }
00036 };
00037 
00040 template<typename Iterator> std::vector<PFCandidatePtr>
00041 filterPFCandidates(const Iterator& begin, const Iterator& end,
00042     int particleId, bool sort=true) {
00043   std::vector<PFCandidatePtr> output;
00044   for(Iterator iter = begin; iter != end; ++iter) {
00045     reco::PFCandidatePtr ptr(*iter);
00046     if (ptr->particleId() == particleId)
00047       output.push_back(ptr);
00048   }
00049   if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
00050   return output;
00051 }
00052 
00055 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
00056                                          int particleId, bool sort=true);
00057 
00059 std::vector<PFCandidatePtr> pfCandidates(const PFJet& jet,
00060                                          const std::vector<int>& particleIds,
00061                                          bool sort=true);
00062 
00064 std::vector<PFCandidatePtr> pfChargedCands(const PFJet& jet, bool sort=true);
00065 
00067 std::vector<PFCandidatePtr> pfGammas(const PFJet& jet, bool sort=true);
00068 
00070 std::vector<PFCandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>&);
00071 
00073 template<typename RefVectorType, typename BaseView>
00074 RefVectorType castView(const edm::Handle<BaseView>& view) {
00075   typedef typename RefVectorType::value_type OutputRef;
00076   // Double check at compile time that the inheritance is okay.  It can still
00077   // fail at runtime if you pass it the wrong collection.
00078   BOOST_STATIC_ASSERT(
00079       (boost::is_base_of<typename BaseView::value_type,
00080                          typename RefVectorType::member_type>::value));
00081   RefVectorType output;
00082   size_t nElements = view->size();
00083   output.reserve(nElements);
00084   // Cast each of our Refs
00085   for (size_t i = 0; i < nElements; ++i) {
00086     output.push_back(view->refAt(i).template castTo<OutputRef>());
00087   }
00088   return output;
00089 }
00090 
00091 /*
00092  *Given a range over a container of type C, return a new 'end' iterator such
00093  *that at max <N> elements are taken.  If there are less than N elements in the
00094  *array, leave the <end>  as it is
00095  */
00096 template<typename InputIterator> InputIterator takeNElements(
00097     const InputIterator& begin, const InputIterator& end, size_t N) {
00098   size_t input_size = end - begin;
00099   return (N > input_size) ? end : begin + N;
00100 }
00101 
00103 template<typename InputIterator, typename FunctionPtr, typename ReturnType>
00104   ReturnType sumPFVector(InputIterator begin, InputIterator end,
00105       FunctionPtr func, ReturnType init) {
00106     ReturnType output = init;
00107     for(InputIterator cand = begin; cand != end; ++cand) {
00108       //#define CALL_MEMBER_FN(object,ptrToMember)  ((object).*(ptrToMember))
00109       output += ((**cand).*(func))();
00110     }
00111     return output;
00112   }
00113 
00114 template<typename InputIterator> reco::Candidate::LorentzVector sumPFCandP4(
00115     InputIterator begin, InputIterator end) {
00116   return sumPFVector(begin, end, &PFCandidate::p4,
00117       reco::Candidate::LorentzVector());
00118 }
00119 
00121 template<typename InputIterator> double sumPFCandPt(InputIterator begin,
00122     InputIterator end) {
00123     return sumPFVector(begin, end, &PFCandidate::pt, 0.0);
00124   }
00125 
00127 template<typename InputIterator> int sumPFCandCharge(InputIterator begin,
00128     InputIterator end) {
00129     return sumPFVector(begin, end, &PFCandidate::charge, 0);
00130   }
00131 
00132 template<typename InputIterator> InputIterator leadPFCand(InputIterator begin,
00133     InputIterator end) {
00134     double max_pt = 0;
00135     InputIterator max_cand = begin;
00136     for(InputIterator cand = begin; cand != end; ++cand) {
00137       if( (*cand)->pt() > max_pt ) {
00138         max_pt = (*cand)->pt();
00139         max_cand = cand;
00140       }
00141     }
00142     return max_cand;
00143   }
00144 }}  // end namespace reco::tau
00145 #endif