CMS 3D CMS Logo

Namespaces | Classes | Functions

reco::tau Namespace Reference

Namespaces

namespace  cone
namespace  disc
namespace  helpers
namespace  qcuts

Classes

class  AssociationMatchRefSelector
class  Combinatoric
class  CombinatoricGenerator
class  CombinatoricIterator
class  RecoTauBuilderCombinatoricPlugin
class  RecoTauBuilderConePlugin
class  RecoTauBuilderPlugin
class  RecoTauCleanerPlugin
class  RecoTauConstructor
class  RecoTauDiscriminantCleanerPlugin
class  RecoTauDiscriminantFunctionPlugin
class  RecoTauDiscriminantPlugin
class  RecoTauDiscriminantVectorFunctionPlugin
class  RecoTauElectronRejectionPlugin
class  RecoTauEventHolderPlugin
class  RecoTauImpactParameterSignificancePlugin
class  RecoTauLexicographicalRanking
class  RecoTauModifierPlugin
class  RecoTauMVAHelper
class  RecoTauNamedPlugin
class  RecoTauPhotonFilter
class  RecoTauPiZeroBuilderPlugin
class  RecoTauPiZeroCombinatoricPlugin
class  RecoTauPiZeroQualityPlugin
class  RecoTauPiZeroStringQuality
class  RecoTauPiZeroStripPlugin
class  RecoTauPiZeroTrivialPlugin
class  RecoTauQualityCuts
class  RecoTauRandomCleanerPlugin
class  RecoTauStringCleanerPlugin
class  RecoTauTagInfoWorkaroundModifer
class  RecoTauTruthEmbedder
class  RecoTauTwoProngFilter
class  SortByDescendingPt
class  SortPFCandsDescendingPt

Functions

template<typename RefVectorType , typename BaseView >
RefVectorType castView (const edm::Handle< BaseView > &view)
 Convert a BaseView (View<T>) to a TRefVector.
unsigned int chargedHadronsInDecayMode (PFTau::hadronicDecayMode mode)
 Reverse mapping of decay modes into multiplicities.
template<typename Container , class OverlapFunction >
Container cleanOverlaps (const Container &dirty)
std::string discPluginName (const std::string &mvaName)
std::vector< PFCandidatePtrflattenPiZeros (const std::vector< RecoTauPiZero > &)
 Flatten a list of pi zeros into a list of there constituent PFCandidates.
template<typename InputIterator >
InputIterator leadPFCand (InputIterator begin, InputIterator end)
std::vector< PFCandidatePtrpfCandidates (const PFJet &jet, int particleId, bool sort=true)
std::vector< PFCandidatePtrpfCandidates (const PFJet &jet, const std::vector< int > &particleIds, bool sort=true)
 Extract pfCandidates of a that match a list of particle Ids from a PFJet.
std::vector< PFCandidatePtrpfChargedCands (const PFJet &jet, bool sort=true)
 Extract all non-neutral candidates from a PFJet.
std::vector< PFCandidatePtrpfGammas (const PFJet &jet, bool sort=true)
 Extract all pfGammas from a PFJet.
unsigned int piZerosInDecayMode (PFTau::hadronicDecayMode mode)
template<typename InputIterator >
int sumPFCandCharge (InputIterator begin, InputIterator end)
 Sum the PT of a collection of PFCandidates.
template<typename InputIterator >
reco::Candidate::LorentzVector sumPFCandP4 (InputIterator begin, InputIterator end)
template<typename InputIterator >
double sumPFCandPt (InputIterator begin, InputIterator end)
 Sum the PT of a collection of PFCandidates.
template<typename InputIterator , typename FunctionPtr , typename ReturnType >
ReturnType sumPFVector (InputIterator begin, InputIterator end, FunctionPtr func, ReturnType init)
 Sum the four vectors in a collection of PFCandidates.
template<typename InputIterator >
InputIterator takeNElements (const InputIterator &begin, const InputIterator &end, size_t N)
PFTau::hadronicDecayMode translateDecayMode (unsigned int nCharged, unsigned int nPiZero)

Function Documentation

template<typename RefVectorType , typename BaseView >
RefVectorType reco::tau::castView ( const edm::Handle< BaseView > &  view)

Convert a BaseView (View<T>) to a TRefVector.

Definition at line 51 of file RecoTauCommonUtilities.h.

References i, convertSQLitetoXML_cfg::output, and relativeConstraints::value.

                                                        {
  typedef typename RefVectorType::value_type OutputRef;
  // Double check at compile time that the inheritance is okay.  It can still
  // fail at runtime if you pass it the wrong collection.
  BOOST_STATIC_ASSERT(
      (boost::is_base_of<typename BaseView::value_type,
                         typename RefVectorType::member_type>::value));
  RefVectorType output;
  size_t nElements = view->size();
  output.reserve(nElements);
  // Cast each of our Refs
  for (size_t i = 0; i < nElements; ++i) {
    output.push_back(view->refAt(i).template castTo<OutputRef>());
  }
  return output;
}
unsigned int reco::tau::chargedHadronsInDecayMode ( PFTau::hadronicDecayMode  mode)

