CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
reco::SuperCluster Class Reference

#include <DataFormats/EgammaReco/interface/SuperCluster.h>

Inheritance diagram for reco::SuperCluster:
reco::CaloCluster

Public Types

typedef math::XYZPoint Point
 
- 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,
  hfnose = 9, undefined = 1000
}
 
typedef AlgoId AlgoID
 
enum  HCalFlags { badHcalMarker = 1 }
 
enum  SCFlags { cleanOnly = 0, common = 100, uncleanOnly = 200 }
 

Public Member Functions

void addCluster (const CaloClusterPtr &r)
 add reference to constituent BasicCluster More...
 
void addPreshowerCluster (const CaloClusterPtr &r)
 add reference to constituent BasicCluster More...
 
void clearHitsAndFractions ()
 
const CaloClusterPtrVectorclusters () const
 const access to the cluster list itself More...
 
CaloCluster_iterator clustersBegin () const
 fist iterator over BasicCluster constituents More...
 
CaloCluster_iterator clustersEnd () const
 last iterator over BasicCluster constituents More...
 
size_t clustersSize () const
 number of BasicCluster constituents More...
 
double etaWidth () const
 
const int getPreshowerPlanesStatus () const
 
double phiWidth () const
 obtain phi and eta width of the Super Cluster More...
 
const CaloClusterPtrVectorpreshowerClusters () const
 const access to the preshower cluster list itself More...
 
CaloCluster_iterator preshowerClustersBegin () const
 fist iterator over PreshowerCluster constituents More...
 
CaloCluster_iterator preshowerClustersEnd () const
 last iterator over PreshowerCluster constituents More...
 
size_t preshowerClustersSize () const
 number of BasicCluster PreShower constituents More...
 
double preshowerEnergy () const
 energy deposited in preshower More...
 
double preshowerEnergyPlane1 () const
 
double preshowerEnergyPlane2 () const
 
double rawEnergy () const
 raw uncorrected energy (sum of energies of component BasicClusters) More...
 
const CaloClusterPtrseed () const
 seed BasicCluster More...
 
const int seedCrysIEtaOrIx () const
 
const int seedCrysIPhiOrIy () const
 
void setClusters (const CaloClusterPtrVector &clusters)
 
void setEtaWidth (double ew)
 
void setPhiWidth (double pw)
 
void setPreshowerClusters (const CaloClusterPtrVector &clusters)
 
void setPreshowerEnergy (double preshowerEnergy)
 
void setPreshowerEnergyPlane1 (double preshowerEnergy1)
 
void setPreshowerEnergyPlane2 (double preshowerEnergy2)
 
void setPreshowerPlanesStatus (const uint32_t &status)
 
void setSeed (const CaloClusterPtr &r)
 list of used xtals by DetId // now inherited by CaloCluster More...
 
 SuperCluster ()
 default constructor More...
 
 SuperCluster (double energy, const Point &position)
 constructor defined by CaloCluster - will have to use setSeed and add() separately More...
 
 SuperCluster (double energy, const Point &position, const CaloClusterPtr &seed, const CaloClusterPtrVector &clusters, double Epreshower=0., double phiWidth=0., double etaWidth=0., double Epreshower1=0., double Epreshower2=0.)
 
 SuperCluster (double energy, const Point &position, const CaloClusterPtr &seed, const CaloClusterPtrVector &clusters, const CaloClusterPtrVector &preshowerClusters, double Epreshower=0., double phiWidth=0., double etaWidth=0., double Epreshower1=0., double Epreshower2=0.)
 
- 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, const CaloID &caloID)
 
 CaloCluster (double energy, const math::XYZPoint &position)
 constructor from values More...
 
 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 Member Functions

void computeRawEnergy ()
 

Private Attributes

CaloClusterPtrVector clusters_
 references to BasicCluster constitunets More...
 
double etaWidth_
 
double phiWidth_
 
CaloClusterPtrVector preshowerClusters_
 references to BasicCluster constitunets More...
 
double preshowerEnergy1_
 
