CMS 3D CMS Logo

Public Member Functions | Private Attributes | Friends

CSCRadialStripTopology Class Reference

#include <CSCRadialStripTopology.h>

Inheritance diagram for CSCRadialStripTopology:
RadialStripTopology StripTopology Topology OffsetRadialStripTopology CSCStripTopology CSCGangedStripTopology CSCUngangedStripTopology

List of all members.

Public Member Functions

float angularWidth () const
float centreToIntersection () const
virtual int channel (const LocalPoint &) const
 CSCRadialStripTopology (int ns, float aw, float dh, float r, float yAx=1.f, float yMid=0.)
float detHeight () const
virtual LocalError localError (float strip, float stripErr2) const
virtual LocalError localError (const MeasurementPoint &, const MeasurementError &) const
virtual float localPitch (const LocalPoint &) const
virtual LocalPoint localPosition (const MeasurementPoint &) const
virtual LocalPoint localPosition (float strip) const
virtual float localStripLength (const LocalPoint &) const
virtual MeasurementError measurementError (const LocalPoint &, const LocalError &) const
virtual MeasurementPoint measurementPosition (const LocalPoint &) const
virtual int nearestStrip (const LocalPoint &) const
virtual int nstrips () const
float originToIntersection () const
float phiOfOneEdge () const
virtual float phiPitch (void) const
virtual float strip (const LocalPoint &) const
virtual float stripAngle (float strip) const
virtual float stripLength () const
float xOfStrip (int strip, float y) const
float yAxisOrientation () const
float yCentreOfStripPlane () const
float yDistanceToIntersection (float y) const
float yExtentOfStripPlane () const
virtual ~CSCRadialStripTopology ()

Private Attributes

float theAngularWidth
float theCentreToIntersection
float theDetHeight
int theNumberOfStrips
float thePhiOfOneEdge
float theYAxisOrientation
float yCentre

Friends

std::ostream & operator<< (std::ostream &, const RadialStripTopology &)

Detailed Description

Definition at line 28 of file CSCRadialStripTopology.h.


Constructor & Destructor Documentation

CSCRadialStripTopology::CSCRadialStripTopology ( int  ns,
float  aw,
float  dh,
float  r,
float  yAx = 1.f,
float  yMid = 0. 
)

Constructor from:

Parameters:
nsnumber of strips
awangular width of a strip
dhdetector height (usually 2 x apothem of TrapezoidalPlaneBounds)
rradial distance from symmetry centre of detector to the point at which the outer edges of the two extreme strips (projected) intersect.
yAxorientation of local y axis: 1 means pointing from the smaller side of the module to the larger side (along apothem), and -1 means in the opposite direction, i.e. from the larger side along the apothem to the smaller side. Default value is 1.
yMidlocal y offset if mid-point of detector (strip plane) does not coincide with local origin. This decouples the extent of strip plane from the boundary of the detector in which the RST is embedded.

Definition at line 7 of file CSCRadialStripTopology.cc.

References LogTrace, theAngularWidth, theNumberOfStrips, thePhiOfOneEdge, theYAxisOrientation, and yCentre.

                                                                                                         :
  theNumberOfStrips(ns), theAngularWidth(aw), 
  theDetHeight(dh), theCentreToIntersection(r),
  theYAxisOrientation(yAx), yCentre( yMid) {   
  // Angular offset of extreme edge of detector, so that angle is
  // zero for a strip lying along local y axis = long symmetry axis of plane of strips
  thePhiOfOneEdge = -(0.5*theNumberOfStrips) * theAngularWidth * yAx;
  
  LogTrace("CSCRadialStripTopology") << "CSCRadialStripTopology: constructed with"
        << " strips = " << ns
        << " width = " << aw << " rad "
        << " det_height = " << dh
        << " ctoi = " << r 
        << " phi_edge = " << thePhiOfOneEdge << " rad "
        << " y_ax_ori = " << theYAxisOrientation
        << " y_det_centre = " << yCentre 
        << "\n";
}    
virtual CSCRadialStripTopology::~CSCRadialStripTopology ( ) [inline, virtual]

Destructor

Definition at line 50 of file CSCRadialStripTopology.h.

{}

Member Function Documentation

float CSCRadialStripTopology::angularWidth ( ) const [inline, virtual]
float CSCRadialStripTopology::centreToIntersection ( ) const [inline, virtual]

Distance from the intersection of the projections of the extreme edges of the two extreme strips to the symmetry centre of the plane of strips.

Implements RadialStripTopology.

Definition at line 185 of file CSCRadialStripTopology.h.

