CMS 3D CMS Logo

Public Member Functions | Private Attributes | Friends

CSCLayerGeometry Class Reference

#include <CSCLayerGeometry.h>

Inheritance diagram for CSCLayerGeometry:
TrapezoidalPlaneBounds Bounds

List of all members.

Public Member Functions

int channel (int strip) const
virtual Boundsclone () const
 CSCLayerGeometry (const CSCLayerGeometry &)
 CSCLayerGeometry (const CSCGeometry *geom, int iChamberType, const TrapezoidalPlaneBounds &bounds, int nstrips, float stripOffset, float stripPhiPitch, float whereStripsMeet, float extentOfStripPlane, float yCentreOfStripPlane, const CSCWireGroupPackage &wg, float wireAngleInDegrees, double yOfFirstWire, float hThickness)
bool inside (const Local3DPoint &, const LocalError &, float scale=1.f) const
bool inside (const Local3DPoint &) const
 Determine if the point is inside the bounds.
bool inside (const Local2DPoint &) const
LocalPoint intersectionOfStripAndWire (float s, int w) const
LocalPoint intersectionOfTwoLines (std::pair< float, float > p1, std::pair< float, float > p2) const
float lengthOfWireGroup (int wireGroup) const
LocalPoint localCenterOfWireGroup (int wireGroup) const
LocalError localError (int strip, float sigmaStrip, float sigmaWire) const
float middleWireOfGroup (int wireGroup) const
int nearestStrip (const LocalPoint &lp) const
int nearestWire (const LocalPoint &lp) const
int numberOfStrips () const
int numberOfWireGroups () const
int numberOfWires () const
int numberOfWiresPerGroup (int wireGroup) const
CSCLayerGeometryoperator= (const CSCLayerGeometry &)
std::pair< LocalPoint, float > possibleRecHitPosition (float s, int w1, int w2) const
void setTopology (CSCStripTopology *topology)
int stagger () const
float strip (const LocalPoint &lp) const
float stripAngle (int strip) const
float stripOffset (void) const
float stripPhiPitch () const
float stripPitch (const LocalPoint &lp) const
float stripPitch () const
LocalPoint stripWireGroupIntersection (int strip, int wireGroup) const
LocalPoint stripWireIntersection (int strip, float wire) const
const CSCStripTopologytopology () const
float wireAngle () const
int wireGroup (int wire) const
float wirePitch () const
const CSCWireTopologywireTopology () const
float xOfStrip (int strip, float y=0.) const
std::pair< float, float > yLimitsOfStripPlane () const
float yOfWire (float wire, float x=0.) const
float yOfWireGroup (int wireGroup, float x=0.) const
float yResolution (int wireGroup=1) const
virtual ~CSCLayerGeometry ()

Private Attributes

float apothem
int chamberType
float hBottomEdge
float hTopEdge
const std::string myName
CSCStripTopologytheStripTopology
CSCWireTopologytheWireTopology

Friends

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

Detailed Description

Encapsulates the geometrical details of a CSCLayer in a WireTopology for the wires and in a StripTopology for the strips. Note that it does not have the capability of calculating global values, so all values are in local coordinates.

Author:
Tim Cox

Definition at line 25 of file CSCLayerGeometry.h.


Constructor & Destructor Documentation

CSCLayerGeometry::CSCLayerGeometry ( const CSCGeometry geom,
int  iChamberType,
const TrapezoidalPlaneBounds bounds,
int  nstrips,
float  stripOffset,
float  stripPhiPitch,
float  whereStripsMeet,
float  extentOfStripPlane,
float  yCentreOfStripPlane,
const CSCWireGroupPackage wg,
float  wireAngleInDegrees,
double  yOfFirstWire,
float  hThickness 
)

Ctor from basic trapezoidal Chamber parameters.