double preshowerEnergy2_
 
double preshowerEnergy_
 used hits by detId - retrieved from BC constituents – now inherited from CaloCluster More...
 
double rawEnergy_
 
CaloClusterPtr seed_
 reference to BasicCluster seed More...
 

Additional Inherited Members

- 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

A SuperCluster reconstructed in the Electromagnetic Calorimeter contains references to seed and constituent BasicClusters

Author
Luca Lista, INFN

Definition at line 20 of file SuperCluster.h.

Member Typedef Documentation

◆ Point

Definition at line 22 of file SuperCluster.h.

Constructor & Destructor Documentation

◆ SuperCluster() [1/4]

reco::SuperCluster::SuperCluster ( )
inline

default constructor

Definition at line 25 of file SuperCluster.h.

26  : CaloCluster(0., Point(0., 0., 0.)),
28  rawEnergy_(-1.),
29  phiWidth_(0),
30  etaWidth_(0),
32  preshowerEnergy2_(0) {}
double preshowerEnergy_
used hits by detId - retrieved from BC constituents – now inherited from CaloCluster ...
Definition: SuperCluster.h:191
CaloCluster()
default constructor. Sets energy and position to zero
Definition: CaloCluster.h:56
math::XYZPoint Point
Definition: SuperCluster.h:22

◆ SuperCluster() [2/4]

reco::SuperCluster::SuperCluster ( double  energy,
const Point position 
)

constructor defined by CaloCluster - will have to use setSeed and add() separately

◆ SuperCluster() [3/4]

reco::SuperCluster::SuperCluster ( double  energy,
const Point position,
const CaloClusterPtr seed,
const CaloClusterPtrVector clusters,
double  Epreshower = 0.,
double  phiWidth = 0.,
double  etaWidth = 0.,
double  Epreshower1 = 0.,
double  Epreshower2 = 0. 
)

◆ SuperCluster() [4/4]

reco::SuperCluster::SuperCluster ( double  energy,
const Point position,
const CaloClusterPtr seed,
const CaloClusterPtrVector clusters,
const CaloClusterPtrVector preshowerClusters,
double  Epreshower = 0.,
double  phiWidth = 0.,
double  etaWidth = 0.,
double  Epreshower1 = 0.,
double  Epreshower2 = 0. 
)

Member Function Documentation

◆ addCluster()

void reco::SuperCluster::addCluster ( const CaloClusterPtr r)
inline

add reference to constituent BasicCluster

Definition at line 124 of file SuperCluster.h.

References clusters_, computeRawEnergy(), and edm::PtrVector< T >::push_back().

Referenced by PFEGammaAlgo::buildRefinedSuperCluster(), PFElectronTranslator::createSuperClusters(), PFPhotonTranslator::createSuperClusters(), and PFECALSuperClusterAlgo::finalizeSuperCluster().

124  {
127  }
CaloClusterPtrVector clusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:183
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149

◆ addPreshowerCluster()

void reco::SuperCluster::addPreshowerCluster ( const CaloClusterPtr r)
inline

add reference to constituent BasicCluster

Definition at line 130 of file SuperCluster.h.

References preshowerClusters_, and edm::PtrVector< T >::push_back().

Referenced by PFEGammaAlgo::buildRefinedSuperCluster(), PFElectronTranslator::createSuperClusters(), PFPhotonTranslator::createSuperClusters(), and PFECALSuperClusterAlgo::finalizeSuperCluster().

void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149
CaloClusterPtrVector preshowerClusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:186

◆ clearHitsAndFractions()

void reco::SuperCluster::clearHitsAndFractions ( )
inline

Definition at line 121 of file SuperCluster.h.

References reco::CaloCluster::hitsAndFractions_.

121 { hitsAndFractions_.clear(); }
std::vector< std::pair< DetId, float > > hitsAndFractions_
Definition: CaloCluster.h:233

◆ clusters()

const CaloClusterPtrVector& reco::SuperCluster::clusters ( ) const
inline

◆ clustersBegin()

