CMS 3D CMS Logo

reco::PFCluster Class Reference

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

#include <DataFormats/ParticleFlowReco/interface/PFCluster.h>

Inheritance diagram for reco::PFCluster:

reco::CaloCluster

List of all members.

Public Types

typedef
ROOT::Math::PositionVector3D
< ROOT::Math::CylindricalEta3D
< Double32_t > > 
REPPoint

Public Member Functions

void addRecHitFraction (const reco::PFRecHitFraction &frac)
 add a given fraction of the rechit
void calculatePositionREP ()
 calculates posrep_ once and for all
int color () const
 
Returns:
color

double energy () const
 cluster energy
PFLayer::Layer layer () const
 cluster layer, see PFLayer.h in this directory
PFClusteroperator= (const PFCluster &)
 PFCluster (PFLayer::Layer layer, double energy, double x, double y, double z)
 constructor
 PFCluster ()
 default constructor
const REPPointpositionREP () const
 cluster position: rho, eta, phi
const std::vector
< reco::PFRecHitFraction > & 
recHitFractions () const
 vector of rechit fractions
void reset ()
 resets clusters parameters
void setColor (int color)
 set cluster color (for the PFRootEventManager display)
void setLayer (PFLayer::Layer layer)
 set layer

Static Public Member Functions

static double getDepthCorrection (double energy, bool isBelowPS=false, bool isHadron=false)
static void setDepthCorParameters (int mode, double a, double b, double ap, double bp)

Static Public Attributes

static unsigned instanceCounter_ = 0
 counter

Private Attributes

int color_
 color (transient)
REPPoint posrep_
 cluster position: rho, eta, phi (transient)
std::vector
< reco::PFRecHitFraction
rechits_
 vector of rechit fractions (transient)

Static Private Attributes

static double depthCorA_ = 0.89
static double depthCorAp_ = 0.89
static double depthCorB_ = 7.3
static double depthCorBp_ = 4.0
static int depthCorMode_ = 0

Friends

class ::PFClusterAlgo
std::ostream & operator<< (std::ostream &out, const PFCluster &cluster)


Detailed Description

Particle flow cluster, see clustering algorithm in PFClusterAlgo.

A particle flow cluster is defined by its energy and position, which are calculated from a vector of PFRecHitFraction. This calculation is performed in PFClusterAlgo.

