CMS 3D CMS Logo

SignPFSpecificAlgo.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: METAlgorithms
4 // Class: SignPFSpecificAlgo
5 //
6 // Authors: A. Khukhunaishvili (Cornell), L. Gibbons (Cornell)
7 // First Implementation: November 11, 2011
8 //
9 //
10 
11 //____________________________________________________________________________||
13 
15 
16 //____________________________________________________________________________||
18 : resolutions_(nullptr), algo_()
19 {
20  clusteredParticlePtrs_.clear();
21 }
22 
23 //____________________________________________________________________________||
25 {
27 }
28 
29 //____________________________________________________________________________||
31 {
32  std::vector<metsig::SigInputObj> vobj;
33  for(edm::View<reco::PFJet>::const_iterator jet = PFJets->begin(); jet != PFJets->end(); ++jet)
34  {
35  vobj.push_back(resolutions_->evalPFJet(&(*jet)));
36  std::vector<reco::PFCandidatePtr> pfs = jet->getPFConstituents();
37  for(std::vector<reco::PFCandidatePtr>::const_iterator it=pfs.begin(); it!=pfs.end(); ++it)
38  {
39  reco::CandidatePtr ptr(*it);
40  clusteredParticlePtrs_.insert(ptr);
41  }
42  }
43  algo_.addObjects(vobj);
44 }
45 
46 //____________________________________________________________________________||
48 {
49  std::set<reco::CandidatePtr>::const_iterator it=clusteredParticlePtrs_.begin();
50  reco::CandidatePtr ptr(*it);
51  if(ptr.id()==productID) return; //If the first element is from the right product, return
52 
53  std::set<reco::CandidatePtr> temp;
54  for(; it!=clusteredParticlePtrs_.end(); ++it)
55  {
56  reco::CandidatePtr ptr(*it);
57  while(ptr.id()!=productID)
58  {
59  ptr = ptr->sourceCandidatePtr(0);
60  if(ptr.isNull()) return; //if it does not get to the correct product, return
61  }
62  temp.insert(ptr);
63  }
64  clusteredParticlePtrs_.clear();
66 }
67 
68 //____________________________________________________________________________||
70 {
72  {
73  return; //pf candidate already added in jet collection
74  }
75  std::vector<metsig::SigInputObj> vobj;
76  vobj.push_back(resolutions_->evalPF(&(*pf)));
77  algo_.addObjects(vobj);
78 }
79 
80 //____________________________________________________________________________||
82 {
84  for(edm::View<reco::Candidate>::const_iterator iParticle = (PFCandidates.product())->begin(); iParticle != (PFCandidates.product())->end(); ++iParticle )
85  {
86  const reco::PFCandidate* pfCandidate = dynamic_cast<const reco::PFCandidate*> (&(*iParticle));
87  if (!pfCandidate) continue;
88  reco::CandidatePtr dau(PFCandidates, iParticle - PFCandidates->begin());
89  if(dau.isNull()) continue;
90  if(!dau.isAvailable()) continue;
91  reco::PFCandidatePtr pf(dau.id(), pfCandidate, dau.key());
93  }
94  return getSignifMatrix();
95 }
96 
97 //____________________________________________________________________________||
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
metsig::SigInputObj evalPF(const reco::PFCandidate *candidate) const
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:40
#define nullptr
reco::METCovMatrix getSignifMatrix() const
const_iterator begin() const
void addPFJets(const edm::View< reco::PFJet > *PFJets)
void setResolutions(metsig::SignAlgoResolutions *resolutions)
bool isNull() const
Checks for null.
Definition: Ptr.h:164
#define end
Definition: vmac.h:37
reco::METCovMatrix mkSignifMatrix(edm::Handle< edm::View< reco::Candidate > > &PFCandidates)
ProductID id() const
Accessor for product ID.
Definition: Ptr.h:180
void addPFCandidate(reco::PFCandidatePtr pf)
metsig::SigInputObj evalPFJet(const reco::PFJet *jet) const
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
#define begin
Definition: vmac.h:30
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
const_iterator end() const
void useOriginalPtrs(const edm::ProductID &productID)
metsig::significanceAlgo algo_
metsig::SignAlgoResolutions * resolutions_
std::set< reco::CandidatePtr > clusteredParticlePtrs_