CaloCluster_iterator reco::SuperCluster::clustersBegin ( ) const
inline

◆ clustersEnd()

CaloCluster_iterator reco::SuperCluster::clustersEnd ( ) const
inline

◆ clustersSize()

size_t reco::SuperCluster::clustersSize ( ) const
inline

number of BasicCluster constituents

Definition at line 100 of file SuperCluster.h.

References clusters_, and edm::PtrVectorBase::size().

Referenced by SuperClusterHelper::clustersSize(), EcalRegressionData::fill(), lowptgsfeleid::findEnergeticClusters(), and lowptgsfeleid::trackClusterMatching().

100 { return clusters_.size(); }
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
CaloClusterPtrVector clusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:183

◆ computeRawEnergy()

void SuperCluster::computeRawEnergy ( )
private

Definition at line 82 of file SuperCluster.cc.

References clustersBegin(), clustersEnd(), and rawEnergy_.

Referenced by addCluster(), and setClusters().

82  {
83  rawEnergy_ = 0.;
84  for (CaloClusterPtrVector::const_iterator bcItr = clustersBegin(); bcItr != clustersEnd(); bcItr++) {
85  rawEnergy_ += (*bcItr)->energy();
86  }
87 }
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:88
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:91

◆ etaWidth()

double reco::SuperCluster::etaWidth ( ) const
inline

◆ getPreshowerPlanesStatus()

const int reco::SuperCluster::getPreshowerPlanesStatus ( ) const
inline

Get preshower planes status : 0 : both planes working 1 : only first plane working 2 : only second plane working 3 : both planes dead

Definition at line 148 of file SuperCluster.h.

References reco::CaloCluster::flags_, and reco::CaloCluster::flagsOffset_.

148 { return (flags_ >> flagsOffset_); }
static const uint32_t flagsOffset_
Definition: CaloCluster.h:247

◆ phiWidth()

double reco::SuperCluster::phiWidth ( ) const
inline

◆ preshowerClusters()

const CaloClusterPtrVector& reco::SuperCluster::preshowerClusters ( ) const
inline

const access to the preshower cluster list itself

Definition at line 85 of file SuperCluster.h.

References preshowerClusters_.

Referenced by ReducedEGProducer::linkCaloClusters(), and ReducedEGProducer::relinkCaloClusters().

85 { return preshowerClusters_; }
CaloClusterPtrVector preshowerClusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:186

◆ preshowerClustersBegin()

CaloCluster_iterator reco::SuperCluster::preshowerClustersBegin ( ) const
inline

fist iterator over PreshowerCluster constituents

Definition at line 94 of file SuperCluster.h.

References edm::PtrVector< T >::begin(), and preshowerClusters_.

Referenced by PFEGammaAlgo::buildRefinedSuperCluster(), PFECALSuperClusterAlgo::finalizeSuperCluster(), and SuperClusterHelper::SuperClusterHelper().

94 { return preshowerClusters_.begin(); }
CaloClusterPtrVector preshowerClusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:186
const_iterator begin() const
Definition: PtrVector.h:144

◆ preshowerClustersEnd()

CaloCluster_iterator reco::SuperCluster::preshowerClustersEnd ( ) const
inline

last iterator over PreshowerCluster constituents

Definition at line 97 of file SuperCluster.h.

References edm::PtrVector< T >::end(), and preshowerClusters_.

Referenced by PFEGammaAlgo::buildRefinedSuperCluster(), PFECALSuperClusterAlgo::finalizeSuperCluster(), and SuperClusterHelper::SuperClusterHelper().

97 { return preshowerClusters_.end(); }
CaloClusterPtrVector preshowerClusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:186
const_iterator end() const
Definition: PtrVector.h:146

◆ preshowerClustersSize()

size_t reco::SuperCluster::preshowerClustersSize ( ) const
inline

number of BasicCluster PreShower constituents

Definition at line 103 of file SuperCluster.h.

References preshowerClusters_, and edm::PtrVectorBase::size().

103 { return preshowerClusters_.size(); }
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
CaloClusterPtrVector preshowerClusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:186

