CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 (edm::Handle< edm::View< reco::PFJet > > PFJets)
 
TMatrixD getSignifMatrix () const
 
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 17 of file SignPFSpecificAlgo.h.

Constructor & Destructor Documentation

metsig::SignPFSpecificAlgo::SignPFSpecificAlgo ( )

Definition at line 10 of file SignPFSpecificAlgo.cc.

References clusteredParticlePtrs_.

10  :
11 resolutions_(0),
12 algo_()
13 {clusteredParticlePtrs_.clear();}
metsig::significanceAlgo algo_
metsig::SignAlgoResolutions * resolutions_
std::set< reco::CandidatePtr > clusteredParticlePtrs_
metsig::SignPFSpecificAlgo::~SignPFSpecificAlgo ( )

Definition at line 15 of file SignPFSpecificAlgo.cc.

15 {}

Member Function Documentation

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

Definition at line 57 of file SignPFSpecificAlgo.cc.

57  {
58  if(clusteredParticlePtrs_.find(pf) != clusteredParticlePtrs_.end()) {
59  return; //pf candidate already added in jet collection
60  }
61  std::vector<metsig::SigInputObj> vobj;
62  vobj.push_back(resolutions_->evalPF(&(*pf)));
63  algo_.addObjects(vobj);
64 }
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 metsig::SignPFSpecificAlgo::addPFJets ( edm::Handle< edm::View< reco::PFJet > >  PFJets)

Definition at line 24 of file SignPFSpecificAlgo.cc.

References metsig::jet.

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

Definition at line 28 of file SignPFSpecificAlgo.h.

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

28 {return algo_.getSignifMatrix();}
TMatrixD getSignifMatrix() const
metsig::significanceAlgo algo_
void metsig::SignPFSpecificAlgo::setResolutions ( metsig::SignAlgoResolutions resolutions)

Definition at line 18 of file SignPFSpecificAlgo.cc.

18  {
19  resolutions_ = resolutions;
20 }
metsig::SignAlgoResolutions * resolutions_
void metsig::SignPFSpecificAlgo::useOriginalPtrs ( const edm::ProductID productID)

Definition at line 38 of file SignPFSpecificAlgo.cc.

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

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

Member Data Documentation

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

Definition at line 34 of file SignPFSpecificAlgo.h.

Referenced by getSignifMatrix().

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

Definition at line 33 of file SignPFSpecificAlgo.h.

Referenced by SignPFSpecificAlgo().

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

Definition at line 32 of file SignPFSpecificAlgo.h.