References theCentreToIntersection.

Referenced by FWRecoGeometryESProducer::addCSCGeometry(), FWTGeoRecoGeometryESProducer::addCSCGeometry(), localError(), and ValidateGeometry::validateCSCLayerGeometry().

int CSCRadialStripTopology::channel ( const LocalPoint lp) const [virtual]

Channel number corresponding to a given LocalPoint.
This is effectively an integer version of strip(), with range 0 to nstrips-1.
LocalPoints outside the detector strip plane will be considered as contributing to the edge channels 0 or nstrips-1.

Implements RadialStripTopology.

Reimplemented in OffsetRadialStripTopology, CSCGangedStripTopology, and CSCUngangedStripTopology.

Definition at line 26 of file CSCRadialStripTopology.cc.

References min, strip(), and theNumberOfStrips.

{ return   std::min( int( strip(lp) ), theNumberOfStrips-1 ) ;}
float CSCRadialStripTopology::detHeight ( ) const [inline, virtual]

Length of long symmetry axis of plane of strips

Implements RadialStripTopology.

Definition at line 173 of file CSCRadialStripTopology.h.

References theDetHeight.

Referenced by localError(), OffsetRadialStripTopology::localPosition(), localPosition(), localStripLength(), measurementError(), and measurementPosition().

{ return theDetHeight;}
LocalError CSCRadialStripTopology::localError ( const MeasurementPoint mp,
const MeasurementError me 
) const [virtual]

LocalError for a given MeasurementPoint with known MeasurementError. This may be used in Kalman filtering and hence must allow possible correlations between the components.

Implements RadialStripTopology.

Definition at line 84 of file CSCRadialStripTopology.cc.

References angularWidth(), alignmentValidation::c1, centreToIntersection(), funct::cos(), fwrapper::cs, detHeight(), phi, dttmaxenums::R, findQualityFiles::rr, indexGen::s2, funct::sin(), stripAngle(), groupFilesInBlocks::tt, MeasurementError::uu(), MeasurementError::uv(), MeasurementError::vv(), PV2DBase< T, PVType, FrameType >::x(), create_public_lumi_plots::xy, PV2DBase< T, PVType, FrameType >::y(), and yAxisOrientation().

                                                                                                          {
  const double
    phi(stripAngle(mp.x())), s1(std::sin(phi)), c1(std::cos(phi)),
    cs(s1*c1), s2(s1*s1), c2(1-s2), // rotation matrix

    T( angularWidth() * ( centreToIntersection() + yAxisOrientation()*mp.y()*detHeight()) / c1 ), // tangential measurement unit (local pitch)
    R( detHeight()/ c1 ),  // radial measurement unit (strip length)
    tt( me.uu() * T*T ),   // tangential sigma^2
    rr( me.vv() * R*R   ), // radial sigma^2
    tr( me.uv() * T*R ),  
    
    xx(  c2*tt  +  2*cs*tr  +  s2*rr      ),
    yy(  s2*tt  -  2*cs*tr  +  c2*rr      ),
    xy( cs*( rr - tt )  +  tr*( c2 - s2 ) );

  return LocalError( xx, xy, yy );
}
LocalError CSCRadialStripTopology::localError ( float  strip,
float  stripErr2 
) const [virtual]

LocalError for a pure strip measurement, where 'strip' is the (float) position (a 'phi' angle wrt y axis) and stripErr2 is the sigma-squared. Both quantities are expressed in units of theAngularWidth of a strip.

Implements RadialStripTopology.

Definition at line 68 of file CSCRadialStripTopology.cc.

References angularWidth(), centreToIntersection(), detHeight(), f, phi, funct::pow(), findQualityFiles::rr, stripAngle(), funct::tan(), groupFilesInBlocks::tt, and create_public_lumi_plots::xy.

                                                                                {
  const double
    phi(stripAngle(strip)), t1(std::tan(phi)), t2(t1*t1),
    // s1(std::sin(phi)), c1(std::cos(phi)),
    // cs(s1*c1), s2(s1*s1), c2(1-s2), // rotation matrix

    tt( stripErr2 * std::pow( centreToIntersection()*angularWidth() ,2.f) ), // tangential sigma^2   *c2
    rr( std::pow(detHeight(), 2.f) * (1.f/12.f) ),                                   // radial sigma^2( uniform prob density along strip)  *c2

    xx( tt + t2*rr  ),
    yy( t2*tt + rr  ),
    xy( t1*( rr - tt ) );
  
  return LocalError( xx, xy, yy );
}
float CSCRadialStripTopology::localPitch ( const LocalPoint lp) const [virtual]

