CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

reco::tau::RecoTauPhotonFilter Class Reference

Inheritance diagram for reco::tau::RecoTauPhotonFilter:
reco::tau::RecoTauModifierPlugin reco::tau::RecoTauEventHolderPlugin reco::tau::RecoTauNamedPlugin

List of all members.

Public Member Functions

void operator() (PFTau &) const
 RecoTauPhotonFilter (const edm::ParameterSet &pset)
virtual ~RecoTauPhotonFilter ()

Private Member Functions

bool filter (const RecoTauPiZero *piZero, const reco::Candidate::LorentzVector &tau) const

Private Attributes

double minPtFractionPiZeroes_
double minPtFractionSinglePhotons_

Detailed Description

Definition at line 20 of file RecoTauPhotonFilter.cc.


Constructor & Destructor Documentation

reco::tau::RecoTauPhotonFilter::RecoTauPhotonFilter ( const edm::ParameterSet pset) [explicit]

Definition at line 32 of file RecoTauPhotonFilter.cc.

References edm::ParameterSet::getParameter(), minPtFractionPiZeroes_, and minPtFractionSinglePhotons_.

                                :RecoTauModifierPlugin(pset) {
  minPtFractionSinglePhotons_ =
      pset.getParameter<double>("minPtFractionSinglePhotons");
  minPtFractionPiZeroes_ =
      pset.getParameter<double>("minPtFractionPiZeroes");
}
virtual reco::tau::RecoTauPhotonFilter::~RecoTauPhotonFilter ( ) [inline, virtual]

Definition at line 23 of file RecoTauPhotonFilter.cc.

{}

Member Function Documentation

bool reco::tau::RecoTauPhotonFilter::filter ( const RecoTauPiZero piZero,
const reco::Candidate::LorentzVector tau 
) const [private]

Definition at line 50 of file RecoTauPhotonFilter.cc.

References minPtFractionPiZeroes_, minPtFractionSinglePhotons_, reco::CompositePtrCandidate::numberOfDaughters(), and reco::LeafCandidate::pt().

Referenced by operator()().

                                                   {
  if (piZero->numberOfDaughters() > 1)
    return piZero->pt()/total.pt() < minPtFractionPiZeroes_;
  return piZero->pt()/total.pt() < minPtFractionSinglePhotons_;
}
void reco::tau::RecoTauPhotonFilter::operator() ( PFTau tau) const [virtual]

Implements reco::tau::RecoTauModifierPlugin.

Definition at line 57 of file RecoTauPhotonFilter.cc.

References filterCSVwithJSON::copy, filter(), reco::tau::flattenPiZeros(), reco::PFTau::isolationPFCands(), reco::PFTau::isolationPFGammaCands(), reco::PFTau::isolationPFGammaCandsEtSum(), reco::PFTau::isolationPiZeroCandidates(), edm::Ref< C, T, F >::key(), edm::Ptr< T >::key(), reco::LeafCandidate::p4(), edm::RefVector< C, T, F >::push_back(), reco::PFTau::setisolationPFCands(), reco::PFTau::setisolationPFGammaCands(), reco::PFTau::setisolationPFGammaCandsEtSum(), reco::PFTau::setisolationPiZeroCandidates(), reco::LeafCandidate::setP4(), reco::PFTau::setsignalPFCands(), reco::PFTau::setsignalPiZeroCandidates(), reco::PFTau::signalPFCands(), reco::PFTau::signalPFGammaCands(), reco::PFTau::signalPiZeroCandidates(), and python::multivaluedict::sort().

                                                     {
  std::vector<const RecoTauPiZero*> signalPiZeros;
  BOOST_FOREACH(const RecoTauPiZero& piZero, tau.signalPiZeroCandidates()) {
    signalPiZeros.push_back(&piZero);
  }
  std::sort(signalPiZeros.begin(), signalPiZeros.end(), PiZeroPtSorter());
  std::vector<const RecoTauPiZero*>::const_iterator wimp =
      signalPiZeros.begin();

  // Loop until we have a sturdy enough pizero
  reco::Candidate::LorentzVector totalP4 = tau.p4();
  while(wimp != signalPiZeros.end() && filter(*wimp, totalP4)) {
    totalP4 -= (*wimp)->p4();
    ++wimp;
  }

  if (wimp != signalPiZeros.begin()) {
    // We filtered stuff, update our tau
    double ptDifference = (totalP4 - tau.p4()).pt();
    tau.setisolationPFGammaCandsEtSum(
        tau.isolationPFGammaCandsEtSum() + ptDifference);

    // Update four vector
    tau.setP4(totalP4);

    std::vector<RecoTauPiZero> toMove;
    std::vector<RecoTauPiZero> newSignal;
    std::vector<RecoTauPiZero> newIsolation;

    // Build our new objects
    for (std::vector<const RecoTauPiZero*>::const_iterator iter =
         signalPiZeros.begin(); iter != wimp; ++iter) {
      toMove.push_back(**iter);
    }
    for (std::vector<const RecoTauPiZero*>::const_iterator iter =
         wimp; iter != signalPiZeros.end(); ++iter) {
      newSignal.push_back(**iter);
    }
    // Build our new isolation collection
    std::copy(toMove.begin(), toMove.end(), std::back_inserter(newIsolation));
    std::copy(tau.isolationPiZeroCandidates().begin(),
              tau.isolationPiZeroCandidates().end(),
              std::back_inserter(newIsolation));

    // Set the collections in the taus.
    tau.setsignalPiZeroCandidates(newSignal);
    tau.setisolationPiZeroCandidates(newIsolation);

    // Now we need to deal with the gamma candidates underlying moved pizeros.
    std::vector<PFCandidatePtr> pfcandsToMove = flattenPiZeros(toMove);

    // Copy the keys to move
    std::set<size_t> keysToMove;
    BOOST_FOREACH(const PFCandidatePtr& ptr, pfcandsToMove) {
      keysToMove.insert(ptr.key());
    }

    PFCandidateRefVector newSignalPFGammas;
    PFCandidateRefVector newSignalPFCands;
    PFCandidateRefVector newIsolationPFGammas = tau.isolationPFGammaCands();
    PFCandidateRefVector newIsolationPFCands = tau.isolationPFCands();

    // Move the necessary signal pizeros - what a mess!
    BOOST_FOREACH(const PFCandidateRef& ref, tau.signalPFCands()) {
      if (keysToMove.count(ref.key()))
        newIsolationPFCands.push_back(ref);
      else
        newSignalPFCands.push_back(ref);
    }

    BOOST_FOREACH(const PFCandidateRef& ref, tau.signalPFGammaCands()) {
      if (keysToMove.count(ref.key()))
        newIsolationPFGammas.push_back(ref);
      else
        newSignalPFGammas.push_back(ref);
    }

    tau.setsignalPFCands(newSignalPFCands);
    tau.setsignalPFCands(newSignalPFGammas);
    tau.setisolationPFGammaCands(newIsolationPFGammas);
    tau.setisolationPFCands(newIsolationPFCands);
  }
}

Member Data Documentation

Definition at line 29 of file RecoTauPhotonFilter.cc.

Referenced by filter(), and RecoTauPhotonFilter().

Definition at line 28 of file RecoTauPhotonFilter.cc.

Referenced by filter(), and RecoTauPhotonFilter().