Parameters:
geomThe pointer to the actual CSCGeometry we're building.
iChamberTypeThe index 1-9 for station/ring combination.
TrapezoidalPlaneBoundsdescribing geometry of face.
nstripsNo. of strips in cathode plane of a Layer.
stripOffsetAlternate strip planes are relatively shifted by +/-0.25 strip widths.
stripPhiPitchDelta-phi width of strips (they're fan-shaped) in radians
whereStripsMeetradial distance from projected intersection of strips to centre of strip plane
extentOfStripPlaneheight of strip plane (along its long symmetry axis)
yCentreOfStripPlanelocal y of symmetry centre of strip plane (before any offset rotation)
wgCSCWireGroupPackage encapsulating wire group info.
wireAngleInDegreesangle of wires w.r.t local x axis.
yOfFirstWirelocal y coordinate of first (lowest) wire in wire plane - nearest narrow edge.
hThicknesshalf-thickness of chamber layer in cm (i.e. half the gas gap).

Definition at line 22 of file CSCLayerGeometry.cc.

References apothem, funct::cos(), CSCGeometry::gangedStrips(), hBottomEdge, LogTrace, myName, CSCGeometry::realWireGeometry(), funct::sin(), theStripTopology, theWireTopology, and CSCWireGroupPackage::wireSpacing.

Referenced by clone().

  :   TrapezoidalPlaneBounds( bounds.widthAtHalfLength() - bounds.width()/2., bounds.width()/2., bounds.length()/2., hThickness ), 
      theWireTopology( 0 ), theStripTopology( 0 ), 
      hBottomEdge( bounds.widthAtHalfLength() - bounds.width()/2. ), 
      hTopEdge( bounds.width()/2. ), apothem( bounds.length()/2. ),  
      myName( "CSCLayerGeometry" ), chamberType( iChamberType ) {

  LogTrace("CSCLayerGeometry|CSC") << myName <<": being constructed, this=" << this;

  // Ganged strips in ME1A?
  bool gangedME1A = ( iChamberType == 1 && geom->gangedStrips() );

  CSCStripTopology* aStripTopology = 
        new CSCUngangedStripTopology(nstrips, stripPhiPitch,
            extentOfStripPlane, whereStripsMeet, stripOffset, yCentreOfStripPlane );

  if ( gangedME1A ) {
    theStripTopology = new CSCGangedStripTopology( *aStripTopology, 16 );
    delete aStripTopology;
  }
  else {
    theStripTopology = aStripTopology;
  }

  if ( ! geom->realWireGeometry() ) {
    // Approximate ORCA_8_8_0 and earlier calculated geometry...
    float wangler = wireAngleInDegrees*degree; // convert angle to radians
    float wireCos = cos(wangler);
    float wireSin = sin(wangler);
    float y2 = apothem * wireCos + hBottomEdge * fabs(wireSin);
    float wireSpacing = wg.wireSpacing/10.; // in cm
    float wireOffset = -y2 + wireSpacing/2.;
    yOfFirstWire = wireOffset/wireCos;
  }

  theWireTopology = new CSCWireTopology( wg, yOfFirstWire, wireAngleInDegrees );
  LogTrace("CSCLayerGeometry|CSC") << myName <<": constructed: "<< *this;
} 
CSCLayerGeometry::CSCLayerGeometry ( const CSCLayerGeometry melg)

Definition at line 65 of file CSCLayerGeometry.cc.

References CSCStripTopology::clone(), theStripTopology, and theWireTopology.

                                                               :
  TrapezoidalPlaneBounds(melg.hBottomEdge, melg.hTopEdge, melg.apothem,
                         0.5 * melg.thickness() ),
  theWireTopology(0), theStripTopology(0), 
  hBottomEdge(melg.hBottomEdge), hTopEdge(melg.hTopEdge),
  apothem(melg.apothem), chamberType(melg.chamberType) 
{
  // CSCStripTopology is abstract, so need clone()
  if (melg.theStripTopology) theStripTopology = melg.theStripTopology->clone();
  // CSCWireTopology is concrete, so direct copy
  if (melg.theWireTopology) theWireTopology = new CSCWireTopology(*(melg.theWireTopology));
}
CSCLayerGeometry::~CSCLayerGeometry ( ) [virtual]

Definition at line 100 of file CSCLayerGeometry.cc.

References LogTrace, myName, theStripTopology, and theWireTopology.

{
  LogTrace("CSCLayerGeometry|CSC") << myName << ": being destroyed, this=" << this << 
    "\nDeleting theStripTopology=" << theStripTopology << 
    " and theWireTopology=" << theWireTopology;
  delete theStripTopology;
  delete theWireTopology;
}

Member Function Documentation

int CSCLayerGeometry::channel ( int  strip) const [inline]

Electronics channel corresponding to a given strip ...sometimes there will be more than one strip OR'ed into a channel

Definition at line 115 of file CSCLayerGeometry.h.

References OffsetRadialStripTopology::channel(), and theStripTopology.

Referenced by MuonTruth::analyze(), MuonTruth::associateCSCHitId(), MuonTruth::associateHitId(), and CSCStripElectronicsSim::readoutElement().

virtual Bounds* CSCLayerGeometry::clone ( void  ) const [inline, virtual]

Utility method to handle proper copying of the class

Reimplemented from TrapezoidalPlaneBounds.

Definition at line 306 of file CSCLayerGeometry.h.

References CSCLayerGeometry().

                                { 
    return new CSCLayerGeometry(*this);
  }
bool CSCLayerGeometry::inside ( const Local3DPoint lp,
const LocalError le,
float  scale = 1.f 
) const [virtual]

Is a supplied LocalPoint inside the strip region?

This is a more reliable fiducial cut for CSCs than the 'Bounds' of the GeomDet(Unit) since those ranges are looser than the sensitive gas region. There are three versions, to parallel those of the TrapezoidalPlaneBounds which a naive user might otherwise employ.

Reimplemented from TrapezoidalPlaneBounds.

Definition at line 283 of file CSCLayerGeometry.cc.

References mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), LocalError::yy(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by CSCWireHitSim::getIonizationClusters(), CSCMake2DRecHit::hitFromStripAndWire(), inside(), CSCMake2DRecHit::isHitInFiducial(), and MuonSimHitProducer::produce().

                                                                                               {
  // Effectively consider that the LocalError components extend the area which is acceptable.
  // Form a little box centered on the point, with x, y diameters defined by the errors
  // and require that ALL four corners of the box fall outside the strip region for failure

  // Note that LocalError is 2-dim x,y and doesn't supply a z error
  float deltaX = scale*sqrt(le.xx());
  float deltaY = scale*sqrt(le.yy());

  LocalPoint lp1( lp.x()-deltaX, lp.y()-deltaY, lp.z() );
  LocalPoint lp2( lp.x()-deltaX, lp.y()+deltaY, lp.z() );
  LocalPoint lp3( lp.x()+deltaX, lp.y()+deltaY, lp.z() );
  LocalPoint lp4( lp.x()+deltaX, lp.y()-deltaY, lp.z() );

  return ( inside(lp1) || inside(lp2) || inside(lp3) || inside(lp4) );
}
bool CSCLayerGeometry::inside ( const Local3DPoint ) const [virtual]

Determine if the point is inside the bounds.

Reimplemented from TrapezoidalPlaneBounds.

Definition at line 263 of file CSCLayerGeometry.cc.

References epsilon, numberOfStrips(), query::result, OffsetRadialStripTopology::strip(), theStripTopology, TrapezoidalPlaneBounds::thickness(), PV3DBase< T, PVType, FrameType >::y(), yLimitsOfStripPlane(), and PV3DBase< T, PVType, FrameType >::z().

                                                            {
  bool result = false;
  const float epsilon = 1.e-06;
  if ( fabs( lp.z() ) < thickness()/2. ) { // thickness of TPB is that of gas layer
    std::pair<float, float> ylims = yLimitsOfStripPlane();
    if ( (lp.y() > ylims.first) && (lp.y() < ylims.second) ) {
      // 'strip' returns float value between 0. and float(Nstrips) and value outside
      // is set to 0. or float(Nstrips)... add a conservative precision of 'epsilon'
      if ( ( theStripTopology->strip(lp) > epsilon ) && 
           ( theStripTopology->strip(lp) < (numberOfStrips() - epsilon) ) ) result = true;
    }
  }
  return result;
}
bool CSCLayerGeometry::inside ( const Local2DPoint lp) const [virtual]

Reimplemented from TrapezoidalPlaneBounds.

Definition at line 278 of file CSCLayerGeometry.cc.

References inside(), PV2DBase< T, PVType, FrameType >::x(), and PV2DBase< T, PVType, FrameType >::y().

                                                            {
  LocalPoint lp2( lp.x(), lp.y(), 0. );
  return inside( lp2 );
}
LocalPoint CSCLayerGeometry::intersectionOfStripAndWire ( float  s,
int  w 
) const

Return 2-dim point at which a strip and a wire intersects.

Input arguments: a (float) strip number, and an (int) wire.
Output: LocalPoint which is at their intersection, or at extreme y of wire plane, as appropriate. (If y is adjusted, x is also adjusted to keep it on same strip.)

Definition at line 197 of file CSCLayerGeometry.cc.

References CSCStripTopology::equationOfStrip(), CSCWireTopology::equationOfWire(), CSCWireTopology::insideYOfWirePlane(), intersectionOfTwoLines(), CSCWireTopology::restrictToYOfWirePlane(), OffsetRadialStripTopology::stripAngle(), funct::tan(), theStripTopology, theWireTopology, PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), and detailsBasic3DVector::y.

