CMS 3D CMS Logo

Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes

TIBRing Class Reference

#include <TIBRing.h>

Inheritance diagram for TIBRing:
GeometricSearchDetWithGroups GeometricSearchDet

List of all members.

Classes

struct  SubRingCrossings

Public Member Functions

virtual const std::vector
< const GeomDet * > & 
basicComponents () const
virtual std::pair< bool,
TrajectoryStateOnSurface
compatible (const TrajectoryStateOnSurface &ts, const Propagator &, const MeasurementEstimator &) const
virtual const std::vector
< const GeometricSearchDet * > & 
components () const
 Returns basic components, if any.
virtual void groupedCompatibleDetsV (const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est, std::vector< DetGroup > &result) const
virtual const BoundCylinderspecificSurface () const
 Return the ring surface as a.
virtual const BoundSurfacesurface () const
 The surface of the GeometricSearchDet.
 TIBRing (std::vector< const GeomDet * > &theGeomDets)
 ~TIBRing ()

Private Types

typedef PeriodicBinFinderInPhi
< double > 
BinFinderType

Private Member Functions

void checkPeriodicity (std::vector< const GeomDet * >::const_iterator first, std::vector< const GeomDet * >::const_iterator last)
void checkRadius (std::vector< const GeomDet * >::const_iterator first, std::vector< const GeomDet * >::const_iterator last)
SubRingCrossings computeCrossings (const TrajectoryStateOnSurface &startingState, PropagationDirection propDir) const
void computeHelicity ()
float computeWindowSize (const GeomDet *det, const TrajectoryStateOnSurface &tsos, const MeasurementEstimator &est) const
void searchNeighbors (const TrajectoryStateOnSurface &tsos, const Propagator &prop, const MeasurementEstimator &est, const SubRingCrossings &crossings, float window, std::vector< DetGroup > &result) const

Private Attributes

BinFinderType theBinFinder
ReferenceCountingPointer
< BoundCylinder
theCylinder
std::vector< const GeomDet * > theDets
int theHelicity

Detailed Description

A concrete implementation for TIB rings

Definition at line 13 of file TIBRing.h.


Member Typedef Documentation

Definition at line 86 of file TIBRing.h.


Constructor & Destructor Documentation

TIBRing::TIBRing ( std::vector< const GeomDet * > &  theGeomDets)

Definition at line 21 of file TIBRing.cc.

References computeHelicity(), i, LogDebug, theBinFinder, theCylinder, and theDets.

                                                   :
  theDets(theGeomDets.begin(),theGeomDets.end())
{
  //checkRadius( first, last);
  //sort( theDets.begin(), theDets.end(), DetLessPhi());
  //checkPeriodicity( theDets.begin(), theDets.end());

  theBinFinder = BinFinderType( theDets.front()->surface().position().phi(),
                                theDets.size());

  theCylinder = CylinderBuilderFromDet()( theDets.begin(), theDets.end());

  computeHelicity();

  
  LogDebug("TkDetLayers") << "==== DEBUG TIBRing =====" ; 
  LogDebug("TkDetLayers") << "radius, thickness, lenght: " 
                          << theCylinder->radius() << " , "
                          << theCylinder->bounds().thickness() << " , "
                          << theCylinder->bounds().length() ;

  for (vector<const GeomDet*>::const_iterator i=theDets.begin();
       i != theDets.end(); i++){
    LogDebug("TkDetLayers") << "Ring's Det pos z,perp,eta,phi: " 
         << (**i).position().z() << " , " 
         << (**i).position().perp() << " , " 
         << (**i).position().eta() << " , " 
         << (**i).position().phi() ;
  }
  LogDebug("TkDetLayers") << "==== end DEBUG TIBRing =====" ; 
 
  
}
TIBRing::~TIBRing ( )

Definition at line 118 of file TIBRing.cc.

                 {

} 

Member Function Documentation

virtual const std::vector<const GeomDet*>& TIBRing::basicComponents ( ) const [inline, virtual]

Implements GeometricSearchDet.

Definition at line 21 of file TIBRing.h.

References theDets.

