CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PFCandidateStripMerger Class Reference

#include <PFCandidateStripMerger.h>

Inheritance diagram for PFCandidateStripMerger:
PFCandidateMergerBase

List of all members.

Public Member Functions

std::vector
< reco::PFCandidateRefVector
mergeCandidates (const reco::PFCandidateRefVector &)
 PFCandidateStripMerger (const edm::ParameterSet &)
 PFCandidateStripMerger ()
 ~PFCandidateStripMerger ()

Private Member Functions

bool candidateMatches (const reco::PFCandidateRef &)

Private Attributes

double etaAssociationDistance_
std::vector< int > inputPdgIds_
double phiAssociationDistance_

Detailed Description

PFCandidateStripMerger

Class that creates strips from Particle Flow Candidates And outputs a Collection of Candidate Lists Michail Bachtis ---------------------- University of Wisconsin bachtis@cern.ch

Definition at line 15 of file PFCandidateStripMerger.h.


Constructor & Destructor Documentation

PFCandidateStripMerger::PFCandidateStripMerger ( )

Definition at line 6 of file PFCandidateStripMerger.cc.

PFCandidateStripMerger::PFCandidateStripMerger ( const edm::ParameterSet config)

Definition at line 10 of file PFCandidateStripMerger.cc.

References etaAssociationDistance_, edm::ParameterSet::getParameter(), inputPdgIds_, and phiAssociationDistance_.

                                                                           :
  PFCandidateMergerBase(config)
{
  inputPdgIds_ = config.getParameter<std::vector<int> >("stripCandidatesPdgIds");
  etaAssociationDistance_ = config.getParameter<double>("stripEtaAssociationDistance");
  phiAssociationDistance_ = config.getParameter<double>("stripPhiAssociationDistance");

}
PFCandidateStripMerger::~PFCandidateStripMerger ( )

Definition at line 20 of file PFCandidateStripMerger.cc.

{}

Member Function Documentation

bool PFCandidateStripMerger::candidateMatches ( const reco::PFCandidateRef cand) [private]

Definition at line 25 of file PFCandidateStripMerger.cc.

References abs, i, and inputPdgIds_.

Referenced by mergeCandidates().

{
  bool matches = false;
  for(unsigned int i=0; i < inputPdgIds_.size(); ++i) {
    if(std::abs(cand->pdgId()) == inputPdgIds_.at(i)) {
      matches = true;
      continue;
    }
  }

  return matches;
}
vector< PFCandidateRefVector > PFCandidateStripMerger::mergeCandidates ( const reco::PFCandidateRefVector candidates) [virtual]

Implements PFCandidateMergerBase.

Definition at line 43 of file PFCandidateStripMerger.cc.

References edm::RefVector< C, T, F >::at(), candidateMatches(), edm::RefVector< C, T, F >::clear(), etaAssociationDistance_, i, phiAssociationDistance_, edm::RefVector< C, T, F >::push_back(), edm::RefVector< C, T, F >::size(), TauTagTools::sortRefVectorByPt(), RecoTauPiZeroBuilderPlugins_cfi::strips, and edm::RefVector< C, T, F >::swap().

{

  //Copy the input getting the relevant candidates and sort by pt 
  PFCandidateRefVector cands;
  for(unsigned int i=0;i<candidates.size();++i)
    if(candidateMatches(candidates.at(i)))
      cands.push_back(candidates.at(i));
 
    if(cands.size()>1)
    TauTagTools::sortRefVectorByPt(cands);

    std::vector<PFCandidateRefVector> strips;


  //Repeat while there are still unclusterized gammas
  while(cands.size()>0) {

    //save the non associated candidates to a different collection
    PFCandidateRefVector notAssociated;
    
    //Create a cluster from the Seed Photon
    PFCandidateRefVector strip;
    math::XYZTLorentzVector stripVector;
    strip.push_back(cands.at(0));
    stripVector=cands.at(0)->p4();

    //Loop and associate
    for(unsigned int i=1;i<cands.size();++i) {
        if(fabs(cands.at(i)->eta()-stripVector.eta())<etaAssociationDistance_ &&  
           fabs(ROOT::Math::VectorUtil::DeltaPhi(cands.at(i)->p4(),stripVector))<phiAssociationDistance_) {
            strip.push_back(cands.at(i));
            stripVector+=cands.at(i)->p4();
          }
          else {
            notAssociated.push_back(cands.at(i));
          }
    }
    //Save the strip 
    strips.push_back(strip);

    //Swap the candidate vector with the non associated vector
    cands.swap(notAssociated);
    //Clear 
    notAssociated.clear();
  }

  return strips;
}

Member Data Documentation

Definition at line 27 of file PFCandidateStripMerger.h.

Referenced by mergeCandidates(), and PFCandidateStripMerger().

std::vector<int> PFCandidateStripMerger::inputPdgIds_ [private]

Definition at line 26 of file PFCandidateStripMerger.h.

Referenced by candidateMatches(), and PFCandidateStripMerger().

Definition at line 28 of file PFCandidateStripMerger.h.

Referenced by mergeCandidates(), and PFCandidateStripMerger().