Referenced by possibleRecHitPosition().

                                                                             {
        
  std::pair<float, float> pw = theWireTopology->equationOfWire( static_cast<float>(w) );
  std::pair<float, float> ps = theStripTopology->equationOfStrip( s );
  LocalPoint sw = intersectionOfTwoLines( ps, pw );
        
  // If point falls outside wire plane, at extremes in local y, 
  // replace its y by that of appropriate edge of wire plane
  if ( !(theWireTopology->insideYOfWirePlane( sw.y() ) ) ) {
     float y  = theWireTopology->restrictToYOfWirePlane( sw.y() );
     // and adjust x to match new y
     float x = sw.x() + (y - sw.y())*tan(theStripTopology->stripAngle(s));
     sw = LocalPoint(x, y);
  }
        
  return sw;
}
LocalPoint CSCLayerGeometry::intersectionOfTwoLines ( std::pair< float, float >  p1,
std::pair< float, float >  p2 
) const

Return the point of intersection of two straight lines (in 2-dim).

Input arguments are pair(m1,c1) and pair(m2,c2) where m=slope, c=intercept (y=mx+c).
BEWARE! Do not call with m1 = m2 ! No trapping !

Definition at line 215 of file CSCLayerGeometry.cc.

References alignmentValidation::c1, x, and detailsBasic3DVector::y.