{return theDets;}
void TIBRing::checkPeriodicity ( std::vector< const GeomDet * >::const_iterator  first,
std::vector< const GeomDet * >::const_iterator  last 
) [private]

Definition at line 78 of file TIBRing.cc.

References first, i, j, pi, stat_mean(), launcher::step, and theDets.

{
  vector<double> adj_diff(last-first-1);
  for (int i=0; i < static_cast<int>(adj_diff.size()); i++) {
    vector<const GeomDet*>::const_iterator curent = first + i;
    adj_diff[i] = (**(curent+1)).surface().position().phi() - 
                  (**curent).surface().position().phi();
  }
  double step = stat_mean( adj_diff);
  double phi_step = 2.*Geom::pi()/(last-first);  

  if ( fabs(step-phi_step)/phi_step > 0.01) {
    int ndets = last-first;
    edm::LogError("TkDetLayers") << "TIBRing Warning: not periodic. ndets=" << ndets ;
    for (int j=0; j<ndets; j++) {
      edm::LogError("TkDetLayers") << "Dets(r,phi): (" << theDets[j]->surface().position().perp() 
                                 << "," << theDets[j]->surface().position().phi() << ") " ;
    }
    throw DetLayerException( "Error: TIBRing is not periodic");
  }
}
void TIBRing::checkRadius ( std::vector< const GeomDet * >::const_iterator  first,
std::vector< const GeomDet * >::const_iterator  last 
) [private]

Definition at line 62 of file TIBRing.cc.

References i, prof2calltree::last, and csvReporter::r.

{
  // check radius range
  float rMin = 10000.;
  float rMax = 0.;
  for (vector<const GeomDet*>::const_iterator i=first; i!=last; i++) {
    float r = (**i).surface().position().perp();
    if (r < rMin) rMin = r;
    if (r > rMax) rMax = r;
  }
  if (rMax - rMin > 0.1) throw DetLayerException(
    "TIBRing construction failed: detectors not at constant radius");
}
pair< bool, TrajectoryStateOnSurface > TIBRing::compatible ( const TrajectoryStateOnSurface ts,
const Propagator ,
const MeasurementEstimator  
) const [virtual]

tests the geometrical compatibility of the Det with the predicted state. The FreeTrajectoryState argument is propagated to the Det surface using the Propagator argument. The resulting TrajectoryStateOnSurface is tested for compatibility with the surface bounds. If compatible, a std::pair< true, propagatedState> is returned. If the propagation fails, or if the state is not compatible, a std::pair< false, propagatedState> is returned.

Implements GeometricSearchDet.

Definition at line 124 of file TIBRing.cc.

                                                    {
  edm::LogError("TkDetLayers") << "temporary dummy implementation of TIBRing::compatible()!!" ;
  return pair<bool,TrajectoryStateOnSurface>();
}
const vector< const GeometricSearchDet * > & TIBRing::components ( ) const [virtual]

Returns basic components, if any.

Returns direct components, if any

Implements GeometricSearchDet.

Definition at line 57 of file TIBRing.cc.

{
  throw DetLayerException("TIBRing doesn't have GeometricSearchDet components");
}
TIBRing::SubRingCrossings TIBRing::computeCrossings ( const TrajectoryStateOnSurface startingState,
PropagationDirection  propDir 
) const [private]

Definition at line 227 of file TIBRing.cc.

References PeriodicBinFinderInPhi< T >::binIndex(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), PV3DBase< T, PVType, FrameType >::phi(), GloballyPositioned< T >::position(), rho, specificSurface(), surface(), theBinFinder, theDets, TrajectoryStateOnSurface::transverseCurvature(), and PV3DBase< T, PVType, FrameType >::x().

Referenced by GeometricSearchDet::groupedCompatibleDetsV().

