CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
metsig::SignPFSpecificAlgo Class Reference

#include <SignPFSpecificAlgo.h>

Public Member Functions

void addPFCandidate (reco::PFCandidatePtr pf)
 
void addPFJets (const edm::View< reco::PFJet > *PFJets)
 
reco::METCovMatrix getSignifMatrix () const
 
reco::METCovMatrix mkSignifMatrix (edm::Handle< edm::View< reco::Candidate > > &PFCandidates)
 
void setResolutions (metsig::SignAlgoResolutions *resolutions)
 
 SignPFSpecificAlgo ()
 
void useOriginalPtrs (const edm::ProductID &productID)
 
 ~SignPFSpecificAlgo ()
 

Private Attributes

metsig::significanceAlgo algo_
 
std::set< reco::CandidatePtrclusteredParticlePtrs_
 
metsig::SignAlgoResolutionsresolutions_
 

Detailed Description

Definition at line 30 of file SignPFSpecificAlgo.h.

Constructor & Destructor Documentation

◆ SignPFSpecificAlgo()

SignPFSpecificAlgo::SignPFSpecificAlgo ( )

Definition at line 17 of file SignPFSpecificAlgo.cc.

References clusteredParticlePtrs_.

17 : resolutions_(nullptr), algo_() { clusteredParticlePtrs_.clear(); }
metsig::significanceAlgo algo_
metsig::SignAlgoResolutions * resolutions_
std::set< reco::CandidatePtr > clusteredParticlePtrs_

◆ ~SignPFSpecificAlgo()

metsig::SignPFSpecificAlgo::~SignPFSpecificAlgo ( )
inline

Definition at line 33 of file SignPFSpecificAlgo.h.

33 {}

Member Function Documentation

◆ addPFCandidate()

void SignPFSpecificAlgo::addPFCandidate ( reco::PFCandidatePtr  pf)

Definition at line 60 of file SignPFSpecificAlgo.cc.

References packedPFCandidateRefMixer_cfi::pf.

60  {
62  return; //pf candidate already added in jet collection
63  }
64  std::vector<metsig::SigInputObj> vobj;
65  vobj.push_back(resolutions_->evalPF(&(*pf)));
66  algo_.addObjects(vobj);
67 }
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
metsig::SigInputObj evalPF(const reco::PFCandidate *candidate) const
metsig::significanceAlgo algo_
metsig::SignAlgoResolutions * resolutions_
std::set< reco::CandidatePtr > clusteredParticlePtrs_

◆ addPFJets()

void SignPFSpecificAlgo::addPFJets ( const edm::View< reco::PFJet > *  PFJets)

Definition at line 25 of file SignPFSpecificAlgo.cc.

References metsig::jet, and JetHTJetPlusHOFilter_cff::PFJets.

25  {
26  std::vector<metsig::SigInputObj> vobj;
27  for (edm::View<reco::PFJet>::const_iterator jet = PFJets->begin(); jet != PFJets->end(); ++jet) {
28  vobj.push_back(resolutions_->evalPFJet(&(*jet)));
29  std::vector<reco::PFCandidatePtr> pfs = jet->getPFConstituents();
30  for (std::vector<reco::PFCandidatePtr>::const_iterator it = pfs.begin(); it != pfs.end(); ++it) {
31  reco::CandidatePtr ptr(*it);
32  clusteredParticlePtrs_.insert(ptr);
33  }
34  }
35  algo_.addObjects(vobj);
36 }
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
metsig::SigInputObj evalPFJet(const reco::Jet *jet) const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
metsig::significanceAlgo algo_
metsig::SignAlgoResolutions * resolutions_
std::set< reco::CandidatePtr > clusteredParticlePtrs_

◆ getSignifMatrix()

reco::METCovMatrix metsig::SignPFSpecificAlgo::getSignifMatrix ( ) const
inline

Definition at line 39 of file SignPFSpecificAlgo.h.

References algo_, and metsig::significanceAlgo::getSignifMatrix().

39 { return algo_.getSignifMatrix(); }
metsig::significanceAlgo algo_
reco::METCovMatrix getSignifMatrix() const

◆ mkSignifMatrix()

reco::METCovMatrix SignPFSpecificAlgo::mkSignifMatrix ( edm::Handle< edm::View< reco::Candidate > > &  PFCandidates)

Definition at line 70 of file SignPFSpecificAlgo.cc.

References mps_fire::end, packedPFCandidateRefMixer_cfi::pf, and pfPileUpJME_cfi::PFCandidates.

70  {
72  for (edm::View<reco::Candidate>::const_iterator iParticle = (PFCandidates.product())->begin();
73  iParticle != (PFCandidates.product())->end();
74  ++iParticle) {
75  const reco::PFCandidate* pfCandidate = dynamic_cast<const reco::PFCandidate*>(&(*iParticle));
76  if (!pfCandidate)
77  continue;
78  reco::CandidatePtr dau(PFCandidates, iParticle - PFCandidates->begin());
79  if (dau.isNull())
80  continue;
81  if (!dau.isAvailable())
82  continue;
83  reco::PFCandidatePtr pf(dau.id(), pfCandidate, dau.key());
85  }
86  return getSignifMatrix();
87 }
reco::METCovMatrix getSignifMatrix() const
void addPFCandidate(reco::PFCandidatePtr pf)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
void useOriginalPtrs(const edm::ProductID &productID)

◆ setResolutions()

void SignPFSpecificAlgo::setResolutions ( metsig::SignAlgoResolutions resolutions)

Definition at line 20 of file SignPFSpecificAlgo.cc.

References electronProducer_cfi::resolutions.

20  {
22 }
metsig::SignAlgoResolutions * resolutions_

◆ useOriginalPtrs()

void SignPFSpecificAlgo::useOriginalPtrs ( const edm::ProductID productID)

Definition at line 39 of file SignPFSpecificAlgo.cc.

References edm::Ptr< T >::id(), edm::Ptr< T >::isNull(), and groupFilesInBlocks::temp.

39  {
40  std::set<reco::CandidatePtr>::const_iterator it = clusteredParticlePtrs_.begin();
41  reco::CandidatePtr ptr(*it);
42  if (ptr.id() == productID)
43  return; //If the first element is from the right product, return
44 
45  std::set<reco::CandidatePtr> temp;
46  for (; it != clusteredParticlePtrs_.end(); ++it) {
47  reco::CandidatePtr ptr(*it);
48  while (ptr.id() != productID) {
49  ptr = ptr->sourceCandidatePtr(0);
50  if (ptr.isNull())
51  return; //if it does not get to the correct product, return
52  }
53  temp.insert(ptr);
54  }
55  clusteredParticlePtrs_.clear();
57 }
std::set< reco::CandidatePtr > clusteredParticlePtrs_

Member Data Documentation

◆ algo_

metsig::significanceAlgo metsig::SignPFSpecificAlgo::algo_
private

Definition at line 45 of file SignPFSpecificAlgo.h.

Referenced by getSignifMatrix().

◆ clusteredParticlePtrs_

std::set<reco::CandidatePtr> metsig::SignPFSpecificAlgo::clusteredParticlePtrs_
private

Definition at line 44 of file SignPFSpecificAlgo.h.

Referenced by SignPFSpecificAlgo().

◆ resolutions_

metsig::SignAlgoResolutions* metsig::SignPFSpecificAlgo::resolutions_
private

Definition at line 43 of file SignPFSpecificAlgo.h.