Referenced by intersectionOfStripAndWire().

                                                                                                              {

  // Calculate the point of intersection of two straight lines (in 2-dim)
  // input arguments are pair(m1,c1) and pair(m2,c2) where m=slope, c=intercept (y=mx+c)
  // BEWARE! Do not call with m1 = m2 ! No trapping !

  float m1 = p1.first;
  float c1 = p1.second;
  float m2 = p2.first;
  float c2 = p2.second;
  float x = (c2-c1)/(m1-m2);
  float y = (m1*c2-m2*c1)/(m1-m2);
  return LocalPoint( x, y );
}
float CSCLayerGeometry::lengthOfWireGroup ( int  wireGroup) const

Length of a wire group (center wire, across chamber face)

Definition at line 175 of file CSCLayerGeometry.cc.

References middleWireOfGroup(), theWireTopology, w(), and CSCWireTopology::wireValues().

                                                               {
  // Return length of 'wire' in the middle of the wire group
   float w = middleWireOfGroup( wireGroup );
   std::vector<float> store = theWireTopology->wireValues( w );
   return store[2];
}
LocalPoint CSCLayerGeometry::localCenterOfWireGroup ( int  wireGroup) const

Local coordinates of center of a wire group. Used to be centerOfWireGroup in ORCA but that version now returns GlobalPoint.

Definition at line 157 of file CSCLayerGeometry.cc.

References middleWireOfGroup(), theWireTopology, w(), wireAngle(), CSCWireTopology::wireValues(), detailsBasic3DVector::y, and yOfWireGroup().

Referenced by ValidateGeometry::validateCSCLayerGeometry().

                                                                         {

  // It can use CSCWireTopology::yOfWireGroup for y,
  // But x requires mixing with 'extent' of wire plane

  // If the wires are NOT tilted, default to simple calculation...
  if ( fabs(wireAngle() ) < 1.E-6 )  {
    float y = yOfWireGroup( wireGroup );
    return LocalPoint( 0., y );
  }
  else {
    // w is "wire" at the center of the wire group
    float w = middleWireOfGroup( wireGroup );
    std::vector<float> store = theWireTopology->wireValues( w );
    return LocalPoint( store[0], store[1] );
  }
}
LocalError CSCLayerGeometry::localError ( int  strip,
float  sigmaStrip,
float  sigmaWire 
) const