Reverse mapping of decay modes into multiplicities.

Definition at line 5 of file PFTauDecayModeTools.cc.

References reco::PFTau::kOneProngNPiZero, and mode.

                                                                    {
   int modeAsInt = static_cast<int>(mode);
   return (modeAsInt / PFTau::kOneProngNPiZero) + 1;
}
template<typename Container , class OverlapFunction >
Container reco::tau::cleanOverlaps ( const Container &  dirty)

Definition at line 34 of file RecoTauCleaningTools.h.

References clean, and analyzePatCleaning_cfg::overlaps.

Referenced by RecoTauCleanerImpl< Prod >::produce().

                                                {
  typedef typename Container::const_iterator Iterator;
  // Output container of clean objects
  Container clean;
  OverlapFunction overlapChecker;
  for (Iterator candidate = dirty.begin(); candidate != dirty.end();
       ++candidate) {
    // Check if this overlaps with a pizero already in the clean list
    bool overlaps = false;
    for (Iterator cleaned = clean.begin();
         cleaned != clean.end() && !overlaps; ++cleaned) {
      overlaps = overlapChecker(*candidate, *cleaned);
    }
    // If it didn't overlap with anything clean, add it to the clean list
    if (!overlaps)
      clean.insert(clean.end(), *candidate);
  }
  return clean;
}
std::string reco::tau::discPluginName ( const std::string &  mvaName) [inline]
std::vector< PFCandidatePtr > reco::tau::flattenPiZeros ( const std::vector< RecoTauPiZero > &  piZeros)

Flatten a list of pi zeros into a list of there constituent PFCandidates.

Definition at line 21 of file RecoTauCommonUtilities.cc.

References convertSQLitetoXML_cfg::output.

Referenced by reco::tau::RecoTauPhotonFilter::operator()().

                                                                                  {
  std::vector<PFCandidatePtr> output;

  for(std::vector<RecoTauPiZero>::const_iterator piZero = piZeros.begin();
      piZero != piZeros.end(); ++piZero)
  {
    for(size_t iDaughter = 0; iDaughter < piZero->numberOfDaughters(); ++iDaughter)
    {
      output.push_back(PFCandidatePtr(piZero->daughterPtr(iDaughter)));
    }
  }
  return output;
}
template<typename InputIterator >
InputIterator reco::tau::leadPFCand ( InputIterator  begin,
InputIterator  end 
)

Definition at line 109 of file RecoTauCommonUtilities.h.

References begin, and end.

Referenced by reco::tau::RecoTauConstructor::get(), and reco::tau::RecoTauBuilderConePlugin::operator()().

                       {
    double max_pt = 0;
    InputIterator max_cand = begin;
    for(InputIterator cand = begin; cand != end; ++cand) {
      if( (*cand)->pt() > max_pt ) {
        max_pt = (*cand)->pt();
        max_cand = cand;
      }
    }
    return max_cand;
  }
std::vector< reco::PFCandidatePtr > reco::tau::pfCandidates ( const PFJet &  jet,
int  particleId,
bool  sort = true 
)

Extract pfCandidates of a given particle Id from a PFJet. If sort is true, candidates will be sorted by descending PT

Definition at line 35 of file RecoTauCommonUtilities.cc.

References reco::PFJet::getPFConstituents(), and python::multivaluedict::sort().

Referenced by PFMETAnalyzer::analyze(), reco::tau::RecoTauBuilderCombinatoricPlugin::operator()(), reco::tau::RecoTauBuilderConePlugin::operator()(), pfCandidates(), pfChargedCands(), pfGammas(), PFElectronTranslator::produce(), PFPileUp::produce(), pat::PATPFParticleProducer::produce(), and PFMETAnalyzer::validateMET().

{
  PFCandPtrs pfCands = jet.getPFConstituents();
  PFCandPtrs selectedPFCands;

  for(PFCandIter cand = pfCands.begin(); cand != pfCands.end(); ++cand)
  {
    if( (**cand).particleId() == particleId )
      selectedPFCands.push_back(*cand);
  }

  if ( sort ) std::sort(selectedPFCands.begin(), selectedPFCands.end(), SortPFCandsDescendingPt());

  return selectedPFCands;
}
std::vector< reco::PFCandidatePtr > reco::tau::pfCandidates ( const PFJet &  jet,
const std::vector< int > &  particleIds,
bool  sort = true 
)

Extract pfCandidates of a that match a list of particle Ids from a PFJet.

Definition at line 51 of file RecoTauCommonUtilities.cc.

References convertSQLitetoXML_cfg::output, pfCandidates(), and python::multivaluedict::sort().

{
  PFCandPtrs output;
  // Get each desired candidate type, unsorted for now
  for(std::vector<int>::const_iterator particleId = particleIds.begin();
      particleId != particleIds.end(); ++particleId) {
    PFCandPtrs selectedPFCands = pfCandidates(jet, *particleId, false);
    output.insert(output.end(), selectedPFCands.begin(), selectedPFCands.end());
  }
  if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
  return output;
}
std::vector< reco::PFCandidatePtr > reco::tau::pfChargedCands ( const PFJet &  jet,
bool  sort = true 
)