◆ preshowerEnergy()

double reco::SuperCluster::preshowerEnergy ( ) const
inline

◆ preshowerEnergyPlane1()

double reco::SuperCluster::preshowerEnergyPlane1 ( ) const
inline

Definition at line 64 of file SuperCluster.h.

References preshowerEnergy1_.

64 { return preshowerEnergy1_; }

◆ preshowerEnergyPlane2()

double reco::SuperCluster::preshowerEnergyPlane2 ( ) const
inline

Definition at line 65 of file SuperCluster.h.

References preshowerEnergy2_.

65 { return preshowerEnergy2_; }

◆ rawEnergy()

double reco::SuperCluster::rawEnergy ( ) const
inline

◆ seed()

const CaloClusterPtr& reco::SuperCluster::seed ( ) const
inline

◆ seedCrysIEtaOrIx()

const int reco::SuperCluster::seedCrysIEtaOrIx ( ) const
inline

Definition at line 150 of file SuperCluster.h.

References EcalBarrel, EcalEndcap, EBDetId::ieta(), EEDetId::ix(), and seed_.

150  {
151  auto detid = seed_->seed();
152  int ietaorix = 0;
153  if (detid.subdetId() == EcalBarrel) {
154  EBDetId ebdetid(detid);
155  ietaorix = ebdetid.ieta();
156  } else if (detid.subdetId() == EcalEndcap) {
157  EEDetId eedetid(detid);
158  ietaorix = eedetid.ix();
159  }
160  return ietaorix;
161  }
CaloClusterPtr seed_
reference to BasicCluster seed
Definition: SuperCluster.h:180

◆ seedCrysIPhiOrIy()

const int reco::SuperCluster::seedCrysIPhiOrIy ( ) const
inline

Definition at line 163 of file SuperCluster.h.

References EcalBarrel, EcalEndcap, EBDetId::iphi(), EEDetId::iy(), and seed_.

163  {
164  auto detid = seed_->seed();
165  int iphioriy = 0;
166  if (detid.subdetId() == EcalBarrel) {
167  EBDetId ebdetid(detid);
168  iphioriy = ebdetid.iphi();
169  } else if (detid.subdetId() == EcalEndcap) {
170  EEDetId eedetid(detid);
171  iphioriy = eedetid.iy();
172  }
173  return iphioriy;
174  }
CaloClusterPtr seed_
reference to BasicCluster seed
Definition: SuperCluster.h:180

◆ setClusters()

void reco::SuperCluster::setClusters ( const CaloClusterPtrVector clusters)
inline

Definition at line 112 of file SuperCluster.h.

References clusters(), clusters_, and computeRawEnergy().

Referenced by LowPtGsfElectronSCProducer::produce(), and ReducedEGProducer::relinkCaloClusters().

112  {
115  }
CaloClusterPtrVector clusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:183
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:82

◆ setEtaWidth()

void reco::SuperCluster::setEtaWidth ( double  ew)
inline

◆ setPhiWidth()

void reco::SuperCluster::setPhiWidth ( double  pw)
inline

◆ setPreshowerClusters()

void reco::SuperCluster::setPreshowerClusters ( const CaloClusterPtrVector clusters)
inline

Definition at line 118 of file SuperCluster.h.

References clusters(), and preshowerClusters_.

Referenced by ReducedEGProducer::relinkCaloClusters().

CaloClusterPtrVector preshowerClusters_
references to BasicCluster constitunets
Definition: SuperCluster.h:186
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:82

◆ setPreshowerEnergy()

void reco::SuperCluster::setPreshowerEnergy ( double  preshowerEnergy)
inline

Definition at line 72 of file SuperCluster.h.

References preshowerEnergy(), and preshowerEnergy_.

Referenced by PFEGammaAlgo::buildRefinedSuperCluster(), PFElectronTranslator::createSuperClusters(), PFPhotonTranslator::createSuperClusters(), and PFECALSuperClusterAlgo::finalizeSuperCluster().