Transform strip and wire errors to local x, y frame. Need to supply (central) strip of the hit. The sigma's are in distance units.

Definition at line 231 of file CSCLayerGeometry.cc.

References funct::cos(), funct::sin(), stripAngle(), and wireAngle().

Referenced by CSCMake2DRecHit::hitFromStripAndWire().

                                                                                            {
  // Input sigmas are expected to be in _distance units_
  // - uncertainty in strip measurement (typically from Gatti fit, value is in local x units)
  // - uncertainty in wire measurement (along direction perpendicular to wires)

  float wangle   = this->wireAngle();
  float strangle = this->stripAngle( strip );

  float sinAngdif  = sin(strangle-wangle);
  float sinAngdif2 = sinAngdif * sinAngdif;
  
  float du = sigmaStrip/sin(strangle); // sigmaStrip is just x-component of strip error
  float dv = sigmaWire;

  // The notation is
  // wsins = wire resol  * sin(strip angle)
  // wcoss = wire resol  * cos(strip angle)
  // ssinw = strip resol * sin(wire angle)
  // scosw = strip resol * cos(wire angle)

  float wsins = dv * sin(strangle);
  float wcoss = dv * cos(strangle);
  float ssinw = du * sin(wangle);
  float scosw = du * cos(wangle);

  float dx2 = (scosw*scosw + wcoss*wcoss)/sinAngdif2;
  float dy2 = (ssinw*ssinw + wsins*wsins)/sinAngdif2;
  float dxy = (scosw*ssinw + wcoss*wsins)/sinAngdif2;
          
  return LocalError(dx2, dxy, dy2);
}
float CSCLayerGeometry::middleWireOfGroup ( int  wireGroup) const [inline]

Middle of wire-group. This is the central wire no. for a group with an odd no. of wires. This is a pseudo-wire no. for a group with an even no. of wires. Accordingly, it is non-integer.

Definition at line 201 of file CSCLayerGeometry.h.

References CSCWireTopology::middleWireOfGroup(), and theWireTopology.

Referenced by CSCMake2DRecHit::hitFromStripAndWire(), lengthOfWireGroup(), localCenterOfWireGroup(), and stripWireGroupIntersection().

int CSCLayerGeometry::nearestStrip ( const LocalPoint lp) const [inline]
int CSCLayerGeometry::nearestWire ( const LocalPoint lp) const [inline]

Wire nearest a given local point

Definition at line 101 of file CSCLayerGeometry.h.

References CSCWireTopology::nearestWire(), and theWireTopology.

Referenced by CSCWireHitSim::simulate().

                                               {
    return theWireTopology->nearestWire( lp ); }
int CSCLayerGeometry::numberOfStrips ( ) const [inline]
int CSCLayerGeometry::numberOfWireGroups ( ) const [inline]
int CSCLayerGeometry::numberOfWires ( ) const [inline]

How many wires in layer

Definition at line 66 of file CSCLayerGeometry.h.

References CSCWireTopology::numberOfWires(), and theWireTopology.

Referenced by operator<<().

int CSCLayerGeometry::numberOfWiresPerGroup ( int  wireGroup) const [inline]

How many wires in a wiregroup

Definition at line 78 of file CSCLayerGeometry.h.

References CSCWireTopology::numberOfWiresPerGroup(), and theWireTopology.

Referenced by CSCMake2DRecHit::hitFromStripAndWire().

CSCLayerGeometry & CSCLayerGeometry::operator= ( const CSCLayerGeometry melg)

Definition at line 78 of file CSCLayerGeometry.cc.

References apothem, CSCStripTopology::clone(), hBottomEdge, hTopEdge, theStripTopology, and theWireTopology.

