CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes
CaloCellGeometry Class Referenceabstract

#include <CaloCellGeometry.h>

Inheritance diagram for CaloCellGeometry:
FlatHexagon FlatTrd IdealCastorTrapezoid IdealObliquePrism IdealZDCTrapezoid IdealZPrism PreshowerStrip TruncatedPyramid

Public Types

typedef float CCGFloat
 
typedef EZMgrFL< GlobalPointCornersMgr
 
typedef EZArrayFL< GlobalPointCornersVec
 
typedef EZMgrFL< CCGFloatParMgr
 
typedef EZArrayFL< CCGFloatParVec
 
typedef std::vector< ParVecParVecVec
 
typedef HepGeom::Point3D
< CCGFloat
Pt3D
 
typedef std::vector< Pt3DPt3DVec
 
using RepCorners = std::array< RhoEtaPhi, k_cornerSize >
 
typedef HepGeom::Transform3D Tr3D
 

Public Member Functions

bool emptyCorners () const
 
virtual float etaPos () const
 
virtual float etaSpan () const
 
GlobalPoint const & getBackPoint () const
 
CornersVec const & getCorners () const
 Returns the corner points of this cell's volume. More...
 
RepCorners const & getCornersREP () const
 
virtual const GlobalPointgetPosition () const
 Returns the position of reference for this cell. More...
 
virtual GlobalPoint getPosition (CCGFloat) const
 
virtual GlobalPoint getPosition (const Pt3D &) const
 
virtual void getTransform (Tr3D &tr, Pt3DVec *lptr) const
 --------— only needed by specific utility; overloaded when needed -— More...
 
bool inside (const GlobalPoint &point) const
 Returns true if the specified point is inside this cell. More...
 
const CCGFloatparam () const
 
virtual float phiPos () const
 
virtual float phiSpan () const
 
RhoEtaPhi const & repPos () const
 
virtual float rhoPos () const
 
void setBackPoint (const GlobalPoint &pos)
 
virtual void vocalCorners (Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref) const =0
 
virtual ~CaloCellGeometry ()
 

Static Public Member Functions

static const CCGFloatcheckParmPtr (const std::vector< CCGFloat > &vd, ParVecVec &pvv)
 
static const CCGFloatgetParmPtr (const std::vector< CCGFloat > &vd, ParMgr *mgr, ParVecVec &pvv)
 

Static Public Attributes

static constexpr unsigned int k_cornerSize = 8
 
static const CCGFloat k_ScaleFromDDDtoGeant
 

Protected Member Functions

 CaloCellGeometry (CornersVec::const_reference gp, CornersMgr *mgr, const CCGFloat *par)
 
 CaloCellGeometry (const CornersVec &cv, const CCGFloat *par)
 
 CaloCellGeometry (void)
 
virtual void initCorners (CornersVec &)=0
 
void initSpan ()
 
void setCornerVec (const std::vector< GlobalPoint > &cor)
 
void setRefPoint (const GlobalPoint &pos)
 

Private Member Functions

void initBack ()
 
void initReps ()
 

Private Attributes

GlobalPoint m_backPoint
 
CornersVec m_corners
 
float m_dEta
 
float m_dPhi
 
const CCGFloatm_parms
 
GlobalPoint m_refPoint
 
RhoEtaPhi m_rep
 
std::array< RhoEtaPhi,
k_cornerSize
m_repCorners
 

Detailed Description

Abstract base class for an individual cell's geometry.

The base class declares a pure virtual function and also writes a definition (body) to force conscious acceptance of default behaviour.

If a derived class doesn't choose to override a normal virtual, it just inherits the base version's behaviour by default. If you want to provide a default behaviour but not let derived classes just inherit it "silently" like this, you can make it pure virtual but still provide a default that the derived class author has to call deliberately if he wants it:

class B {
public:
virtual bool f() = 0;
};
bool B::f() {
return true; // this is a good default, but
} // shouldn't be used blindly
class D : public B {
public:
bool f() {
return B::f(); // if D wants the default
} // behaviour, it has to say so
};
Author
J. Mans, P. Meridiani

Definition at line 50 of file CaloCellGeometry.h.

Member Typedef Documentation

Definition at line 52 of file CaloCellGeometry.h.