Pitch (strip width) at a given LocalPoint.
BEWARE: are you sure you really want to call this for a RadialStripTopology?

Implements RadialStripTopology.

Definition at line 119 of file CSCRadialStripTopology.cc.

References angularWidth(), funct::cos(), f, min, nstrips(), funct::pow(), funct::sin(), strip(), stripAngle(), PV3DBase< T, PVType, FrameType >::y(), and yDistanceToIntersection().

Referenced by CSCLayerGeometry::stripPitch().

                                                                   { 
 // The local pitch is the local x width of the strip at the local (x,y)
  const int istrip = std::min(nstrips(), static_cast<int>(strip(lp)) + 1); // which strip number
  const float fangle = stripAngle(static_cast<float>(istrip) - 0.5); // angle of strip centre
  return
    yDistanceToIntersection( lp.y() ) * std::sin(angularWidth()) /
    std::pow( std::cos(fangle-0.5f*angularWidth()), 2.f);
}
LocalPoint CSCRadialStripTopology::localPosition ( const MeasurementPoint mp) const [virtual]

LocalPoint for a given MeasurementPoint
What's a MeasurementPoint?
In analogy with that used with TrapezoidalStripTopology objects, a MeasurementPoint is a 2-dim object.
The first dimension measures the angular position wrt central line of symmetry of detector, in units of strip (angular) widths (range 0 to total angle subtended by a detector).
The second dimension measures the fractional position along the strip (range -0.5 to +0.5).
BEWARE! The components are not Cartesian.

Implements RadialStripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 54 of file CSCRadialStripTopology.cc.

References detHeight(), stripAngle(), funct::tan(), PV2DBase< T, PVType, FrameType >::x(), x, PV2DBase< T, PVType, FrameType >::y(), detailsBasic3DVector::y, yAxisOrientation(), yCentreOfStripPlane(), and yDistanceToIntersection().

                                                                                 {
  const float  // y = (L/cos(phi))*mp.y()*cos(phi) 
    y( mp.y()*detHeight()  +  yCentreOfStripPlane() ),
    x( yAxisOrientation() * yDistanceToIntersection( y ) * std::tan ( stripAngle( mp.x() ) ) );
  return LocalPoint( x, y );
}
LocalPoint CSCRadialStripTopology::localPosition ( float  strip) const [virtual]

LocalPoint on x axis for given 'strip' 'strip' is a float in units of the strip (angular) width

Implements RadialStripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 50 of file CSCRadialStripTopology.cc.

References originToIntersection(), stripAngle(), funct::tan(), and yAxisOrientation().

float CSCRadialStripTopology::localStripLength ( const LocalPoint lp) const [virtual]

Length of a strip passing through a given LocalPpoint

Implements RadialStripTopology.

Definition at line 34 of file CSCRadialStripTopology.cc.

References detHeight(), f, funct::pow(), mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and yDistanceToIntersection().

                                                                         {  
  return detHeight() * std::sqrt(1.f + std::pow( lp.x()/yDistanceToIntersection(lp.y()), 2.f) );
}
MeasurementError CSCRadialStripTopology::measurementError ( const LocalPoint p,
const LocalError e 
) const [virtual]

Implements RadialStripTopology.

Definition at line 102 of file CSCRadialStripTopology.cc.

References angularWidth(), fwrapper::cs, detHeight(), f, funct::pow(), indexGen::s2, mathSSE::sqrt(), lumiQTWidget::t, PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), LocalError::xy(), PV3DBase< T, PVType, FrameType >::y(), yAxisOrientation(), yDistanceToIntersection(), and LocalError::yy().

                                                                                                         {
  const double
    yHitToInter(yDistanceToIntersection(p.y())),
    t(yAxisOrientation() * p.x() / yHitToInter),   // tan(strip angle) 
    cs(t/(1+t*t)), s2(t*cs), c2(1-s2),             // rotation matrix

    T2( 1./(std::pow(angularWidth(),2.f) * ( std::pow(p.x(),2.f) + std::pow(yHitToInter,2)) )), // 1./tangential measurement unit (local pitch) ^2
    R2( c2/std::pow(detHeight(),2.f) ),                                    // 1./ radial measurement unit (strip length) ^2

    uu(       ( c2*e.xx() - 2*cs*e.xy() + s2*e.yy() )   * T2 ),
    vv(       ( s2*e.xx() + 2*cs*e.xy() + c2*e.yy() )   * R2 ),
    uv( ( cs*( e.xx() - e.yy() ) + e.xy()*( c2 - s2 ) )  * std::sqrt (T2*R2) );
  
  return MeasurementError(uu, uv, vv);
}
MeasurementPoint CSCRadialStripTopology::measurementPosition ( const LocalPoint lp) const [virtual]