{
  if ( &melg != this ) {
    delete theStripTopology;
    if ( melg.theStripTopology )
      theStripTopology=melg.theStripTopology->clone();
    else
      theStripTopology=0;

    delete theWireTopology;
    if ( melg.theWireTopology )
      theWireTopology=new CSCWireTopology(*(melg.theWireTopology));
    else
      theWireTopology=0;

    hBottomEdge     = melg.hBottomEdge;
    hTopEdge        = melg.hTopEdge;
    apothem         = melg.apothem;
  }
  return *this;
}
std::pair< LocalPoint, float > CSCLayerGeometry::possibleRecHitPosition ( float  s,
int  w1,
int  w2 
) const

Return estimate of the 2-dim point of intersection of a strip and a cluster of wires.

Input arguments: a (float) strip number, and the wires which delimit a cluster of wires. The wires are expected to be real wire numbers, and not wire-group numbers.

Returned: pair, with members:
first: LocalPoint which is midway along "the" strip between the wire limits,
or the chamber edges, as appropriate. <bf> second: length of the strip between the wires (or edges as appropriate).

Definition at line 182 of file CSCLayerGeometry.cc.

References intersectionOfStripAndWire(), TrapezoidalPlaneBounds::length(), mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

                                                                                                   {
        
  LocalPoint sw1 = intersectionOfStripAndWire( s, w1 );
  LocalPoint sw2 = intersectionOfStripAndWire( s, w2 );
                
  // Average the two points
  LocalPoint midpt( (sw1.x()+sw2.x())/2., (sw1.y()+sw2.y())/2 );
        
  // Length of strip crossing this group of wires
  float length = sqrt( (sw1.x()-sw2.x())*(sw1.x()-sw2.x()) + 
                     (sw1.y()-sw2.y())*(sw1.y()-sw2.y()) );
        
  return std::pair<LocalPoint,float>( midpt, length );
}
void CSCLayerGeometry::setTopology ( CSCStripTopology topology)

This class takes ownership of the pointer, and will destroy it

Definition at line 300 of file CSCLayerGeometry.cc.

References theStripTopology.

                                                                   {
   delete theStripTopology;
   theStripTopology = newTopology;
}
int CSCLayerGeometry::stagger ( ) const [inline]

Return +1 or -1 for a stripOffset of +0.25 or -0.25 respectively. Requested by trigger people.

Definition at line 129 of file CSCLayerGeometry.h.

References stripOffset().

Referenced by CSCCathodeLCTProcessor::run().

{ return static_cast<int>( 4.1*stripOffset() ); }
float CSCLayerGeometry::strip ( const LocalPoint lp) const [inline]

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.

Definition at line 193 of file CSCLayerGeometry.h.

References OffsetRadialStripTopology::strip(), and theStripTopology.

Referenced by CSCRecHit2DValidation::plotResolution().

{ return theStripTopology->strip(lp); }
float CSCLayerGeometry::stripAngle ( int  strip) const

The angle (in radians) of a strip wrt local x-axis.

Definition at line 140 of file CSCLayerGeometry.cc.

References M_PI_2, OffsetRadialStripTopology::stripAngle(), and theStripTopology.

Referenced by MuonCSCChamberResidual::addResidual(), CSCPairResidualsConstraint::calculatePhi(), CSCMake2DRecHit::hitFromStripAndWire(), localError(), CSCEfficiency::recHitSegment_Efficiencies(), CSCStripHitSim::simulate(), and stripWireIntersection().

{
  // Cleverly subtly change meaning of stripAngle once more.
  // In TrapezoidalStripTopology it is angle measured
  // counter-clockwise from y axis. 
  // In APTST and RST it is angle measured 
  // clockwise from y axis.
  // Output of this function is measured counter-clockwise 
  // from x-axis, so it is a conventional 2-dim azimuthal angle
  // in the (x,y) local coordinates

  // We want angle at centre of strip (strip N covers
  // *float* range N-1 to N-epsilon)

  return M_PI_2 - theStripTopology->stripAngle(strip-0.5);
}
float CSCLayerGeometry::stripOffset ( void  ) const [inline]

Offset of strips from symmetrical distribution about local y axis as a fraction of a strip (0 default, but usually +0.25 or -0.25)

Definition at line 123 of file CSCLayerGeometry.h.

References OffsetRadialStripTopology::stripOffset(), and theStripTopology.

Referenced by stagger().

float CSCLayerGeometry::stripPhiPitch ( ) const [inline]

The phi width of the strips (radians)

