CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateStripMerger.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
3 
4 using namespace reco;
5 using namespace std;
8 {}
9 
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 }
18 
19 
21 {}
22 
23 
24 bool
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 }
37 
38 
39 
40 
41 
42 vector<PFCandidateRefVector>
44 {
45 
46  //Copy the input getting the relevant candidates and sort by pt
48  for(unsigned int i=0;i<candidates.size();++i)
49  if(candidateMatches(candidates.at(i)))
50  cands.push_back(candidates.at(i));
51 
52  if(cands.size()>1)
54 
55  std::vector<PFCandidateRefVector> strips;
56 
57 
58  //Repeat while there are still unclusterized gammas
59  while(cands.size()>0) {
60 
61  //save the non associated candidates to a different collection
62  PFCandidateRefVector notAssociated;
63 
64  //Create a cluster from the Seed Photon
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 }
T getParameter(std::string const &) const
std::vector< reco::PFCandidateRefVector > mergeCandidates(const reco::PFCandidateRefVector &)
int i
Definition: DBlmapReader.cc:9
void strip(std::string &input, const std::string &blanks=" \n\t")
Definition: stringTools.cc:16
void sortRefVectorByPt(reco::PFCandidateRefVector &)
Definition: TauTagTools.cc:242
#define abs(x)
Definition: mlp_lapack.h:159
std::vector< int > inputPdgIds_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
void swap(RefVector< C, T, F > &other) noexcept
Swap two vectors.
Definition: RefVector.h:162
bool candidateMatches(const reco::PFCandidateRef &)
void clear()
Clear the vector.
Definition: RefVector.h:133
value_type const at(size_type idx) const
Retrieve an element of the RefVector.
Definition: RefVector.h:76
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89