Definition at line 58 of file CaloCellGeometry.h.

Definition at line 57 of file CaloCellGeometry.h.

Definition at line 62 of file CaloCellGeometry.h.

Definition at line 60 of file CaloCellGeometry.h.

typedef std::vector<ParVec> CaloCellGeometry::ParVecVec

Definition at line 61 of file CaloCellGeometry.h.

Definition at line 54 of file CaloCellGeometry.h.

typedef std::vector<Pt3D> CaloCellGeometry::Pt3DVec

Definition at line 55 of file CaloCellGeometry.h.

Definition at line 66 of file CaloCellGeometry.h.

typedef HepGeom::Transform3D CaloCellGeometry::Tr3D

Definition at line 53 of file CaloCellGeometry.h.

Constructor & Destructor Documentation

CaloCellGeometry::~CaloCellGeometry ( )
virtual

Definition at line 22 of file CaloCellGeometry.cc.

22 {}
CaloCellGeometry::CaloCellGeometry ( CornersVec::const_reference  gp,
CornersMgr mgr,
const CCGFloat par 
)
protected

Definition at line 24 of file CaloCellGeometry.cc.

25  : m_refPoint(gp), m_corners(mgr), m_parms(par), m_rep(gp.perp(), gp.eta(), gp.barePhi()), m_dEta(0), m_dPhi(0) {}
const CCGFloat * m_parms
GlobalPoint m_refPoint
CaloCellGeometry::CaloCellGeometry ( const CornersVec cv,
const CCGFloat par 
)
protected

Definition at line 27 of file CaloCellGeometry.cc.

28  : m_refPoint(0.25 * (cv[0].x() + cv[1].x() + cv[2].x() + cv[3].x()),
29  0.25 * (cv[0].y() + cv[1].y() + cv[2].y() + cv[3].y()),
30  0.25 * (cv[0].z() + cv[1].z() + cv[2].z() + cv[3].z())),
31  m_corners(cv),
32  m_parms(par),
34  m_dEta(0),
35  m_dPhi(0) {}
T perp() const
Definition: PV3DBase.h:69
const CCGFloat * m_parms
T barePhi() const
Definition: PV3DBase.h:65
GlobalPoint m_refPoint
dictionary cv
Definition: cuy.py:363
T eta() const
Definition: PV3DBase.h:73
CaloCellGeometry::CaloCellGeometry ( void  )
protected

Definition at line 19 of file CaloCellGeometry.cc.

20  : m_refPoint(0., 0., 0.), m_corners(), m_parms((CCGFloat*)nullptr), m_dEta(0), m_dPhi(0) {}
const CCGFloat * m_parms
CaloCellGeometry::CCGFloat CCGFloat
GlobalPoint m_refPoint

Member Function Documentation

const float * CaloCellGeometry::checkParmPtr ( const std::vector< CCGFloat > &  vd,
CaloCellGeometry::ParVecVec pvv 
)
static

Definition at line 96 of file CaloCellGeometry.cc.

References cms::cuda::assert(), EZArrayFL< T >::begin(), cuy::ii, dqmiolumiharvest::j, EZArrayFL< T >::size(), and findQualityFiles::v.

Referenced by getParmPtr().

96  {
97  const float* pP(nullptr);
98 
99  for (unsigned int ii(0); ii != pvv.size(); ++ii) {
100  const ParVec& v(pvv[ii]);
101  assert(v.size() == vv.size());
102 
103  bool same(true);
104  for (unsigned int j(0); j != vv.size(); ++j) {
105  same = same && (fabs(vv[j] - v[j]) < 1.e-6);
106  if (!same)
107  break;
108  }
109  if (same) {
110  pP = &(*v.begin());
111  break;
112  }
113  }
114  return pP;
115 }
EZArrayFL< CCGFloat > ParVec
int ii
Definition: cuy.py:589
assert(be >=bs)
bool CaloCellGeometry::emptyCorners ( ) const
inline

Definition at line 97 of file CaloCellGeometry.h.

References EZArrayFL< T >::empty(), and m_corners.

Referenced by operator<<().

97 { return m_corners.empty(); }
bool empty() const
Definition: EZArrayFL.h:65
virtual float CaloCellGeometry::etaPos ( ) const
inlinevirtual
virtual float CaloCellGeometry::etaSpan ( ) const
inlinevirtual