Definition at line 165 of file CSCLayerGeometry.h.

References CSCRadialStripTopology::phiPitch(), and theStripTopology.

                              { 
    return theStripTopology->phiPitch(); }
float CSCLayerGeometry::stripPitch ( ) const [inline]

The width of the strips (in middle)

Definition at line 171 of file CSCLayerGeometry.h.

Referenced by CSCMake2DRecHit::hitFromStripAndWire(), operator<<(), and CSCStripHitSim::simulate().

                           { 
    //    return theStripTopology->pitch(); }
    return stripPitch( LocalPoint(0.,0.) ); }
float CSCLayerGeometry::stripPitch ( const LocalPoint lp) const [inline]

The width of the strip at a given local point

Definition at line 178 of file CSCLayerGeometry.h.

References CSCRadialStripTopology::localPitch(), and theStripTopology.

                                                { 
    return theStripTopology->localPitch(lp); }
LocalPoint CSCLayerGeometry::stripWireGroupIntersection ( int  strip,
int  wireGroup 
) const

Local point at which strip and centre of wire group intersect

Definition at line 133 of file CSCLayerGeometry.cc.

References middleWireOfGroup(), and stripWireIntersection().

Referenced by CSCSectorReceiverLUT::getGlobalEtaValue().

{
  // middleWire is only an actual wire for a group with an odd no. of wires
  float middleWire = middleWireOfGroup( wireGroup );
  return stripWireIntersection(strip, middleWire);
}
LocalPoint CSCLayerGeometry::stripWireIntersection ( int  strip,
float  wire 
) const

Local point at which strip and wire intersect

Definition at line 111 of file CSCLayerGeometry.cc.

References stripAngle(), funct::tan(), wireAngle(), xOfStrip(), and yOfWire().

Referenced by CSCMake2DRecHit::hitFromStripAndWire(), and stripWireGroupIntersection().

{
  // This allows _float_ wire no. so that we can calculate the
  // intersection of a strip with the mid point of a wire group 
  // containing an even no. of wires (which is not an actual wire),
  // as well as for a group containing an odd no. of wires.

  // Equation of wire and strip as straight lines in local xy
  // y = mx + c where m = tan(angle w.r.t. x axis)
  // At the intersection x = -(cs-cw)/(ms-mw)
  // At y=0, 0 = ms * xOfStrip(strip) + cs => cs = -ms*xOfStrip
  // At x=0, yOfWire(wire) = 0 + cw => cw = yOfWire

  float ms = tan( stripAngle(strip) );
  float mw = tan( wireAngle() );
  float xs = xOfStrip(strip);
  float xi = ( ms * xs + yOfWire(wire) ) / ( ms - mw );
  float yi = ms * (xi - xs );

  return LocalPoint(xi, yi);
}
const CSCStripTopology* CSCLayerGeometry::topology ( ) const [inline]
float CSCLayerGeometry::wireAngle ( ) const [inline]

The angle (in radians) of (any) wire wrt local x-axis.

Definition at line 139 of file CSCLayerGeometry.h.

References theWireTopology, and CSCWireTopology::wireAngle().

Referenced by CSCDriftSim::getWireHit(), localCenterOfWireGroup(), localError(), operator<<(), CSCStripHitSim::simulate(), and stripWireIntersection().

                          {
    return theWireTopology->wireAngle(); }
int CSCLayerGeometry::wireGroup ( int  wire) const [inline]

Wire group containing a given wire

Definition at line 107 of file CSCLayerGeometry.h.

References theWireTopology, and CSCWireTopology::wireGroup().

Referenced by CSCWireElectronicsSim::readoutElement(), and yResolution().

                                {
    return theWireTopology->wireGroup( wire );
  }
float CSCLayerGeometry::wirePitch ( ) const [inline]

The distance (in cm) between anode wires

Definition at line 145 of file CSCLayerGeometry.h.

References theWireTopology, and CSCWireTopology::wirePitch().

Referenced by operator<<(), and CSCChamberSpecs::wireSpacing().

                          {
    return theWireTopology->wirePitch(); }
const CSCWireTopology* CSCLayerGeometry::wireTopology ( ) const [inline]
float CSCLayerGeometry::xOfStrip ( int  strip,
float  y = 0. 
) const [inline]
std::pair<float, float> CSCLayerGeometry::yLimitsOfStripPlane ( ) const [inline]