Todo:
Clean up this class to a common base (talk to Paolo Meridiani) the extra internal stuff (like the vector of PFRecHitFraction's) could be moved to a PFClusterExtra.
Todo:
Now that PFRecHitFraction's hold a reference to the PFRecHit's, put back the calculation of energy and position to PFCluster.
Todo:
Add an operator+=
Author:
Colin Bernet
Date:
July 2006

Definition at line 42 of file PFCluster.h.


Member Typedef Documentation

typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<Double32_t> > reco::PFCluster::REPPoint

Definition at line 46 of file PFCluster.h.


Constructor & Destructor Documentation

reco::PFCluster::PFCluster (  )  [inline]

default constructor

Definition at line 49 of file PFCluster.h.

00049 :  color_(1) {}

PFCluster::PFCluster ( PFLayer::Layer  layer,
double  energy,
double  x,
double  y,
double  z 
)

constructor

Definition at line 17 of file PFCluster.cc.

00018                                                     : 
00019   CaloCluster( energy, math::XYZPoint(x,y,z), PFLayer::toCaloID(layer) ),
00020   posrep_( position_.Rho(), position_.Eta(), position_.Phi() ),
00021   color_(2)
00022 {  }


Member Function Documentation

void PFCluster::addRecHitFraction ( const reco::PFRecHitFraction frac  ) 

add a given fraction of the rechit

Definition at line 35 of file PFCluster.cc.

References rechits_.

Referenced by PFClusterAlgo::buildPFClusters().

00035                                                                     {
00036 
00037   rechits_.push_back( frac );
00038 }

void reco::PFCluster::calculatePositionREP (  )  [inline]

calculates posrep_ once and for all

Definition at line 80 of file PFCluster.h.

References reco::CaloCluster::position_, and posrep_.

Referenced by PFEnergyCalibration::energyEm().

00080                                 {
00081       posrep_.SetCoordinates( position_.Rho(), position_.Eta(), position_.Phi() ); 
00082     }

int reco::PFCluster::color (  )  const [inline]

Returns:
color

Definition at line 92 of file PFCluster.h.

References color_.

00092 {return color_;}

double reco::PFCluster::energy (  )  const [inline]

cluster energy

Reimplemented from reco::CaloCluster.

Definition at line 74 of file PFCluster.h.

References reco::CaloCluster::energy_.

Referenced by PFEnergyCalibration::energyEm(), GPFCluster::GPFCluster(), reco::operator<<(), PFConversionAlgo::setCandidates(), PFBlockAlgo::testECALAndHCAL(), PFBlockAlgo::testLinkByRecHit(), PFBlockAlgo::testPSAndECAL(), PFBlockAlgo::testTrackAndECAL(), and PFBlockAlgo::testTrackAndHCAL().

00074 {return energy_;}

double PFCluster::getDepthCorrection ( double  energy,
bool  isBelowPS = false,
bool  isHadron = false 
) [static]

Todo:
move to PFClusterTools

Definition at line 41 of file PFCluster.cc.

References depthCorA_, depthCorAp_, depthCorB_, depthCorBp_, and funct::log().

Referenced by PFBlockAlgo::testTrackAndECAL().

00043 {
00044   double corrA = depthCorA_;
00045   double corrB = depthCorB_;
00046   if (isBelowPS) {
00047     corrA = depthCorAp_;
00048     corrB = depthCorBp_;
00049   }
00050   double depth = 0;
00051   switch(isHadron) {
00052   case 0: // e/gamma
00053     depth = corrA*(corrB + log(energy)); 
00054     break;
00055   case 1: // hadrons
00056     depth = corrA;
00057     break;
00058   default:
00059     assert(0);
00060     //     edm::LogError("PFCluster") << "unknown function for depth correction!"
00061     //                         << std::endl;
00062   }
00063   return depth;
00064 }

PFLayer::Layer PFCluster::layer (  )  const

cluster layer, see PFLayer.h in this directory

Definition at line 73 of file PFCluster.cc.

References reco::CaloCluster::caloID(), and PFLayer::fromCaloID().

Referenced by PFClusterAlgo::calculateClusterPosition(), DisplayManager::createGCluster(), PFAlgo::isSatelliteCluster(), reco::operator<<(), PFAlgo::reconstructCluster(), PFBlockAlgo::testLinkByRecHit(), PFBlockAlgo::testPSAndECAL(), and PFBlockAlgo::testTrackAndPS().

00073                                      {
00074   
00075   // cout<<"calling PFCluster::layer "<<caloID()<<" "<<PFLayer::fromCaloID( caloID() )<<endl;
00076   return PFLayer::fromCaloID( caloID() );
00077 }     

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

Definition at line 80 of file PFCluster.cc.

References color_, reco::CaloCluster::energy_, reco::CaloCluster::operator=(), reco::CaloCluster::position_, posrep_, and rechits_.

00080                                                       {
00081 
00082   CaloCluster::operator=(other); 
00083   rechits_ = other.rechits_;
00084   energy_ = other.energy_;
00085   position_ = other.position_;
00086   posrep_ = other.posrep_;
00087   color_ = other.color_;
00088 
00089   return *this;
00090 }

const REPPoint& reco::PFCluster::positionREP (  )  const [inline]

cluster position: rho, eta, phi

Definition at line 77 of file PFCluster.h.

References posrep_.

Referenced by PFEnergyCalibration::energyEm(), PFAlgo::isSatelliteCluster(), reco::operator<<(), PFBlockAlgo::testECALAndHCAL(), PFBlockAlgo::testLinkByRecHit(), PFBlockAlgo::testPSAndECAL(), PFBlockAlgo::testTrackAndECAL(), and PFBlockAlgo::testTrackAndHCAL().

00077 {return posrep_;}

const std::vector< reco::PFRecHitFraction >& reco::PFCluster::recHitFractions (  )  const [inline]

vector of rechit fractions

Definition at line 64 of file PFCluster.h.

References rechits_.

Referenced by reco::operator<<(), and PFBlockAlgo::testLinkByRecHit().

00065       { return rechits_; }

void PFCluster::reset ( void   ) 

resets clusters parameters

Definition at line 25 of file PFCluster.cc.

References reco::CaloCluster::energy_, reco::CaloCluster::position_, posrep_, and rechits_.

00025                       {
00026   
00027   energy_ = 0;
00028   position_ *= 0;
00029   posrep_ *= 0;
00030   
00031   rechits_.clear();
00032 }

void reco::PFCluster::setColor ( int  color  )  [inline]

set cluster color (for the PFRootEventManager display)

Definition at line 89 of file PFCluster.h.

References color_.

00089 {color_ = color;}

static void reco::PFCluster::setDepthCorParameters ( int  mode,
double  a,
double  b,
double  ap,
double  bp 
) [inline, static]

Todo:
move to PFClusterTools

Definition at line 103 of file PFCluster.h.

References depthCorA_, depthCorAp_, depthCorB_, depthCorBp_, and depthCorMode_.

Referenced by PFClusterProducer::PFClusterProducer(), and PFRootEventManager::readOptions().

00105                                                              {
00106       depthCorMode_ = mode;
00107       depthCorA_ = a; 
00108       depthCorB_ = b; 
00109       depthCorAp_ = ap; 
00110       depthCorBp_ = bp; 
00111     } 

void PFCluster::setLayer ( PFLayer::Layer  layer  ) 

set layer

Definition at line 66 of file PFCluster.cc.

References reco::CaloCluster::caloID_, and PFLayer::toCaloID().

Referenced by PFClusterAlgo::calculateClusterPosition().

00066                                             {
00067   // cout<<"calling PFCluster::setLayer "<<layer<<endl;
00068   caloID_ = PFLayer::toCaloID( layer );
00069   // cout<<"done "<<caloID_<<endl;
00070 }


Friends And Related Function Documentation

friend class ::PFClusterAlgo [friend]

Definition at line 138 of file PFCluster.h.

std::ostream& operator<< ( std::ostream &  out,
const PFCluster cluster 
) [friend]

Definition at line 93 of file PFCluster.cc.

00094                                                          {
00095   
00096   if(!out) return out;
00097   
00098   const math::XYZPoint&  pos = cluster.position();
00099   const PFCluster::REPPoint&  posrep = cluster.positionREP();
00100   const std::vector< reco::PFRecHitFraction >& fracs = 
00101     cluster.recHitFractions();
00102   
00103   out<<"cluster "
00104      <<"\tlayer: "<<cluster.layer()
00105      <<"\tenergy: "<<cluster.energy()
00106      <<"\tXYZ: "
00107      <<pos.X()<<","<<pos.Y()<<","<<pos.Z()<<" | "
00108      <<"\tREP: "
00109      <<posrep.Rho()<<","<<posrep.Eta()<<","<<posrep.Phi()<<" | "
00110      <<fracs.size()<<" rechits: ";
00111   for(unsigned i=0; i<fracs.size(); i++) {
00112     out<<fracs[i]<<", ";
00113   }
00114   
00115   return out;
00116 }


Member Data Documentation

int reco::PFCluster::color_ [private]

color (transient)

Definition at line 136 of file PFCluster.h.

Referenced by color(), operator=(), and setColor().

double PFCluster::depthCorA_ = 0.89 [static, private]

Todo:
move to PFClusterTools

Definition at line 126 of file PFCluster.h.

Referenced by PFClusterAlgo::calculateClusterPosition(), getDepthCorrection(), and setDepthCorParameters().

double PFCluster::depthCorAp_ = 0.89 [static, private]

Todo:
move to PFClusterTools

Definition at line 130 of file PFCluster.h.

Referenced by PFClusterAlgo::calculateClusterPosition(), getDepthCorrection(), and setDepthCorParameters().

double PFCluster::depthCorB_ = 7.3 [static, private]

Todo:
move to PFClusterTools

Definition at line 128 of file PFCluster.h.

Referenced by PFClusterAlgo::calculateClusterPosition(), getDepthCorrection(), and setDepthCorParameters().

double PFCluster::depthCorBp_ = 4.0 [static, private]

Todo:
move to PFClusterTools

Definition at line 132 of file PFCluster.h.

Referenced by PFClusterAlgo::calculateClusterPosition(), getDepthCorrection(), and setDepthCorParameters().

int PFCluster::depthCorMode_ = 0 [static, private]

Todo:
move to PFClusterTools

Definition at line 124 of file PFCluster.h.

Referenced by PFClusterAlgo::calculateClusterPosition(), and setDepthCorParameters().

unsigned PFCluster::instanceCounter_ = 0 [static]

counter

Definition at line 100 of file PFCluster.h.

REPPoint reco::PFCluster::posrep_ [private]

cluster position: rho, eta, phi (transient)

Definition at line 120 of file PFCluster.h.

Referenced by PFClusterAlgo::calculateClusterPosition(), calculatePositionREP(), operator=(), positionREP(), and reset().

std::vector< reco::PFRecHitFraction > reco::PFCluster::rechits_ [private]

vector of rechit fractions (transient)

Definition at line 117 of file PFCluster.h.

Referenced by addRecHitFraction(), PFClusterAlgo::calculateClusterPosition(), operator=(), recHitFractions(), and reset().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:51:21 2009 for CMSSW by  doxygen 1.5.4