Reimplemented in FlatTrd, and FlatHexagon.

Definition at line 91 of file CaloCellGeometry.h.

References m_dEta.

91 { return m_dEta; }
GlobalPoint const& CaloCellGeometry::getBackPoint ( ) const
inline

Definition at line 84 of file CaloCellGeometry.h.

References m_backPoint.

84 { return m_backPoint; }
GlobalPoint m_backPoint
CornersVec const& CaloCellGeometry::getCorners ( ) const
inline
RepCorners const& CaloCellGeometry::getCornersREP ( ) const
inline

Definition at line 77 of file CaloCellGeometry.h.

References m_repCorners.

Referenced by reco::PFRecHit::getCornersREP().

77 { return m_repCorners; }
std::array< RhoEtaPhi, k_cornerSize > m_repCorners
const float * CaloCellGeometry::getParmPtr ( const std::vector< CCGFloat > &  vd,
CaloCellGeometry::ParMgr mgr,
CaloCellGeometry::ParVecVec pvv 
)
static

Definition at line 117 of file CaloCellGeometry.cc.

References checkParmPtr(), and mps_fire::i.

Referenced by FastTimeGeometryLoader::buildGeom(), HGCalGeometryLoader::buildGeom(), HcalHardcodeGeometryLoader::fillHBHO(), HcalFlexiHardcodeGeometryLoader::fillHBHO(), HcalHardcodeGeometryLoader::fillHE(), HcalFlexiHardcodeGeometryLoader::fillHE(), HcalHardcodeGeometryLoader::fillHF(), HcalFlexiHardcodeGeometryLoader::fillHF(), CaloTowerHardcodeGeometryLoader::makeCell(), CastorHardcodeGeometryLoader::makeCell(), ZdcHardcodeGeometryLoader::makeCell(), HcalDDDGeometryLoader::makeCell(), EcalTBHodoscopeGeometryLoaderFromDDD::makeGeometry(), and CaloGeometryDBEP< T, U >::produceAligned().

119  {
120  const float* pP(checkParmPtr(vv, pvv));
121 
122  if (nullptr == pP) {
123  pvv.emplace_back(ParVec(mgr));
124  ParVec& back(pvv.back());
125  for (unsigned int i(0); i != vv.size(); ++i) {
126  back[i] = vv[i];
127  }
128  pP = &(*back.begin());
129  }
130  return pP;
131 }
EZArrayFL< CCGFloat > ParVec
static const CCGFloat * checkParmPtr(const std::vector< CCGFloat > &vd, ParVecVec &pvv)
virtual const GlobalPoint& CaloCellGeometry::getPosition ( ) const
inlinevirtual
virtual GlobalPoint CaloCellGeometry::getPosition ( CCGFloat  ) const
inlinevirtual

Reimplemented in TruncatedPyramid.

Definition at line 81 of file CaloCellGeometry.h.

References m_refPoint.

81 { return m_refPoint; }
GlobalPoint m_refPoint
virtual GlobalPoint CaloCellGeometry::getPosition ( const Pt3D ) const
inlinevirtual

Reimplemented in FlatTrd, and FlatHexagon.

Definition at line 82 of file CaloCellGeometry.h.

References m_refPoint.

82 { return m_refPoint; }
GlobalPoint m_refPoint
void CaloCellGeometry::getTransform ( Tr3D tr,
Pt3DVec lptr 
) const
virtual

--------— only needed by specific utility; overloaded when needed -—

Reimplemented in FlatTrd, TruncatedPyramid, FlatHexagon, and PreshowerStrip.

Definition at line 51 of file CaloCellGeometry.cc.

References angle(), cms::cuda::assert(), getCorners(), getPosition(), mps_fire::i, mag(), AlCaHLTBitMon_ParallelJobs::p, param(), unit(), vocalCorners(), x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