Local y limits of the strip plane

Definition at line 231 of file CSCLayerGeometry.h.

References theStripTopology, and CSCStripTopology::yLimitsOfStripPlane().

Referenced by inside().

float CSCLayerGeometry::yOfWire ( float  wire,
float  x = 0. 
) const [inline]

Local y of a given wire 'number' (float) at given x

Definition at line 207 of file CSCLayerGeometry.h.

References theWireTopology, x, and CSCWireTopology::yOfWire().

Referenced by CSCDriftSim::getWireHit(), CSCMake2DRecHit::hitFromStripAndWire(), CSCStripHitSim::simulate(), and stripWireIntersection().

                                              {
    return theWireTopology->yOfWire( wire, x ); }
float CSCLayerGeometry::yOfWireGroup ( int  wireGroup,
float  x = 0. 
) const [inline]

Local y of a given wire group at given x

Definition at line 213 of file CSCLayerGeometry.h.

References theWireTopology, x, and CSCWireTopology::yOfWireGroup().

Referenced by CSCEfficiency::fillWG_info(), localCenterOfWireGroup(), and CSCWireDigiValidation::plotResolution().

float CSCLayerGeometry::yResolution ( int  wireGroup = 1) const [inline]

The measurement resolution from wire groups (in cm.) This approximates the measurement resolution in the local y direction but may be too small by a factor of up to 1.26 due to stripAngle contributions which are neglected here. The last wiregroup may have more wires than others. The other wiregroups, including the first, are the same. One day the wiregroups will be matched to the hardware by using the DDD.

Definition at line 159 of file CSCLayerGeometry.h.

References theWireTopology, wireGroup(), and CSCWireTopology::yResolution().

Referenced by CSCMake2DRecHit::hitFromStripAndWire().


Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  stream,
const CSCLayerGeometry lg 
) [friend]

Output operator for members of class.

Definition at line 305 of file CSCLayerGeometry.cc.

                                                                          {
  stream << "LayerGeometry " << std::endl
         << "------------- " << std::endl
         << "numberOfStrips               " << lg.numberOfStrips() << std::endl
         << "numberOfWires                " << lg.numberOfWires() << std::endl
         << "numberOfWireGroups           " << lg.numberOfWireGroups() << std::endl
         << "wireAngle  (rad)             " << lg.wireAngle() << std::endl
    //         << "wireAngle  (deg)      " << lg.theWireAngle << std::endl
    //         << "sin(wireAngle)        " << lg.theWireSin << std::endl
    //         << "cos(wireAngle)        " << lg.theWireCos << std::endl
         << "wirePitch                    " << lg.wirePitch() << std::endl
         << "stripPitch                   " << lg.stripPitch() << std::endl
    //         << "numberOfWiresPerGroup " << lg.theNumberOfWiresPerGroup << std::endl
    //         << "numberOfWiresInLastGroup " << lg.theNumberOfWiresInLastGroup << std::endl
    //         << "wireOffset            " << lg.theWireOffset << std::endl
    //         << "whereStripsMeet       " << lg.whereStripsMeet << std::endl;
         << "hBottomEdge                  " << lg.hBottomEdge << std::endl
         << "hTopEdge                     " << lg.hTopEdge << std::endl
         << "apothem                      " << lg.apothem << std::endl
         << "length (should be 2xapothem) " << lg.length() << std::endl
         << "thickness                    " << lg.thickness() << std::endl;
    return stream;
}

Member Data Documentation

float CSCLayerGeometry::apothem [private]

Definition at line 334 of file CSCLayerGeometry.h.

Referenced by CSCLayerGeometry(), operator<<(), and operator=().

Definition at line 337 of file CSCLayerGeometry.h.

Definition at line 332 of file CSCLayerGeometry.h.

Referenced by CSCLayerGeometry(), operator<<(), and operator=().

float CSCLayerGeometry::hTopEdge [private]

Definition at line 333 of file CSCLayerGeometry.h.

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

const std::string CSCLayerGeometry::myName [private]

Definition at line 336 of file CSCLayerGeometry.h.

Referenced by CSCLayerGeometry(), and ~CSCLayerGeometry().