CMS 3D CMS Logo

RadialStripTopology Class Reference

A StripTopology in which the component strips subtend a constant angular width, and, if projected, intersect at a point. More...

#include <Geometry/CommonTopologies/interface/RadialStripTopology.h>

Inheritance diagram for RadialStripTopology:

StripTopology Topology OffsetRadialStripTopology CSCStripTopology CSCGangedStripTopology CSCUngangedStripTopology

List of all members.

Public Member Functions

float angularWidth () const
 Angular width of a each strip.
float centreToIntersection () const
 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.
virtual int channel (const LocalPoint &) const
 Channel number corresponding to a given LocalPoint.
float detHeight () const
 Length of long symmetry axis of plane of strips.
virtual LocalError localError (const MeasurementPoint &, const MeasurementError &) const
 LocalError for a given MeasurementPoint with known MeasurementError.
virtual LocalError localError (float strip, float stripErr2) const
 LocalError for a pure strip measurement, where 'strip' is the (float) position (a 'phi' angle wrt y axis) and stripErr2 is the sigma-squared.
virtual float localPitch (const LocalPoint &) const
 Pitch (strip width) at a given LocalPoint.
virtual LocalPoint localPosition (const MeasurementPoint &) const
 LocalPoint for a given MeasurementPoint
What's a MeasurementPoint?
In analogy with that used with TrapezoidalStripTopology objects, a MeasurementPoint is a 2-dim object.
virtual LocalPoint localPosition (float strip) const
 LocalPoint on x axis for given 'strip' 'strip' is a float in units of the strip (angular) width.
virtual float localStripLength (const LocalPoint &) const
 Length of a strip passing through a given LocalPpoint.
virtual MeasurementError measurementError (const LocalPoint &, const LocalError &) const
virtual MeasurementPoint measurementPosition (const LocalPoint &) const
virtual int nearestStrip (const LocalPoint &) const
 Nearest strip to given LocalPoint.
virtual int nstrips () const
 Total number of strips.
float originToIntersection () const
 (y) distance from intersection of the projections of the strips to the local coordinate origin.
float phiOfOneEdge () const
 Convenience function to access azimuthal angle of extreme edge of first strip measured relative to long symmetry axis of the plane of strips.
virtual float phiPitch (void) const
 Phi pitch of each strip (= angular width!).
virtual float pitch () const
 BEWARE: calling pitch() throws an exception.
 RadialStripTopology (int ns, float aw, float dh, float r, int yAx=1, float yMid=0.)
 Constructor from:.
virtual float strip (const LocalPoint &) const
 Strip in which a given LocalPoint lies.
virtual float stripAngle (float strip) const
 Angle between strip and symmetry axis (=local y axis) for given strip.
virtual float stripLength () const
 Height of detector (= length of long symmetry axis of the plane of strips).
float xOfStrip (int strip, float y) const
 Local x where centre of strip intersects input local y
'strip' should be in range 1 to nstrips()
.
int yAxisOrientation () const
 y axis orientation, 1 means detector width increases with local y
float yCentreOfStripPlane () const
 Offset in local y between midpoint of detector (strip plane) extent and local origin.
float yDistanceToIntersection (float y) const
 Distance in local y from a hit to the point of intersection of projected strips.
float yExtentOfStripPlane () const
 y extent of strip plane
virtual ~RadialStripTopology ()
 Destructor.

Private Attributes

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

Friends

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


Detailed Description

A StripTopology in which the component strips subtend a constant angular width, and, if projected, intersect at a point.

Author:
Tim Cox
WARNING! Wherever 'float strip' is used the units of 'strip' are angular widths of each strip. The range is from 0.0 at the extreme edge of the 'first' strip at one edge of the detector, to nstrip*angular width at the other edge.
The centre of the first strip is at strip = 0.5
The centre of the last strip is at strip = 0.5 + (nstrip-1)
This is for consistency with CommonDet usage of 'float strip' (but where units are strip pitch rather than strip angular width.)

WARNING! If the mid-point along local y of the plane of strips does not correspond to the local coordinate origin, set the final ctor argument appropriately.

Definition at line 26 of file RadialStripTopology.h.


Constructor & Destructor Documentation

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