51  {
53  const Pt3D gFront(p.x(), p.y(), p.z());
54  const DPt3D dgFront(p.x(), p.y(), p.z());
55 
56  Pt3D lFront;
57  assert(nullptr != param());
58  Pt3DVec lc(8, Pt3D(0, 0, 0));
59  vocalCorners(lc, param(), lFront);
60 
61  DPt3D dlFront(lFront.x(), lFront.y(), lFront.z());
62 
63  const Pt3D lBack(0.25 * (lc[4] + lc[5] + lc[6] + lc[7]));
64  const DPt3D dlBack(lBack.x(), lBack.y(), lBack.z());
65 
66  const Pt3D dlOne(lc[0].x(), lc[0].y(), lc[0].z());
67 
68  const CornersVec& cor(getCorners());
69  DPt3DVec dkor(8, DPt3D(0, 0, 0));
70  for (unsigned int i(0); i != 8; ++i) {
71  dkor[i] = DPt3D(cor[i].x(), cor[i].y(), cor[i].z());
72  }
73 
74  DPt3D dgBack(0.25 * (dkor[4] + dkor[5] + dkor[6] + dkor[7]));
75 
76  const DVec3D dgAxis((dgBack - dgFront).unit());
77 
78  dgBack = (dgFront + (dlBack - dlFront).mag() * dgAxis);
79  const DPt3D dgOneT(dgFront + (dlOne - dlFront).mag() * (dkor[0] - dgFront).unit());
80 
81  const double dlangle((dlBack - dlFront).angle(dlOne - dlFront));
82  const double dgangle((dgBack - dgFront).angle(dgOneT - dgFront));
83  const double ddangle(dlangle - dgangle);
84 
85  const DPlane3D dgPl(dgFront, dgOneT, dgBack);
86  const DPt3D dp2(dgFront + dgPl.normal().unit());
87 
88  const DPt3D dgOne(dgFront + HepGeom::Rotate3D(-ddangle, dgFront, dp2) * DVec3D(dgOneT - dgFront));
89 
90  tr = Tr3D(dlFront, dlBack, dlOne, dgFront, dgBack, dgOne);
91 
92  if (nullptr != lptr)
93  (*lptr) = lc;
94 }
EZArrayFL< GlobalPoint > CornersVec
HepGeom::Point3D< double > DPt3D
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
virtual const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
HepGeom::Transform3D Tr3D
assert(be >=bs)
const CCGFloat * param() const
HepGeom::Vector3D< double > DVec3D
virtual void vocalCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref) const =0
CaloCellGeometry::Pt3D Pt3D
HepGeom::Plane3D< double > DPlane3D
CaloCellGeometry::Pt3DVec Pt3DVec
HepGeom::Point3D< CCGFloat > Pt3D
CornersVec const & getCorners() const
Returns the corner points of this cell&#39;s volume.
std::vector< DPt3D > DPt3DVec
Basic3DVector unit() const
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
void CaloCellGeometry::initBack ( )
inlineprivate

Definition at line 137 of file CaloCellGeometry.h.

References cuy::cv, getCorners(), m_backPoint, x, y, and z.

Referenced by initSpan().

137  {
138  // from CaloTower code
139  CornersVec const& cv = getCorners();
140  m_backPoint = GlobalPoint(0.25 * (cv[4].x() + cv[5].x() + cv[6].x() + cv[7].x()),
141  0.25 * (cv[4].y() + cv[5].y() + cv[6].y() + cv[7].y()),
142  0.25 * (cv[4].z() + cv[5].z() + cv[6].z() + cv[7].z()));
143  }
EZArrayFL< GlobalPoint > CornersVec
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
GlobalPoint m_backPoint
dictionary cv
Definition: cuy.py:363
CornersVec const & getCorners() const
Returns the corner points of this cell&#39;s volume.
virtual void CaloCellGeometry::initCorners ( CornersVec )
protectedpure virtual
void CaloCellGeometry::initReps ( )
inlineprivate

Definition at line 144 of file CaloCellGeometry.h.

References getCorners(), mps_fire::i, k_cornerSize, and m_repCorners.

Referenced by initSpan().

144  {
145  for (auto i = 0U; i < k_cornerSize; ++i)
146  m_repCorners[i] = {getCorners()[i].perp(), getCorners()[i].eta(), getCorners()[i].barePhi()};
147  }
static constexpr unsigned int k_cornerSize
std::array< RhoEtaPhi, k_cornerSize > m_repCorners
CornersVec const & getCorners() const
Returns the corner points of this cell&#39;s volume.
void CaloCellGeometry::initSpan ( )
inlineprotected
bool CaloCellGeometry::inside ( const GlobalPoint point) const