{
  typedef HelixBarrelPlaneCrossing2OrderLocal    Crossing;
  typedef MeasurementEstimator::Local2DVector Local2DVector;

  GlobalPoint startPos( startingState.globalPosition());
  GlobalVector startDir( startingState.globalMomentum());
  double rho( startingState.transverseCurvature());

  HelixBarrelCylinderCrossing cylCrossing( startPos, startDir, rho,
                                           propDir,specificSurface());

  if (!cylCrossing.hasSolution()) return SubRingCrossings();

  GlobalPoint  cylPoint( cylCrossing.position());
  GlobalVector cylDir( cylCrossing.direction());
  int closestIndex = theBinFinder.binIndex(cylPoint.phi());

  const BoundPlane& closestPlane( theDets[closestIndex]->surface());

  LocalPoint closestPos = Crossing( cylPoint, cylDir, rho, closestPlane).position();
  float closestDist = closestPos.x(); // use fact that local X perp to global Z 

  //int next = cylPoint.phi() - closestPlane.position().phi() > 0 ? closest+1 : closest-1;
  int nextIndex = PhiLess()( closestPlane.position().phi(), cylPoint.phi()) ? 
    closestIndex+1 : closestIndex-1;

  const BoundPlane& nextPlane( theDets[ theBinFinder.binIndex(nextIndex)]->surface());
  LocalPoint nextPos = Crossing( cylPoint, cylDir, rho, nextPlane).position();
  float nextDist = nextPos.x();

  if (fabs(closestDist) < fabs(nextDist)) {
    return SubRingCrossings( closestIndex, nextIndex, nextDist);
  }
  else {
    return SubRingCrossings( nextIndex, closestIndex, closestDist);
  }
}
void TIBRing::computeHelicity ( ) [private]

Definition at line 101 of file TIBRing.cc.

References Vector3DBase< T, FrameTag >::dot(), PV3DBase< T, PVType, FrameType >::phi(), GloballyPositioned< T >::position(), GeomDet::surface(), theDets, theHelicity, and Surface::toGlobal().

Referenced by TIBRing().

                              {

  const GeomDet& det = *theDets.front();
  GlobalVector radial = det.surface().position() - GlobalPoint(0,0,0);
  GlobalVector normal = det.surface().toGlobal( LocalVector(0,0,1));
  if(normal.dot(radial)<=0)normal*=-1;
//   edm::LogInfo(TkDetLayers) << "BarrelDetRing::computeHelicity: phi(normal) " << normal.phi()
//        << " phi(radial) " << radial.phi() ;
  if (PhiLess()( normal.phi(), radial.phi())) {
    theHelicity = 1;  // smaller phi angles mean "inner" group
  }
  else {
    theHelicity = 0;  // smaller phi angles mean "outer" group
  }
}
float TIBRing::computeWindowSize ( const GeomDet det,
const TrajectoryStateOnSurface tsos,
const MeasurementEstimator est 
) const [private]
virtual void TIBRing::groupedCompatibleDetsV ( const TrajectoryStateOnSurface startingState,
const Propagator prop,
const MeasurementEstimator est,
std::vector< DetGroup > &  result 
) const [virtual]

Reimplemented from GeometricSearchDet.

void TIBRing::searchNeighbors ( const TrajectoryStateOnSurface tsos,
const Propagator prop,
const MeasurementEstimator est,
const SubRingCrossings crossings,
float  window,
std::vector< DetGroup > &  result 
) const [private]
virtual const BoundCylinder& TIBRing::specificSurface ( ) const [inline, virtual]

Return the ring surface as a.

Definition at line 41 of file TIBRing.h.

References theCylinder.

Referenced by computeCrossings().

{return *theCylinder;}
virtual const BoundSurface& TIBRing::surface ( ) const [inline, virtual]

The surface of the GeometricSearchDet.

Implements GeometricSearchDet.

Definition at line 19 of file TIBRing.h.

References theCylinder.

Referenced by computeCrossings().

{return *theCylinder;}  

Member Data Documentation

Definition at line 87 of file TIBRing.h.

Referenced by computeCrossings(), GeometricSearchDet::groupedCompatibleDetsV(), and TIBRing().

Definition at line 90 of file TIBRing.h.

Referenced by specificSurface(), surface(), and TIBRing().

std::vector<const GeomDet*> TIBRing::theDets [private]
int TIBRing::theHelicity [private]

Definition at line 91 of file TIBRing.h.

Referenced by computeHelicity(), and GeometricSearchDet::groupedCompatibleDetsV().