Constructor from:.

Parameters:
ns number of strips
aw angular width of a strip
dh detector height (usually 2 x apothem of TrapezoidalPlaneBounds)
r radial distance from symmetry centre of detector to the point at which the outer edges of the two extreme strips (projected) intersect.
yAx orientation 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.
yMid local 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 9 of file RadialStripTopology.cc.

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

00009                                                                                                  :
00010   theNumberOfStrips(ns), theAngularWidth(aw), 
00011   theDetHeight(dh), theCentreToIntersection(r),
00012   theYAxisOrientation(yAx), yCentre( yMid) {   
00013   // Angular offset of extreme edge of detector, so that angle is
00014   // zero for a strip lying along local y axis = long symmetry axis of plane of strips
00015   thePhiOfOneEdge = -(theNumberOfStrips/2.) * theAngularWidth * yAx;
00016   
00017   LogTrace("RadialStripTopology") << "RadialStripTopology: constructed with"
00018         << " strips = " << ns
00019         << " width = " << aw << " rad "
00020         << " det_height = " << dh
00021         << " ctoi = " << r 
00022         << " phi_edge = " << thePhiOfOneEdge << " rad "
00023         << " y_ax_ori = " << theYAxisOrientation
00024         << " y_det_centre = " << yCentre 
00025         << "\n";
00026 }    

virtual RadialStripTopology::~RadialStripTopology (  )  [inline, virtual]

Destructor.

Definition at line 48 of file RadialStripTopology.h.

00048 {}


Member Function Documentation

float RadialStripTopology::angularWidth (  )  const [inline]

Angular width of a each strip.

Definition at line 169 of file RadialStripTopology.h.

References theAngularWidth.

Referenced by localError(), localPitch(), OffsetRadialStripTopology::localPosition(), measurementError(), measurementPosition(), OffsetRadialStripTopology::OffsetRadialStripTopology(), phiPitch(), OffsetRadialStripTopology::strip(), strip(), OffsetRadialStripTopology::stripAngle(), and stripAngle().

00169 { return theAngularWidth;}

float RadialStripTopology::centreToIntersection (  )  const [inline]

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.

Definition at line 191 of file RadialStripTopology.h.

References theCentreToIntersection.

Referenced by localError(), and operator<<().

00191 { return theCentreToIntersection; }

int RadialStripTopology::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 Topology.

Reimplemented in OffsetRadialStripTopology, CSCGangedStripTopology, and CSCUngangedStripTopology.

Definition at line 182 of file RadialStripTopology.cc.

References min, strip(), and theNumberOfStrips.

Referenced by CSCUngangedStripTopology::channel().

00182                                                        {
00183   return std::min( int( strip(lp) ), theNumberOfStrips-1 );
00184 }

float RadialStripTopology::detHeight (  )  const [inline]

Length of long symmetry axis of plane of strips.

Definition at line 179 of file RadialStripTopology.h.

References theDetHeight.

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

00179 { return theDetHeight;}

LocalError RadialStripTopology::localError ( const MeasurementPoint mp,
const MeasurementError merr 
) 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 StripTopology.

Definition at line 83 of file RadialStripTopology.cc.

References funct::A, angularWidth(), c1, c2, centreToIntersection(), funct::D, detHeight(), L, phi, s1, s2, funct::sin(), funct::sqrt(), stripAngle(), MeasurementError::uu(), MeasurementError::uv(), MeasurementError::vv(), PV2DBase< T, PVType, FrameType >::x(), PV2DBase< T, PVType, FrameType >::y(), and yAxisOrientation().

