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 31 of file SignPFSpecificAlgo.h.

Constructor & Destructor Documentation

SignPFSpecificAlgo::SignPFSpecificAlgo ( )

Definition at line 17 of file SignPFSpecificAlgo.cc.

References clusteredParticlePtrs_.

18 : resolutions_(nullptr), algo_()
19 {
20  clusteredParticlePtrs_.clear();
21 }
metsig::significanceAlgo algo_
metsig::SignAlgoResolutions * resolutions_
std::set< reco::CandidatePtr > clusteredParticlePtrs_
metsig::SignPFSpecificAlgo::~SignPFSpecificAlgo ( )
inline

Member Function Documentation

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

Definition at line 69 of file SignPFSpecificAlgo.cc.

References metsig::significanceAlgo::addObjects(), algo_, clusteredParticlePtrs_, metsig::SignAlgoResolutions::evalPF(), and resolutions_.

Referenced by mkSignifMatrix(), and ~SignPFSpecificAlgo().

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 }
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_
void SignPFSpecificAlgo::addPFJets ( const edm::View< reco::PFJet > *  PFJets)

Definition at line 30 of file SignPFSpecificAlgo.cc.

References metsig::significanceAlgo::addObjects(), algo_, edm::View< T >::begin(), clusteredParticlePtrs_, edm::View< T >::end(), metsig::SignAlgoResolutions::evalPFJet(), metsig::jet, and resolutions_.

Referenced by ~SignPFSpecificAlgo().

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 }
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
const_iterator begin() const
metsig::SigInputObj evalPFJet(const reco::PFJet *jet) const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
const_iterator end() const
metsig::significanceAlgo algo_
metsig::SignAlgoResolutions * resolutions_
std::set< reco::CandidatePtr > clusteredParticlePtrs_
reco::METCovMatrix metsig::SignPFSpecificAlgo::getSignifMatrix ( ) const
inline

Definition at line 42 of file SignPFSpecificAlgo.h.

References algo_, metsig::significanceAlgo::getSignifMatrix(), mkSignifMatrix(), and valueMapTest_cfi::PFCandidates.

Referenced by mkSignifMatrix().

42 {return algo_.getSignifMatrix();}
reco::METCovMatrix getSignifMatrix() const
metsig::significanceAlgo algo_
reco::METCovMatrix SignPFSpecificAlgo::mkSignifMatrix ( edm::Handle< edm::View< reco::Candidate > > &  PFCandidates)

Definition at line 81 of file SignPFSpecificAlgo.cc.

References addPFCandidate(), begin, end, getSignifMatrix(), packedPFCandidateRefMixer_cfi::pf, valueMapTest_cfi::PFCandidates, and useOriginalPtrs().

Referenced by getSignifMatrix().

82 {
83  useOriginalPtrs(PFCandidates.id());
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 }
ProductID id() const
Definition: HandleBase.cc:15
reco::METCovMatrix getSignifMatrix() const
const_iterator begin() const
#define end
Definition: vmac.h:37
T const * product() const
Definition: Handle.h:81
void addPFCandidate(reco::PFCandidatePtr pf)
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
void useOriginalPtrs(const edm::ProductID &productID)
void SignPFSpecificAlgo::setResolutions ( metsig::SignAlgoResolutions resolutions)

Definition at line 24 of file SignPFSpecificAlgo.cc.

References electronProducer_cfi::resolutions, and resolutions_.

Referenced by ~SignPFSpecificAlgo().

25 {
27 }
metsig::SignAlgoResolutions * resolutions_
void SignPFSpecificAlgo::useOriginalPtrs ( const edm::ProductID productID)

Definition at line 47 of file SignPFSpecificAlgo.cc.

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

Referenced by mkSignifMatrix(), and ~SignPFSpecificAlgo().

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 }
std::set< reco::CandidatePtr > clusteredParticlePtrs_

Member Data Documentation

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

Definition at line 48 of file SignPFSpecificAlgo.h.

Referenced by addPFCandidate(), addPFJets(), and getSignifMatrix().

std::set<reco::CandidatePtr> metsig::SignPFSpecificAlgo::clusteredParticlePtrs_
private
metsig::SignAlgoResolutions* metsig::SignPFSpecificAlgo::resolutions_
private

Definition at line 46 of file SignPFSpecificAlgo.h.

Referenced by addPFCandidate(), addPFJets(), and setResolutions().