Implements RadialStripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 61 of file CSCRadialStripTopology.cc.

References angularWidth(), detHeight(), phi, phiOfOneEdge(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), yAxisOrientation(), yCentreOfStripPlane(), and yDistanceToIntersection().

                                                                                       {
  const float // phi is [pi/2 - conventional local phi], use atan2(x,y) rather than atan2(y,x)
    phi( yAxisOrientation() * std::atan2( lp.x(), yDistanceToIntersection( lp.y() ) ));
  return MeasurementPoint( yAxisOrientation()*( phi-phiOfOneEdge() ) / angularWidth(),
                           ( lp.y() - yCentreOfStripPlane() )        / detHeight() );
}
int CSCRadialStripTopology::nearestStrip ( const LocalPoint lp) const [virtual]

Nearest strip to given LocalPoint

Implements RadialStripTopology.

Definition at line 28 of file CSCRadialStripTopology.cc.

References max(), min, nstrips(), and strip().

Referenced by CSCLayerGeometry::nearestStrip().

{   return std::min( nstrips(), static_cast<int>( std::max(float(0), strip(lp)) ) + 1);}
virtual int CSCRadialStripTopology::nstrips ( ) const [inline, virtual]

Total number of strips

Implements RadialStripTopology.

Definition at line 124 of file CSCRadialStripTopology.h.

References theNumberOfStrips.

Referenced by localPitch(), nearestStrip(), CSCLayerGeometry::numberOfStrips(), OffsetRadialStripTopology::strip(), and strip().

{ return theNumberOfStrips; }
float CSCRadialStripTopology::originToIntersection ( ) const [inline, virtual]

(y) distance from intersection of the projections of the strips to the local coordinate origin. Same as centreToIntersection() if symmetry centre of strip plane coincides with local origin.

Implements RadialStripTopology.

Definition at line 192 of file CSCRadialStripTopology.h.

References theCentreToIntersection, and yCentre.

Referenced by CSCStripTopology::equationOfStrip(), OffsetRadialStripTopology::localPosition(), localPosition(), OffsetRadialStripTopology::strip(), OffsetRadialStripTopology::toLocal(), OffsetRadialStripTopology::toPrime(), and yDistanceToIntersection().

float CSCRadialStripTopology::phiOfOneEdge ( ) const [inline, virtual]

Convenience function to access azimuthal angle of extreme edge of first strip measured relative to long symmetry axis of the plane of strips.

WARNING! This angle is measured clockwise from the local y axis which means it is in the conventional azimuthal phi plane, but azimuth is of course measured from local x axis not y. The range of this angle is -(full angle)/2 to +(full angle)/2.
where (full angle) = nstrips() * angularWidth().

Implements RadialStripTopology.

Definition at line 206 of file CSCRadialStripTopology.h.

References thePhiOfOneEdge.

Referenced by FWRecoGeometryESProducer::addCSCGeometry(), FWTGeoRecoGeometryESProducer::addCSCGeometry(), OffsetRadialStripTopology::localPosition(), measurementPosition(), OffsetRadialStripTopology::strip(), strip(), stripAngle(), OffsetRadialStripTopology::stripAngle(), and ValidateGeometry::validateCSCLayerGeometry().

{ return thePhiOfOneEdge; }
virtual float CSCRadialStripTopology::phiPitch ( void  ) const [inline, virtual]

Phi pitch of each strip (= angular width!)

Implements RadialStripTopology.

Definition at line 168 of file CSCRadialStripTopology.h.

References angularWidth().

Referenced by CSCLayerGeometry::stripPhiPitch().

{ return angularWidth(); }
float CSCRadialStripTopology::strip ( const LocalPoint lp) const [virtual]

Strip in which a given LocalPoint lies. This is a float which represents the fractional strip position within the detector.
Returns zero if the LocalPoint falls at the extreme low edge of the detector or BELOW, and float(nstrips) if it falls at the extreme high edge or ABOVE.

Implements RadialStripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 43 of file CSCRadialStripTopology.cc.

References angularWidth(), max(), min, nstrips(), phi, phiOfOneEdge(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), yAxisOrientation(), and yDistanceToIntersection().