00084                                                                     {
00085   // Here we need to allow the possibility of correlated errors, since
00086   // that may happen during Kalman filtering
00087 
00088   //  float phi = phiOfOneEdge() + yAxisOrientation() * mp.x() * angularWidth();
00089   float phi = stripAngle( mp.x() );
00090 
00091   //  float t = tan( phi );                        // tan(angle between strip and y)
00092   //float c2 = 1./(1. + t*t);                    // cos(angle)**2
00093   // float cs = t*c2;                             // sin(angle)*cos(angle); tan carries sign of sin!
00094   // float s2 = t*t * c2;                         // sin(angle)**2
00095   
00096 
00097   float s1 = sin(phi);
00098   float s2 = s1*s1;
00099   float c2 = 1.-s2;
00100   float c1 = std::sqrt(c2);
00101   float cs = c1*s1;
00102 
00103   float A  = angularWidth();
00104   float A2 = A * A;
00105 
00106   // D is distance from intersection of edges to hit on strip
00107   float D = (centreToIntersection() + yAxisOrientation() * mp.y() * detHeight()) /c1;
00108   float D2 = D * D;
00109 
00110   // L is length of strip across face of chamber
00111   float L = detHeight()/ c1;  
00112   float L2  = L*L; 
00113 
00114   // MeasurementError elements are already squared
00115   // but they're normalized to products of A and L 
00116   // (N.B. uses L=length of strip, and not D=distance to intersection)
00117   // Remember to ensure measurement error components are indeed normalized like this!
00118 
00119   float SA2 = merr.uu() * A2;
00120   float SD2 = merr.vv() * L2; // Note this norm uses stripLength**2
00121   float RHOSASR = merr.uv() * A * L;
00122 
00123   float sx2 = SA2 * D2 * c2  +  2. * RHOSASR * D * cs  +  SD2 * s2;
00124   float sy2 = SA2 * D2 * s2  -  2. * RHOSASR * D * cs  +  SD2 * c2;
00125   float rhosxsy = cs * ( SD2 - D2 * SA2 )  +  RHOSASR * D * ( c2 - s2 );
00126 
00127   return LocalError(sx2, rhosxsy, sy2);
00128 }

LocalError RadialStripTopology::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 StripTopology.

Definition at line 55 of file RadialStripTopology.cc.

References angularWidth(), c2, centreToIntersection(), detHeight(), s2, stripAngle(), t, and funct::tan().

00055                                                                   {
00056   // Consider measurement as strip (phi) position and mid-point of strip
00057   // Since 'strip' is in units of angular strip-widths, stripErr2 is
00058   // required to be in corresponding units.
00059 
00060   float t = std::tan( stripAngle( strip ) );        // tan(angle between strip and y)
00061   float c2 = 1./(1. + t*t);                    // cos(angle)**2
00062   float cs = t*c2;                             // sin(angle)*cos(angle); tan carries sign of sin!
00063   float s2 = t*t * c2;                         // sin(angle)**2
00064 
00065   float D2 = centreToIntersection()*centreToIntersection() / c2; // (mid pt of strip to intersection)**2
00066   float L2 = detHeight()*detHeight() / c2;   // length**2 of strip across detector
00067   float A2 = angularWidth()*angularWidth();
00068 
00069   // The error**2 we're assigning is L2/12 wherever we measure the position along the strip
00070   // from... we just know the measurement is somewhere ON the strip. 
00071 
00072   float SD2 = L2 / 12.;       // SD = Sigma-Distance-along-strip
00073   float SA2 = A2 * stripErr2; // SA = Sigma-Angle-of-strip
00074 
00075   float sx2 = c2 * D2 * SA2 + s2 * SD2;
00076   float sy2 = s2 * D2 * SA2 + c2 * SD2;
00077   float rhosxsy = cs * ( SD2 - D2 * SA2 );
00078 
00079   return LocalError(sx2, rhosxsy, sy2);
00080 }

float RadialStripTopology::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 StripTopology.

Definition at line 197 of file RadialStripTopology.cc.

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

Referenced by TrackerValidationVariables::fillHitQuantities(), and CSCLayerGeometry::stripPitch().

00197                                                           {
00198   // The local pitch is the local x width of the strip at the local (x,y)
00199 
00200   // Calculating it is a nightmare...deriving the expression below is left
00201   // as a exercise for the reader. 
00202 
00203   float fstrip = strip(lp); // position in strip units
00204   int istrip = static_cast<int>(fstrip) + 1; // which strip number
00205   if (istrip>nstrips() ) istrip = nstrips(); // enforce maximum
00206   float fangle = stripAngle(static_cast<float>(istrip) - 0.5); // angle of strip centre
00207   float localp = yDistanceToIntersection( lp.y() ) * sin(angularWidth()) /
00208     ( cos(fangle-0.5*angularWidth())*cos(fangle+0.5*angularWidth()) );
00209   return localp;
00210 }

