CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes | Friends
reco::PFSuperCluster Class Reference

Particle flow cluster, see clustering algorithm in PFSuperClusterAlgo. More...

#include <PFSuperCluster.h>

Inheritance diagram for reco::PFSuperCluster:
reco::PFCluster reco::CaloCluster

Public Member Functions

const edm::PtrVector< reco::PFCluster > & clusters () const
 vector of clusters More...
 
PFSuperClusteroperator= (const PFSuperCluster &)
 
 PFSuperCluster ()
 
 PFSuperCluster (const edm::PtrVector< reco::PFCluster > &clusters)
 constructor More...
 
void reset ()
 resets clusters parameters More...
 
- Public Member Functions inherited from reco::PFCluster
void addRecHitFraction (const reco::PFRecHitFraction &frac)
 add a given fraction of the rechit More...
 
void calculatePositionREP ()
 computes posrep_ once and for all More...
 
double charge () const
 dummy charge More...
 
double depth () const
 cluster depth More...
 
double energy () const
 cluster energy More...
 
PFLayer::Layer layer () const
 cluster layer, see PFLayer.h in this directory More...
 
PFClusteroperator= (const PFCluster &)
 
 PFCluster ()
 
 PFCluster (PFLayer::Layer layer, double energy, double x, double y, double z)
 constructor More...
 
const REPPointpositionREP () const
 cluster position: rho, eta, phi More...
 
template<typename pruner >
void pruneUsing (pruner prune)
 
double pt () const
 transverse momentum, massless approximation More...
 
const std::vector< reco::PFRecHitFraction > & recHitFractions () const
 vector of rechit fractions More...
 
void reset ()
 resets clusters parameters More...
 
void resetHitsAndFractions ()
 reset only hits and fractions More...
 
void setDepth (double depth)
 
void setLayer (PFLayer::Layer layer)
 set layer More...
 
void setTime (float time, float timeError=0)
 
void setTimeError (float timeError)
 
double theta () const
 angle More...
 
float time () const
 
float timeError () const
 
math::XYZPoint const & vertex () const
 dummy vertex access More...
 
double vx () const
 
double vy () const
 
double vz () const
 
- Public Member Functions inherited from reco::CaloCluster
void addHitAndFraction (DetId id, float fraction)
 
AlgoId algo () const
 algorithm identifier More...
 
AlgoID algoID () const
 
 CaloCluster ()
 default constructor. Sets energy and position to zero More...
 
 CaloCluster (AlgoID algoID)
 constructor with algoId, to be used in all child classes More...
 
 CaloCluster (double energy, const math::XYZPoint &position)
 constructor from values More...
 
 CaloCluster (double energy, const math::XYZPoint &position, const CaloID &caloID)
 
 CaloCluster (double energy, const math::XYZPoint &position, const CaloID &caloID, const AlgoID &algoID, uint32_t flags=0)
 
 CaloCluster (double energy, const math::XYZPoint &position, const CaloID &caloID, const std::vector< std::pair< DetId, float > > &usedHitsAndFractions, const AlgoId algoId, const DetId seedId=DetId(0), uint32_t flags=0)
 
 CaloCluster (double energy, const math::XYZPoint &position, float chi2, const std::vector< DetId > &usedHits, const AlgoId algoId, uint32_t flags=0)
 temporary compatibility constructor More...
 
const CaloIDcaloID () const
 
double correctedEnergy () const
 
float correctedEnergyUncertainty () const
 
double energy () const
 cluster energy More...
 
double eta () const
 pseudorapidity of cluster centroid More...
 
uint32_t flags () const
 
const std::vector< std::pair< DetId, float > > & hitsAndFractions () const
 
bool isInClean () const
 
bool isInUnclean () const
 
bool operator< (const CaloCluster &rhs) const
 comparison < operator More...
 
bool operator<= (const CaloCluster &rhs) const
 comparison <= operator More...
 
bool operator== (const CaloCluster &rhs) const
 comparison == operator More...
 
bool operator> (const CaloCluster &rhs) const
 comparison > operator More...
 
bool operator>= (const CaloCluster &rhs) const
 comparison >= operator More...
 
double phi () const
 azimuthal angle of cluster centroid More...
 
const math::XYZPointposition () const
 cluster centroid position More...
 
std::string printHitAndFraction (unsigned i) const
 print hitAndFraction More...
 
void reset ()
 resets the CaloCluster (position, energy, hitsAndFractions) More...
 
DetId seed () const
 return DetId of seed More...
 
void setAlgoId (const AlgoId &id)
 
void setCaloId (const CaloID &id)
 
void setCorrectedEnergy (double cenergy)
 
void setCorrectedEnergyUncertainty (float energyerr)
 
void setEnergy (double energy)
 
void setFlags (uint32_t flags)
 
void setPosition (const math::XYZPoint &p)
 
