CMS 3D CMS Logo

Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes | Static Private Attributes | Friends

reco::PFCluster Class Reference

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

#include <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 ()
 computes posrep_ once and for all
double charge () const
 dummy charge
int color () const
double energy () const
 cluster energy
PFLayer::Layer layer () const
 cluster layer, see PFLayer.h in this directory
PFClusteroperator= (const PFCluster &)
 PFCluster ()
 PFCluster (PFLayer::Layer layer, double energy, double x, double y, double z)
 constructor
const REPPointpositionREP () const
 cluster position: rho, eta, phi
double pt () const
 transverse momentum, massless approximation
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
double theta () const
 angle
math::XYZPoint const & vertex () const
 dummy vertex access
double vx () const
double vy () const
double vz () const

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.

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]

Definition at line 48 of file PFCluster.h.

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

constructor

Definition at line 16 of file PFCluster.cc.


Member Function Documentation

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

add a given fraction of the rechit

Definition at line 40 of file PFCluster.cc.

References reco::CaloCluster::addHitAndFraction(), reco::PFRecHitFraction::fraction(), reco::PFRecHitFraction::recHitRef(), and rechits_.

Referenced by PFClusterAlgo::buildPFClusters().

                                                                    {

  rechits_.push_back( frac );

  addHitAndFraction( frac.recHitRef()->detId(), 
                     frac.fraction() );
}
void reco::PFCluster::calculatePositionREP ( ) [inline]

computes posrep_ once and for all

Definition at line 79 of file PFCluster.h.

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

Referenced by PFEnergyCalibration::energyEm(), PFElecTkProducer::isSharingEcalEnergyWithEgSC(), and ConvBremPFTrackFinder::runConvBremFinder().

                                {
      posrep_.SetCoordinates( position_.Rho(), 
                              position_.Eta(), 
                              position_.Phi() ); 
    }
double reco::PFCluster::charge ( ) const [inline]

dummy charge

some classes to make this fit into a template footprint for RecoPFClusterRefCandidate so we can make jets and MET out of PFClusters.

Definition at line 120 of file PFCluster.h.

{ return 0;}
int reco::PFCluster::color ( ) const [inline]
Returns:
color

Definition at line 93 of file PFCluster.h.

References color_.

Referenced by setColor().

{return color_;}
double reco::PFCluster::energy ( ) const [inline]
double PFCluster::getDepthCorrection ( double  energy,
bool  isBelowPS = false,
bool  isHadron = false 
) [static]

Definition at line 49 of file PFCluster.cc.

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

{
  double corrA = depthCorA_;
  double corrB = depthCorB_;
  if (isBelowPS) {
    corrA = depthCorAp_;
    corrB = depthCorBp_;
  }
  double depth = 0;
  switch(isHadron) {
  case 0: // e/gamma
    depth = corrA*(corrB + log(energy)); 
    break;
  case 1: // hadrons
    depth = corrA;
    break;
  default:
    assert(0);
    //     edm::LogError("PFCluster") << "unknown function for depth correction!"
    //                         << std::endl;
  }
  return depth;
}
PFLayer::Layer PFCluster::layer ( ) const
PFCluster & PFCluster::operator= ( const PFCluster other)

Definition at line 88 of file PFCluster.cc.

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

                                                      {

  CaloCluster::operator=(other); 
  rechits_ = other.rechits_;
  energy_ = other.energy_;
  position_ = other.position_;
  posrep_ = other.posrep_;
  color_ = other.color_;

  return *this;
}
const REPPoint& reco::PFCluster::positionREP ( ) const [inline]
double reco::PFCluster::pt ( ) const [inline]

transverse momentum, massless approximation

Definition at line 123 of file PFCluster.h.

References energy(), reco::CaloCluster::position_, and funct::sin().

                      { 
      return (energy() * sin(position_.theta()));
    }
const std::vector< reco::PFRecHitFraction >& reco::PFCluster::recHitFractions ( ) const [inline]

vector of rechit fractions

Definition at line 63 of file PFCluster.h.

References rechits_.

Referenced by LinkByRecHit::testECALAndPSByRecHit(), and LinkByRecHit::testTrackAndClusterByRecHit().

      { return rechits_; }
void PFCluster::reset ( void  )

resets clusters parameters

Reimplemented from reco::CaloCluster.

Definition at line 27 of file PFCluster.cc.

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

                      {
  
  energy_ = 0;
  position_ *= 0;
  posrep_ *= 0;
  
  rechits_.clear();

  CaloCluster::reset();
  
}
void reco::PFCluster::setColor ( int  color) [inline]

set cluster color (for the PFRootEventManager display)

Definition at line 90 of file PFCluster.h.

References color(), and color_.

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

set layer

Definition at line 74 of file PFCluster.cc.

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

Referenced by PFClusterAlgo::calculateClusterPosition().

                                            {
  // cout<<"calling PFCluster::setLayer "<<layer<<endl;
  caloID_ = PFLayer::toCaloID( layer );
  // cout<<"done "<<caloID_<<endl;
}
double reco::PFCluster::theta ( ) const [inline]

angle

Definition at line 128 of file PFCluster.h.

References posrep_.

                         { 
      return posrep_.theta();
    }
math::XYZPoint const& reco::PFCluster::vertex ( ) const [inline]

dummy vertex access

Definition at line 133 of file PFCluster.h.

Referenced by vx(), vy(), and vz().

                                        { 
      static math::XYZPoint dummyVtx(0,0,0);
      return dummyVtx;      
    }
double reco::PFCluster::vx ( ) const [inline]

Definition at line 137 of file PFCluster.h.

References vertex().

{ return vertex().x(); }
double reco::PFCluster::vy ( ) const [inline]

Definition at line 138 of file PFCluster.h.

References vertex().

{ return vertex().y(); }
double reco::PFCluster::vz ( ) const [inline]

Definition at line 139 of file PFCluster.h.

References vertex().

{ return vertex().z(); }    

Friends And Related Function Documentation

friend class ::PFClusterAlgo [friend]

Definition at line 169 of file PFCluster.h.

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

Member Data Documentation

int reco::PFCluster::color_ [private]

color (transient)

Definition at line 167 of file PFCluster.h.

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

double PFCluster::depthCorA_ = 0.89 [static, private]
double PFCluster::depthCorAp_ = 0.89 [static, private]
double PFCluster::depthCorB_ = 7.3 [static, private]
double PFCluster::depthCorBp_ = 4.0 [static, private]
int PFCluster::depthCorMode_ = 0 [static, private]

Definition at line 151 of file PFCluster.h.

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

unsigned PFCluster::instanceCounter_ = 0 [static]

counter

Definition at line 101 of file PFCluster.h.

cluster position: rho, eta, phi (transient)

Definition at line 147 of file PFCluster.h.

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

vector of rechit fractions (transient)

Definition at line 144 of file PFCluster.h.

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