Extract all non-neutral candidates from a PFJet.

Definition at line 68 of file RecoTauCommonUtilities.cc.

References reco::PFCandidate::e, reco::PFCandidate::h, reco::PFCandidate::mu, convertSQLitetoXML_cfg::output, pfCandidates(), and python::multivaluedict::sort().

Referenced by reco::tau::RecoTauBuilderCombinatoricPlugin::operator()(), and reco::tau::RecoTauBuilderConePlugin::operator()().

                                                            {
  PFCandPtrs output;
  PFCandPtrs mus = pfCandidates(jet, reco::PFCandidate::mu, false);
  PFCandPtrs es = pfCandidates(jet, reco::PFCandidate::e, false);
  PFCandPtrs chs = pfCandidates(jet, reco::PFCandidate::h, false);
  output.reserve(mus.size() + es.size() + chs.size());
  output.insert(output.end(), mus.begin(), mus.end());
  output.insert(output.end(), es.begin(), es.end());
  output.insert(output.end(), chs.begin(), chs.end());
  if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
  return output;
}
std::vector< reco::PFCandidatePtr > reco::tau::pfGammas ( const PFJet &  jet,
bool  sort = true 
)
unsigned int reco::tau::piZerosInDecayMode ( PFTau::hadronicDecayMode  mode)

Definition at line 10 of file PFTauDecayModeTools.cc.

References reco::PFTau::kOneProngNPiZero, and mode.

                                                             {
   int modeAsInt = static_cast<int>(mode);
   return (modeAsInt % PFTau::kOneProngNPiZero);
}
template<typename InputIterator >
int reco::tau::sumPFCandCharge ( InputIterator  begin,
InputIterator  end 
)

Sum the PT of a collection of PFCandidates.

Definition at line 104 of file RecoTauCommonUtilities.h.

References reco::LeafCandidate::charge(), and sumPFVector().

Referenced by reco::tau::RecoTauConstructor::get().

template<typename InputIterator >
reco::Candidate::LorentzVector reco::tau::sumPFCandP4 ( InputIterator  begin,
InputIterator  end 
)
template<typename InputIterator >
double reco::tau::sumPFCandPt ( InputIterator  begin,
InputIterator  end 
)

Sum the PT of a collection of PFCandidates.

Definition at line 98 of file RecoTauCommonUtilities.h.

References reco::LeafCandidate::pt(), and sumPFVector().

Referenced by reco::tau::RecoTauConstructor::get().

                       {
    return sumPFVector(begin, end, &PFCandidate::pt, 0.0);
  }
template<typename InputIterator , typename FunctionPtr , typename ReturnType >
ReturnType reco::tau::sumPFVector ( InputIterator  begin,
InputIterator  end,
FunctionPtr  func,
ReturnType  init 
)

Sum the four vectors in a collection of PFCandidates.

Definition at line 81 of file RecoTauCommonUtilities.h.

References end, init, and convertSQLitetoXML_cfg::output.

Referenced by sumPFCandCharge(), sumPFCandP4(), and sumPFCandPt().

                                         {
    ReturnType output = init;
    for(InputIterator cand = begin; cand != end; ++cand) {
      //#define CALL_MEMBER_FN(object,ptrToMember)  ((object).*(ptrToMember))
      output += ((**cand).*(func))();
    }
    return output;
  }
template<typename InputIterator >
InputIterator reco::tau::takeNElements ( const InputIterator &  begin,
const InputIterator &  end,
size_t  N 
)

Definition at line 73 of file RecoTauCommonUtilities.h.

References begin, and MultiGaussianStateTransform::N.

Referenced by reco::tau::RecoTauBuilderCombinatoricPlugin::operator()(), and reco::tau::RecoTauPiZeroCombinatoricPlugin::operator()().

                                                                    {
  size_t input_size = end - begin;
  return (N > input_size) ? end : begin + N;
}
PFTau::hadronicDecayMode reco::tau::translateDecayMode ( unsigned int  nCharged,
unsigned int  nPiZero 
)

Definition at line 15 of file PFTauDecayModeTools.cc.

References kNull, reco::PFTau::kOneProngNPiZero, and reco::PFTau::kRareDecayMode.

Referenced by reco::tau::RecoTauTruthEmbedder::operator()(), RecoTauMVADiscriminator::RecoTauMVADiscriminator(), and RecoTauMVATransform::RecoTauMVATransform().

                                                  {
   // If no tracks exist, this is definitely not a tau!
   if(!nCharged) return PFTau::kNull;
   // Find the maximum number of PiZeros our parameterization can hold
   const unsigned int maxPiZeros = PFTau::kOneProngNPiZero;
   // Determine our track index
   unsigned int trackIndex = (nCharged-1)*(maxPiZeros+1);
   // Check if we handle the given number of tracks
   if(trackIndex >= PFTau::kRareDecayMode) return PFTau::kRareDecayMode;

   nPiZeros = (nPiZeros <= maxPiZeros) ? nPiZeros : maxPiZeros;
   return static_cast<PFTau::hadronicDecayMode>(trackIndex + nPiZeros);
}