LocalPoint RadialStripTopology::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 StripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 41 of file RadialStripTopology.cc.

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

00041                                                                    {
00042   float phi = stripAngle( mp.x() );
00043   // 2nd component of MP is fractional position along strip, range -0.5 to +0.5 in direction of y
00044   // so distance along strip (measured from mid-point of length of strip) is
00045   //     mp.y() * (length of strip). 
00046   // Thus distance along y direction, from mid-point of strip, is 
00047   //    mp.y() * (length of strip) * cos(phi)
00048   // But (length of strip) = theDetHeight/cos(phi), so
00049   float y =  mp.y() * detHeight() + yCentreOfStripPlane();
00050   float x = yAxisOrientation() * yDistanceToIntersection( y ) * tan ( phi );
00051   return LocalPoint( x, y );
00052 }

LocalPoint RadialStripTopology::localPosition ( float  strip  )  const [virtual]

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

Implements StripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 36 of file RadialStripTopology.cc.

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

Referenced by OffsetRadialStripTopology::localPosition().

00036                                                     {
00037   return LocalPoint( yAxisOrientation() * originToIntersection() * tan( stripAngle(strip) ), 0.0 );
00038 }

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

Length of a strip passing through a given LocalPpoint.

Implements StripTopology.

Definition at line 217 of file RadialStripTopology.cc.

References c2, detHeight(), funct::sqrt(), t, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and yDistanceToIntersection().

00217                                                                       {
00218   float yHitToInter = yDistanceToIntersection( lp.y() );
00219   // since we're dealing with magnitudes, sign is unimportant
00220   float t  = lp.x() / yHitToInter;    // tan(angle between strip and y)
00221   float c2 = 1./(1. + t*t);           // cos(angle)**2
00222   return detHeight() / sqrt(c2);
00223 }

MeasurementError RadialStripTopology::measurementError ( const LocalPoint lp,
const LocalError lerr 
) const [virtual]

Implements Topology.

Definition at line 152 of file RadialStripTopology.cc.

References funct::A, angularWidth(), c2, funct::D, detHeight(), L, s2, funct::sqrt(), t, PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), LocalError::xy(), PV3DBase< T, PVType, FrameType >::y(), yAxisOrientation(), yDistanceToIntersection(), and LocalError::yy().

Referenced by TrackerValidationVariables::fillHitQuantities().

00153                                 {
00154 
00155   float yHitToInter = yDistanceToIntersection( lp.y() );
00156   float t  = yAxisOrientation() * lp.x() / yHitToInter; // tan(angle between strip and y) 
00157   float c2 = 1./(1. + t*t);  // cos(angle)**2
00158   float cs = t*c2;           // sin(angle)*cos(angle); tan carries sign of sin!
00159   float s2 = 1. - c2;        // sin(angle)**2
00160 
00161   float A  = angularWidth();
00162 
00163   // L is length of strip across face of chamber
00164   float L2 = detHeight()*detHeight() / c2;  
00165   float L  = std::sqrt(L2); 
00166 
00167   // D is distance from intersection of edges to hit on strip
00168   float D2 = lp.x()*lp.x() + yHitToInter*yHitToInter;
00169   float D = std::sqrt(D2);
00170 
00171   float LP = D * A;
00172   float LP2 = LP * LP;
00173 
00174   float SA2 = ( c2 * lerr.xx() - 2. * cs * lerr.xy() + s2 * lerr.yy() ) / LP2;
00175   float SD2 = ( s2 * lerr.xx() + 2. * cs * lerr.xy() + c2 * lerr.yy() ) / L2;
00176   float RHOSASR = ( cs * ( lerr.xx() - lerr.yy() ) + ( c2 - s2 ) * lerr.xy() ) / (LP*L);
00177 
00178   return MeasurementError(SA2, RHOSASR, SD2);
00179 }

MeasurementPoint RadialStripTopology::measurementPosition ( const LocalPoint lp  )  const [virtual]

Implements Topology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 143 of file RadialStripTopology.cc.

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

Referenced by TrackerValidationVariables::fillHitQuantities(), and OffsetRadialStripTopology::measurementPosition().