void setSeed (const DetId &id)
 
size_t size () const
 size in number of hits (e.g. in crystals for ECAL) More...
 
double x () const
 x coordinate of cluster centroid More...
 
double y () const
 y coordinate of cluster centroid More...
 
double z () const
 z coordinate of cluster centroid More...
 
virtual ~CaloCluster ()
 destructor More...
 

Private Attributes

edm::PtrVector< reco::PFClusterclusters_
 vector of clusters More...
 

Friends

class ::PFSuperClusterAlgo
 

Additional Inherited Members

- Public Types inherited from reco::PFCluster
typedef std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
 
typedef ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
 
- Public Types inherited from reco::CaloCluster
enum  AlgoId {
  island = 0, hybrid = 1, fixedMatrix = 2, dynamicHybrid = 3,
  multi5x5 = 4, particleFlow = 5, hgcal_em = 6, hgcal_had = 7,
  hgcal_mixed = 8, hfnose = 9, undefined = 1000
}
 
typedef AlgoId AlgoID
 
enum  HCalFlags { badHcalMarker = 1 }
 
enum  SCFlags { cleanOnly = 0, common = 100, uncleanOnly = 200 }
 
- Static Public Member Functions inherited from reco::PFCluster
static double getDepthCorrection (double energy, bool isBelowPS=false, bool isHadron=false)
 
- Protected Attributes inherited from reco::CaloCluster
AlgoID algoID_
 
CaloID caloID_
 bitmask for detector information More...
 
double correctedEnergy_
 
float correctedEnergyUncertainty_
 
double energy_
 cluster energy More...
 
uint32_t flags_
 
std::vector< std::pair< DetId, float > > hitsAndFractions_
 
math::XYZPoint position_
 cluster centroid position More...
 
DetId seedId_
 DetId of seed. More...
 
- Static Protected Attributes inherited from reco::CaloCluster
static const uint32_t flagsMask_ = 0x0FFFFFFF
 
static const uint32_t flagsOffset_ = 28
 

Detailed Description

Particle flow cluster, see clustering algorithm in PFSuperClusterAlgo.

A particle flow supercluster is constructed from clusters. This calculation is performed in PFSuperClusterAlgo.

Author
Chris Tully
Date
July 2012

Definition at line 24 of file PFSuperCluster.h.

Constructor & Destructor Documentation

◆ PFSuperCluster() [1/2]

reco::PFSuperCluster::PFSuperCluster ( )
inline

Definition at line 26 of file PFSuperCluster.h.

26 {}

◆ PFSuperCluster() [2/2]

PFSuperCluster::PFSuperCluster ( const edm::PtrVector< reco::PFCluster > &  clusters)

constructor

Definition at line 8 of file PFSuperCluster.cc.

Member Function Documentation

◆ clusters()

const edm::PtrVector<reco::PFCluster>& reco::PFSuperCluster::clusters ( ) const
inline

vector of clusters

Definition at line 35 of file PFSuperCluster.h.

35 { return clusters_; }

References clusters_.

Referenced by reco::operator<<().

◆ operator=()

PFSuperCluster & PFSuperCluster::operator= ( const PFSuperCluster other)

Definition at line 15 of file PFSuperCluster.cc.

15  {
17  clusters_ = other.clusters_;
18 
19  return *this;
20 }

References clusters_, reco::PFCluster::operator=(), and trackingPlots::other.

◆ reset()

void PFSuperCluster::reset ( void  )

resets clusters parameters

Definition at line 10 of file PFSuperCluster.cc.

10  {
12  clusters_.clear();
13 }

References edm::PtrVectorBase::clear(), clusters_, and reco::PFCluster::reset().

Friends And Related Function Documentation

◆ ::PFSuperClusterAlgo

friend class ::PFSuperClusterAlgo
friend

Definition at line 43 of file PFSuperCluster.h.

Member Data Documentation

◆ clusters_

edm::PtrVector<reco::PFCluster> reco::PFSuperCluster::clusters_
private

vector of clusters

Definition at line 41 of file PFSuperCluster.h.

Referenced by clusters(), operator=(), and reset().

reco::PFSuperCluster::clusters
const edm::PtrVector< reco::PFCluster > & clusters() const
vector of clusters
Definition: PFSuperCluster.h:35
reco::PFCluster::operator=
PFCluster & operator=(const PFCluster &)
Definition: PFCluster.cc:63
trackingPlots.other
other
Definition: trackingPlots.py:1467
reco::PFSuperCluster::clusters_
edm::PtrVector< reco::PFCluster > clusters_
vector of clusters
Definition: PFSuperCluster.h:41
reco::PFCluster::reset
void reset()
resets clusters parameters
Definition: PFCluster.cc:17
edm::PtrVectorBase::clear
void clear()
Clear the PtrVector.
Definition: PtrVectorBase.h:79
reco::PFCluster
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42