Referenced by channel(), localPitch(), and nearestStrip().

                                                              {
  const float   // phi is measured from y axis --> sign of angle is sign of x * yAxisOrientation --> use atan2(x,y), not atan2(y,x)
    phi( std::atan2( lp.x(), yDistanceToIntersection( lp.y() ) )),
    aStrip( ( phi - yAxisOrientation() * phiOfOneEdge() )/angularWidth());
  return  std::max(float(0), std::min( (float)nstrips(), aStrip ));
}
float CSCRadialStripTopology::stripAngle ( float  strip) const [virtual]

Angle between strip and symmetry axis (=local y axis) for given strip.
This is like a phi angle but measured clockwise from y axis rather than counter clockwise from x axis. Note that 'strip' is a float with a continuous range from 0 to float(nstrips) to cover the whole detector, and the centres of strips correspond to half-integer values 0.5, 1.5, ..., nstrips-0.5 whereas values 1, 2, ... nstrips correspond to the upper phi edges of the strips.

Implements RadialStripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 30 of file CSCRadialStripTopology.cc.

References angularWidth(), phiOfOneEdge(), and yAxisOrientation().

Referenced by localError(), localPitch(), localPosition(), and xOfStrip().

virtual float CSCRadialStripTopology::stripLength ( ) const [inline, virtual]

Height of detector (= length of long symmetry axis of the plane of strips).

Implements RadialStripTopology.

Definition at line 129 of file CSCRadialStripTopology.h.

References theDetHeight.

{ return theDetHeight; }
float CSCRadialStripTopology::xOfStrip ( int  strip,
float  y 
) const [virtual]

Local x where centre of strip intersects input local y
'strip' should be in range 1 to nstrips()

Implements RadialStripTopology.

Definition at line 38 of file CSCRadialStripTopology.cc.

References stripAngle(), funct::tan(), yAxisOrientation(), and yDistanceToIntersection().

Referenced by CSCLayerGeometry::xOfStrip().

                                                               { 
  return   
    yAxisOrientation() * yDistanceToIntersection( y ) * std::tan( stripAngle(static_cast<float>(strip) - 0.5 ) );
}
float CSCRadialStripTopology::yAxisOrientation ( ) const [inline, virtual]
float CSCRadialStripTopology::yCentreOfStripPlane ( ) const [inline, virtual]
float CSCRadialStripTopology::yDistanceToIntersection ( float  y) const [virtual]

Distance in local y from a hit to the point of intersection of projected strips

Implements RadialStripTopology.

Definition at line 32 of file CSCRadialStripTopology.cc.

References originToIntersection(), and yAxisOrientation().

Referenced by localPitch(), localPosition(), localStripLength(), measurementError(), measurementPosition(), strip(), and xOfStrip().

float CSCRadialStripTopology::yExtentOfStripPlane ( ) const [inline, virtual]

y extent of strip plane

Implements RadialStripTopology.

Definition at line 178 of file CSCRadialStripTopology.h.

References theDetHeight.

Referenced by CSCStripTopology::yLimitsOfStripPlane().

{ return theDetHeight; } // same as detHeight()

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const RadialStripTopology rst 
) [friend]

Reimplemented from RadialStripTopology.

Definition at line 10 of file RadialStripTopology.cc.

                                                                            {
  os  << "RadialStripTopology " << std::endl
      << " " << std::endl
      << "number of strips          " << rst.nstrips() << std::endl
      << "centre to whereStripsMeet " << rst.centreToIntersection() << std::endl
      << "detector height in y      " << rst.detHeight() << std::endl
      << "angular width of strips   " << rst.phiPitch() << std::endl
      << "phi of one edge           " << rst.phiOfOneEdge() << std::endl
      << "y axis orientation        " << rst.yAxisOrientation() << std::endl
      << "y of centre of strip plane " << rst.yCentreOfStripPlane() << std::endl;
  return os;
}

Member Data Documentation

Definition at line 239 of file CSCRadialStripTopology.h.

Referenced by angularWidth(), and CSCRadialStripTopology().

Definition at line 241 of file CSCRadialStripTopology.h.

Referenced by centreToIntersection(), and originToIntersection().

Definition at line 240 of file CSCRadialStripTopology.h.

Referenced by detHeight(), stripLength(), and yExtentOfStripPlane().

Definition at line 238 of file CSCRadialStripTopology.h.

Referenced by channel(), CSCRadialStripTopology(), and nstrips().

Definition at line 242 of file CSCRadialStripTopology.h.

Referenced by CSCRadialStripTopology(), and phiOfOneEdge().

Definition at line 243 of file CSCRadialStripTopology.h.

Referenced by CSCRadialStripTopology(), and yAxisOrientation().