00143                                                                    {
00144   // Note that this phi is (pi/2 - conventional local phi)
00145   // This means use atan2(x,y) rather than more usual atan2(y,x)
00146   float phi = yAxisOrientation() * atan2( lp.x(), yDistanceToIntersection( lp.y() ) );
00147   return MeasurementPoint( yAxisOrientation()*( phi-phiOfOneEdge() )/angularWidth(),
00148                           (lp.y() - yCentreOfStripPlane())/detHeight() );
00149 }

int RadialStripTopology::nearestStrip ( const LocalPoint lp  )  const [virtual]

Nearest strip to given LocalPoint.

Definition at line 225 of file RadialStripTopology.cc.

References nstrips(), and strip().

Referenced by CSCLayerGeometry::nearestStrip().

00226 {
00227   // xxxStripTopology::strip() is expected to have range 0. to 
00228   // float(no of strips), but be extra careful and enforce that
00229 
00230   float fstrip = this->strip(lp);
00231   int n_strips = this->nstrips();
00232   fstrip = ( fstrip>=0. ? fstrip :  0. ); // enforce minimum 0.
00233   int near = static_cast<int>( fstrip ) + 1; // first strip is 1
00234   near = ( near<=n_strips ? near : n_strips ); // enforce maximum at right edge
00235   return near;
00236 }

virtual int RadialStripTopology::nstrips (  )  const [inline, virtual]

Total number of strips.

Implements StripTopology.

Definition at line 130 of file RadialStripTopology.h.

References theNumberOfStrips.

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

00130 { return theNumberOfStrips; }

float RadialStripTopology::originToIntersection (  )  const [inline]

