CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
PFCandidateStripMerger Class Reference

#include <PFCandidateStripMerger.h>

Inheritance diagram for PFCandidateStripMerger:
PFCandidateMergerBase

Public Member Functions

std::vector< std::vector< reco::PFCandidatePtr > > mergeCandidates (const std::vector< reco::PFCandidatePtr > &) override
 
 PFCandidateStripMerger ()
 
 PFCandidateStripMerger (const edm::ParameterSet &)
 
 ~PFCandidateStripMerger () override
 
- Public Member Functions inherited from PFCandidateMergerBase
 PFCandidateMergerBase (const edm::ParameterSet &)
 
 PFCandidateMergerBase ()
 
virtual ~PFCandidateMergerBase ()=0
 

Private Member Functions

bool candidateMatches (const reco::PFCandidatePtr &)
 

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 bacht.nosp@m.is@c.nosp@m.ern.c.nosp@m.h

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_.

10  :
11  PFCandidateMergerBase(config)
12 {
13  inputPdgIds_ = config.getParameter<std::vector<int> >("stripCandidatesPdgIds");
14  etaAssociationDistance_ = config.getParameter<double>("stripEtaAssociationDistance");
15  phiAssociationDistance_ = config.getParameter<double>("stripPhiAssociationDistance");
16 
17 }
T getParameter(std::string const &) const
std::vector< int > inputPdgIds_
PFCandidateStripMerger::~PFCandidateStripMerger ( )
override

Definition at line 20 of file PFCandidateStripMerger.cc.

21 {}

Member Function Documentation

bool PFCandidateStripMerger::candidateMatches ( const reco::PFCandidatePtr cand)
private

Definition at line 25 of file PFCandidateStripMerger.cc.

References funct::abs(), mps_fire::i, inputPdgIds_, patRefSel_triggerMatching_cfi::matches, and reco::LeafCandidate::pdgId().

Referenced by mergeCandidates().

26 {
27  bool matches = false;
28  for(unsigned int i=0; i < inputPdgIds_.size(); ++i) {
29  if(std::abs(cand->pdgId()) == inputPdgIds_.at(i)) {
30  matches = true;
31  continue;
32  }
33  }
34 
35  return matches;
36 }
int pdgId() const final
PDG identifier.
std::vector< int > inputPdgIds_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
vector< vector< PFCandidatePtr > > PFCandidateStripMerger::mergeCandidates ( const std::vector< reco::PFCandidatePtr > &  candidates)
overridevirtual

Implements PFCandidateMergerBase.

Definition at line 43 of file PFCandidateStripMerger.cc.

References candidateMatches(), allConversions_cfi::DeltaPhi, etaAssociationDistance_, mps_fire::i, phiAssociationDistance_, TauTagTools::sortRefVectorByPt(), digi_MixPreMix_cfi::strip, and RecoTauPiZeroBuilderPlugins_cfi::strips.

44 {
45 
46  //Copy the input getting the relevant candidates and sort by pt
47  vector<PFCandidatePtr> cands;
48  for(unsigned int i=0;i<candidates.size();++i)
50  cands.push_back(candidates.at(i));
51 
52  if(cands.size()>1)
54 
55  std::vector<vector<PFCandidatePtr>> strips;
56 
57 
58  //Repeat while there are still unclusterized gammas
59  while(!cands.empty()) {
60 
61  //save the non associated candidates to a different collection
62  vector<PFCandidatePtr> notAssociated;
63 
64  //Create a cluster from the Seed Photon
65  vector<PFCandidatePtr> strip;
66  math::XYZTLorentzVector stripVector;
67  strip.push_back(cands.at(0));
68  stripVector=cands.at(0)->p4();
69 
70  //Loop and associate
71  for(unsigned int i=1;i<cands.size();++i) {
72  if(fabs(cands.at(i)->eta()-stripVector.eta())<etaAssociationDistance_ &&
73  fabs(ROOT::Math::VectorUtil::DeltaPhi(cands.at(i)->p4(),stripVector))<phiAssociationDistance_) {
74  strip.push_back(cands.at(i));
75  stripVector+=cands.at(i)->p4();
76  }
77  else {
78  notAssociated.push_back(cands.at(i));
79  }
80  }
81  //Save the strip
82  strips.push_back(strip);
83 
84  //Swap the candidate vector with the non associated vector
85  cands.swap(notAssociated);
86  //Clear
87  notAssociated.clear();
88  }
89 
90  return strips;
91 }
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool candidateMatches(const reco::PFCandidatePtr &)
void sortRefVectorByPt(std::vector< reco::PFCandidatePtr > &)
Definition: TauTagTools.cc:242

Member Data Documentation

double PFCandidateStripMerger::etaAssociationDistance_
private

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().

double PFCandidateStripMerger::phiAssociationDistance_
private

Definition at line 28 of file PFCandidateStripMerger.h.

Referenced by mergeCandidates(), and PFCandidateStripMerger().