Returns true if the specified point is inside this cell.

Definition at line 133 of file CaloCellGeometry.cc.

References cms::cuda::co, EE, getCorners(), mps_fire::i, AlCaHLTBitMon_ParallelJobs::p, x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

133  {
134  bool ans(false);
135  const Pt3D p(point.x(), point.y(), point.z());
136  const CornersVec& cog(getCorners());
137  Pt3D co[8];
138  for (unsigned int i(0); i != 8; ++i) {
139  co[i] = Pt3D(cog[i].x(), cog[i].y(), cog[i].z());
140  }
141 
142  const Plane3D AA(co[0], co[1], co[2]); // z<0
143  const Plane3D BB(co[6], co[5], co[4]); // z>0
144 
145  if (AA.distance(p) * BB.distance(p) >= 0) {
146  const Plane3D CC(co[0], co[4], co[5]); // x<0
147  const Plane3D DD(co[2], co[6], co[7]); // x>0
148  if (CC.distance(p) * DD.distance(p) >= 0) {
149  const Plane3D EE(co[3], co[7], co[4]); // y<0
150  const Plane3D FF(co[1], co[5], co[6]); // y>0
151  if (EE.distance(p) * FF.distance(p) >= 0) {
152  ans = true;
153  }
154  }
155  }
156  return ans;
157 }
EZArrayFL< GlobalPoint > CornersVec
T y() const
Definition: PV3DBase.h:60
__host__ __device__ VT * co
Definition: prefixScan.h:47
T z() const
Definition: PV3DBase.h:61
CaloCellGeometry::Pt3D Pt3D
HepGeom::Point3D< CCGFloat > Pt3D
CornersVec const & getCorners() const
Returns the corner points of this cell&#39;s volume.
T x() const
Definition: PV3DBase.h:59
ROOT::Math::Plane3D Plane3D
const CCGFloat* CaloCellGeometry::param ( ) const
inline

Definition at line 99 of file CaloCellGeometry.h.

References m_parms.

Referenced by FWTGeoRecoGeometryESProducer::addCaloTowerGeometry(), FWTGeoRecoGeometryESProducer::addEcalCaloGeometry(), FWTGeoRecoGeometryESProducer::addHcalCaloGeometryBarrel(), FWTGeoRecoGeometryESProducer::addHcalCaloGeometryEndcap(), FWTGeoRecoGeometryESProducer::addHcalCaloGeometryForward(), FWTGeoRecoGeometryESProducer::addHcalCaloGeometryOuter(), IdealZDCTrapezoid::an(), IdealCastorTrapezoid::an(), FlatHexagon::backCtr(), FlatTrd::backCtr(), IdealObliquePrism::dEta(), IdealZPrism::dEta(), IdealCastorTrapezoid::dh(), IdealObliquePrism::dPhi(), IdealZPrism::dPhi(), IdealCastorTrapezoid::dR(), PreshowerStrip::dx(), IdealZDCTrapezoid::dx(), IdealCastorTrapezoid::dxh(), IdealCastorTrapezoid::dxl(), PreshowerStrip::dy(), IdealZDCTrapezoid::dy(), PreshowerStrip::dz(), IdealZDCTrapezoid::dz(), IdealObliquePrism::dz(), IdealZPrism::dz(), IdealCastorTrapezoid::dz(), FlatHexagon::dz(), FlatTrd::dz(), IdealObliquePrism::eta(), IdealZPrism::eta(), FlatHexagon::etaSpan(), FlatTrd::etaSpan(), EcalTBHodoscopeGeometry::getGeometryRawPtr(), HcalDDDGeometry::getGeometryRawPtr(), ZdcGeometry::getGeometryRawPtr(), CastorGeometry::getGeometryRawPtr(), CaloTowerGeometry::getGeometryRawPtr(), EcalPreshowerGeometry::getGeometryRawPtr(), EcalEndcapGeometry::getGeometryRawPtr(), FastTimeGeometry::getGeometryRawPtr(), EcalBarrelGeometry::getGeometryRawPtr(), HcalGeometry::getGeometryRawPtr(), HGCalGeometry::getGeometryRawPtr(), FlatHexagon::getTransform(), TruncatedPyramid::getTransform(), FlatTrd::getTransform(), getTransform(), IdealZDCTrapezoid::initCorners(), IdealCastorTrapezoid::initCorners(), makeEcalShape(), operator<<(), FlatHexagon::phiSpan(), FlatTrd::phiSpan(), PreshowerStrip::tilt(), IdealObliquePrism::z(), and IdealZPrism::z().