(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.

Definition at line 198 of file RadialStripTopology.h.

References theCentreToIntersection, and yCentre.

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

00198 { return (theCentreToIntersection - yCentre); }

float RadialStripTopology::phiOfOneEdge (  )  const [inline]

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().

Definition at line 212 of file RadialStripTopology.h.

References thePhiOfOneEdge.

Referenced by OffsetRadialStripTopology::localPosition(), measurementPosition(), operator<<(), OffsetRadialStripTopology::strip(), strip(), OffsetRadialStripTopology::stripAngle(), and stripAngle().

00212 { return thePhiOfOneEdge; }

virtual float RadialStripTopology::phiPitch ( void   )  const [inline, virtual]

Phi pitch of each strip (= angular width!).

Definition at line 174 of file RadialStripTopology.h.

References angularWidth().

Referenced by operator<<(), and CSCLayerGeometry::stripPhiPitch().

00174 { return angularWidth(); }

float RadialStripTopology::pitch (  )  const [virtual]

BEWARE: calling pitch() throws an exception.


Pitch is conventional name for width of something, but this is not sensible for a RadialStripTopology since strip widths vary with local y. Use localPitch(.) instead.

Implements StripTopology.

Definition at line 187 of file RadialStripTopology.cc.

00187                                  { 
00188   // BEWARE: this originally returned the pitch at the local origin 
00189   // but Tracker prefers throwing an exception since the strip is 
00190   // not constant along local x axis for a RadialStripTopology.
00191 
00192   throw Genexception("RadialStripTopology::pitch() called - makes no sense, use localPitch(.) instead.");
00193   return 0.;
00194 }

float RadialStripTopology::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 StripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 131 of file RadialStripTopology.cc.

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

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

00131                                                      {
00132   // Note that this phi is measured from y axis, so sign of angle <-> sign of x * yAxisOrientation
00133   // Use atan2(x,y) rather than more usual atan2(y,x)
00134 
00135   float phi = atan2( lp.x(), yDistanceToIntersection( lp.y() ) );
00136   float aStrip = ( phi - yAxisOrientation() * phiOfOneEdge() )/angularWidth();
00137   if (aStrip < 0. ) aStrip = 0.;
00138   else if (aStrip > theNumberOfStrips)  aStrip = theNumberOfStrips;
00139   return aStrip;
00140 }

float RadialStripTopology::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 StripTopology.

Reimplemented in OffsetRadialStripTopology.

Definition at line 213 of file RadialStripTopology.cc.

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

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

00213                                                  {
00214   return ( phiOfOneEdge() + yAxisOrientation() * strip * angularWidth() );
00215 }

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

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

Implements StripTopology.

Definition at line 135 of file RadialStripTopology.h.

References theDetHeight.

00135 { return theDetHeight; }

float RadialStripTopology::xOfStrip ( int  strip,
float  y 
) const

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

Definition at line 29 of file RadialStripTopology.cc.

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

Referenced by CSCLayerGeometry::xOfStrip().

00029                                                       {
00030   // Expect input 'strip' to be in range 1 to nstrips()
00031   float tanPhi = std::tan( stripAngle(static_cast<float>(strip) - 0.5 ) );
00032   return yAxisOrientation()* yDistanceToIntersection( y ) * tanPhi;
00033 }

int RadialStripTopology::yAxisOrientation (  )  const [inline]

y axis orientation, 1 means detector width increases with local y

Definition at line 228 of file RadialStripTopology.h.

References theYAxisOrientation.

Referenced by localError(), localPosition(), measurementError(), measurementPosition(), operator<<(), strip(), stripAngle(), xOfStrip(), and yDistanceToIntersection().

00228 { return theYAxisOrientation; }

float RadialStripTopology::yCentreOfStripPlane (  )  const [inline]

Offset in local y between midpoint of detector (strip plane) extent and local origin.

Definition at line 233 of file RadialStripTopology.h.

References yCentre.

Referenced by OffsetRadialStripTopology::localPosition(), localPosition(), measurementPosition(), operator<<(), and CSCStripTopology::yLimitsOfStripPlane().

00233 { return yCentre; }

float RadialStripTopology::yDistanceToIntersection ( float  y  )  const

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

Definition at line 238 of file RadialStripTopology.cc.

References originToIntersection(), and yAxisOrientation().

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

00238                                                                   {
00239   return yAxisOrientation() * y + originToIntersection();
00240 }

float RadialStripTopology::yExtentOfStripPlane (  )  const [inline]

y extent of strip plane

Definition at line 184 of file RadialStripTopology.h.

References theDetHeight.

Referenced by CSCStripTopology::yLimitsOfStripPlane().

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


Friends And Related Function Documentation

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

Definition at line 242 of file RadialStripTopology.cc.

00243 {
00244   os  << "RadialStripTopology " << std::endl
00245       << " " << std::endl
00246       << "number of strips          " << rst.nstrips() << std::endl
00247       << "centre to whereStripsMeet " << rst.centreToIntersection() << std::endl
00248       << "detector height in y      " << rst.detHeight() << std::endl
00249       << "angular width of strips   " << rst.phiPitch() << std::endl
00250       << "phi of one edge           " << rst.phiOfOneEdge() << std::endl
00251       << "y axis orientation        " << rst.yAxisOrientation() << std::endl
00252       << "y of centre of strip plane " << rst.yCentreOfStripPlane() << std::endl;
00253   return os;
00254 }


Member Data Documentation

float RadialStripTopology::theAngularWidth [private]

Definition at line 245 of file RadialStripTopology.h.

Referenced by angularWidth(), and RadialStripTopology().

float RadialStripTopology::theCentreToIntersection [private]

Definition at line 247 of file RadialStripTopology.h.

Referenced by centreToIntersection(), and originToIntersection().

float RadialStripTopology::theDetHeight [private]

Definition at line 246 of file RadialStripTopology.h.

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

int RadialStripTopology::theNumberOfStrips [private]

Definition at line 244 of file RadialStripTopology.h.

Referenced by channel(), nstrips(), RadialStripTopology(), and strip().

float RadialStripTopology::thePhiOfOneEdge [private]

Definition at line 248 of file RadialStripTopology.h.

Referenced by phiOfOneEdge(), and RadialStripTopology().

int RadialStripTopology::theYAxisOrientation [private]

Definition at line 249 of file RadialStripTopology.h.

Referenced by RadialStripTopology(), and yAxisOrientation().

float RadialStripTopology::yCentre [private]

Definition at line 250 of file RadialStripTopology.h.

Referenced by originToIntersection(), RadialStripTopology(), and yCentreOfStripPlane().


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