double preshowerEnergy_
used hits by detId - retrieved from BC constituents – now inherited from CaloCluster ...
Definition: SuperCluster.h:191
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:63

◆ setPreshowerEnergyPlane1()

void reco::SuperCluster::setPreshowerEnergyPlane1 ( double  preshowerEnergy1)
inline

◆ setPreshowerEnergyPlane2()

void reco::SuperCluster::setPreshowerEnergyPlane2 ( double  preshowerEnergy2)
inline

Definition at line 74 of file SuperCluster.h.

References preshowerEnergy2_.

Referenced by PFEGammaAlgo::buildRefinedSuperCluster(), and PFECALSuperClusterAlgo::finalizeSuperCluster().

74 { preshowerEnergy2_ = preshowerEnergy2; };

◆ setPreshowerPlanesStatus()

void reco::SuperCluster::setPreshowerPlanesStatus ( const uint32_t &  status)
inline

Set preshower planes status : 0 : both planes working 1 : only first plane working 2 : only second plane working 3 : both planes dead

Definition at line 138 of file SuperCluster.h.

References reco::CaloCluster::flags(), reco::CaloCluster::flags_, reco::CaloCluster::flagsMask_, reco::CaloCluster::flagsOffset_, and mps_update::status.

Referenced by PreshowerClusterProducer::produce().

138  {
139  uint32_t flags = flags_ & flagsMask_;
140  flags_ = flags | (status << flagsOffset_);
141  }
uint32_t flags() const
Definition: CaloCluster.h:192
static const uint32_t flagsOffset_
Definition: CaloCluster.h:247
static const uint32_t flagsMask_
Definition: CaloCluster.h:246

◆ setSeed()

void reco::SuperCluster::setSeed ( const CaloClusterPtr r)
inline

Member Data Documentation

◆ clusters_

CaloClusterPtrVector reco::SuperCluster::clusters_
private

references to BasicCluster constitunets

Definition at line 183 of file SuperCluster.h.

Referenced by addCluster(), clusters(), clustersBegin(), clustersEnd(), clustersSize(), and setClusters().

◆ etaWidth_

double reco::SuperCluster::etaWidth_
private

Definition at line 196 of file SuperCluster.h.

Referenced by etaWidth(), and setEtaWidth().

◆ phiWidth_

double reco::SuperCluster::phiWidth_
private

Definition at line 195 of file SuperCluster.h.

Referenced by phiWidth(), and setPhiWidth().

◆ preshowerClusters_

CaloClusterPtrVector reco::SuperCluster::preshowerClusters_
private

references to BasicCluster constitunets

Definition at line 186 of file SuperCluster.h.

Referenced by addPreshowerCluster(), preshowerClusters(), preshowerClustersBegin(), preshowerClustersEnd(), preshowerClustersSize(), and setPreshowerClusters().

◆ preshowerEnergy1_

double reco::SuperCluster::preshowerEnergy1_
private

Definition at line 198 of file SuperCluster.h.

Referenced by preshowerEnergyPlane1(), and setPreshowerEnergyPlane1().

◆ preshowerEnergy2_

double reco::SuperCluster::preshowerEnergy2_
private

Definition at line 199 of file SuperCluster.h.

Referenced by preshowerEnergyPlane2(), and setPreshowerEnergyPlane2().

◆ preshowerEnergy_

double reco::SuperCluster::preshowerEnergy_
private

used hits by detId - retrieved from BC constituents – now inherited from CaloCluster

Definition at line 191 of file SuperCluster.h.

Referenced by preshowerEnergy(), and setPreshowerEnergy().

◆ rawEnergy_

double reco::SuperCluster::rawEnergy_
private

Definition at line 193 of file SuperCluster.h.

Referenced by computeRawEnergy(), and rawEnergy().

◆ seed_

CaloClusterPtr reco::SuperCluster::seed_
private

reference to BasicCluster seed

Definition at line 180 of file SuperCluster.h.

Referenced by seed(), seedCrysIEtaOrIx(), seedCrysIPhiOrIy(), and setSeed().