99 { return m_parms; }
const CCGFloat * m_parms
virtual float CaloCellGeometry::phiPos ( ) const
inlinevirtual
virtual float CaloCellGeometry::phiSpan ( ) const
inlinevirtual

Reimplemented in FlatTrd, and FlatHexagon.

Definition at line 92 of file CaloCellGeometry.h.

References m_dPhi.

92 { return m_dPhi; }
RhoEtaPhi const& CaloCellGeometry::repPos ( ) const
inline

Definition at line 86 of file CaloCellGeometry.h.

References m_rep.

Referenced by reco::PFRecHit::positionREP().

86 { return m_rep; }
virtual float CaloCellGeometry::rhoPos ( ) const
inlinevirtual

Definition at line 87 of file CaloCellGeometry.h.

References m_rep, and RhoEtaPhi::rho().

87 { return m_rep.rho(); }
float rho() const
transverse momentum
Definition: PtEtaPhiMass.h:50
void CaloCellGeometry::setBackPoint ( const GlobalPoint pos)
inline

Definition at line 108 of file CaloCellGeometry.h.

References m_backPoint.

108 { m_backPoint = pos; }
GlobalPoint m_backPoint
void CaloCellGeometry::setCornerVec ( const std::vector< GlobalPoint > &  cor)
inlineprotected

Definition at line 131 of file CaloCellGeometry.h.

References isotrackApplyRegressor::k, and m_corners.

131  {
132  for (unsigned int k = 0; k < cor.size(); ++k)
133  m_corners[k] = cor[k];
134  }
void CaloCellGeometry::setRefPoint ( const GlobalPoint pos)
inlineprotected

Definition at line 130 of file CaloCellGeometry.h.

References m_refPoint.

Referenced by FlatHexagon::setPosition(), and FlatTrd::setPosition().

130 { m_refPoint = pos; }
GlobalPoint m_refPoint
virtual void CaloCellGeometry::vocalCorners ( Pt3DVec vec,
const CCGFloat pv,
Pt3D ref 
) const
pure virtual

Member Data Documentation

constexpr unsigned int CaloCellGeometry::k_cornerSize = 8
static

Definition at line 64 of file CaloCellGeometry.h.

Referenced by CaloSubdetectorGeometry::allocateCorners(), and initReps().

const float CaloCellGeometry::k_ScaleFromDDDtoGeant
static
GlobalPoint CaloCellGeometry::m_backPoint
private

Definition at line 150 of file CaloCellGeometry.h.

Referenced by getBackPoint(), initBack(), and setBackPoint().

CornersVec CaloCellGeometry::m_corners
private

Definition at line 151 of file CaloCellGeometry.h.

Referenced by emptyCorners(), getCorners(), initSpan(), and setCornerVec().

float CaloCellGeometry::m_dEta
private

Definition at line 154 of file CaloCellGeometry.h.

Referenced by etaSpan(), and initSpan().

float CaloCellGeometry::m_dPhi
private

Definition at line 155 of file CaloCellGeometry.h.

Referenced by initSpan(), and phiSpan().

const CCGFloat* CaloCellGeometry::m_parms
private

Definition at line 152 of file CaloCellGeometry.h.

Referenced by param().

GlobalPoint CaloCellGeometry::m_refPoint
private

Definition at line 149 of file CaloCellGeometry.h.

Referenced by getPosition(), and setRefPoint().

RhoEtaPhi CaloCellGeometry::m_rep
private

Definition at line 153 of file CaloCellGeometry.h.

Referenced by etaPos(), phiPos(), repPos(), and rhoPos().

std::array<RhoEtaPhi, k_cornerSize> CaloCellGeometry::m_repCorners
private

Definition at line 156 of file CaloCellGeometry.h.

Referenced by getCornersREP(), and initReps().