CMS 3D CMS Logo

RoadSearchCircleSeed Class Reference

#include <RecoTracker/RoadSearchSeedFinder/interface/RoadSearchCircleSeed.h>

List of all members.

Public Types

typedef std::pair< double, double > line
enum  type { circle, straightLine }

Public Member Functions

std::vector< const
TrackingRecHit * >
::const_iterator 
begin_hits () const
std::vector< GlobalPoint >
::const_iterator 
begin_points () const
double calculateEta (double theta) const
double calculateImpactParameter (GlobalPoint &center, double radius)
bool calculateInBarrel ()
GlobalPoint Center () const
bool Compare (const RoadSearchCircleSeed *circle, double centerCut, double radiusCut, unsigned int differentHitsCut) const
bool CompareCenter (const RoadSearchCircleSeed *circle, double centerCut) const
bool CompareDifferentHits (const RoadSearchCircleSeed *circle, unsigned int differentHitsCut) const
bool CompareRadius (const RoadSearchCircleSeed *circle, double radiusCut) const
double determinant (double array[][3], unsigned int bins)
std::vector< const
TrackingRecHit * >
::const_iterator 
end_hits () const
std::vector< GlobalPoint >
::const_iterator 
end_points () const
double Eta () const
const Roads::RoadSeedgetSeed ()
const Roads::RoadSetgetSet ()
std::vector< const
TrackingRecHit * > 
Hits () const
double ImpactParameter () const
bool InBarrel () const
double Phi0 () const
std::vector< GlobalPointPoints () const
std::string print () const
double Radius () const
 RoadSearchCircleSeed (const TrackingRecHit *hit1, const TrackingRecHit *hit2, GlobalPoint &point1, GlobalPoint &point2)
 RoadSearchCircleSeed (const TrackingRecHit *hit1, const TrackingRecHit *hit2, const TrackingRecHit *hit3, GlobalPoint &point1, GlobalPoint &point2, GlobalPoint &point3)
void setSeed (const Roads::RoadSeed *input)
void setSet (const Roads::RoadSet *input)
double Theta () const
double Type () const
 ~RoadSearchCircleSeed ()

Private Attributes

GlobalPoint center_
std::vector< const
TrackingRecHit * > 
hits_
double impactParameter_
bool inBarrel_
std::vector< GlobalPointpoints_
double radius_
const Roads::RoadSeedseed_
const Roads::RoadSetset_
type type_


Detailed Description

Definition at line 29 of file RoadSearchCircleSeed.h.


Member Typedef Documentation

typedef std::pair<double,double> RoadSearchCircleSeed::line

Definition at line 34 of file RoadSearchCircleSeed.h.


Member Enumeration Documentation

enum RoadSearchCircleSeed::type

Enumerator:
circle 
straightLine 

Definition at line 36 of file RoadSearchCircleSeed.h.

00036             {
00037     circle,
00038     straightLine
00039   };


Constructor & Destructor Documentation

RoadSearchCircleSeed::RoadSearchCircleSeed ( const TrackingRecHit hit1,
const TrackingRecHit hit2,
const TrackingRecHit hit3,
GlobalPoint point1,
GlobalPoint point2,
GlobalPoint point3 
)

Definition at line 25 of file RoadSearchCircleSeed.cc.

References calculateImpactParameter(), calculateInBarrel(), center_, circle, hits_, impactParameter_, inBarrel_, FastCircle::isValid(), points_, radius_, FastCircle::rho(), straightLine, type_, FastCircle::x0(), and FastCircle::y0().

00030                                                                 { 
00031 
00032   hits_.reserve(3);
00033   hits_.push_back(hit1);
00034   hits_.push_back(hit2);
00035   hits_.push_back(hit3);
00036 
00037   points_.reserve(3);
00038   points_.push_back(point1);
00039   points_.push_back(point2);
00040   points_.push_back(point3);
00041 
00042   FastCircle kreis(point1,
00043                    point2,
00044                    point3);
00045 
00046   if ( !kreis.isValid() ) {
00047     // line
00048     type_ = straightLine;
00049     inBarrel_ = true; // Not used for lines
00050     center_ = GlobalPoint(0,0,0);
00051     radius_ = 0;
00052     impactParameter_ = 0;
00053   } else {
00054     type_ = circle;
00055     inBarrel_        = calculateInBarrel();
00056     radius_          = kreis.rho();
00057     center_          = GlobalPoint(kreis.x0(),kreis.y0(),0);
00058     impactParameter_ = calculateImpactParameter(center_,radius_);
00059   }
00060 
00061 }

RoadSearchCircleSeed::RoadSearchCircleSeed ( const TrackingRecHit hit1,
const TrackingRecHit hit2,
GlobalPoint point1,
GlobalPoint point2 
)

Definition at line 63 of file RoadSearchCircleSeed.cc.

References center_, hits_, impactParameter_, inBarrel_, points_, radius_, straightLine, and type_.

00066                                                                 { 
00067   //
00068   // straight line constructor
00069   //
00070 
00071   hits_.reserve(2);
00072   hits_.push_back(hit1);
00073   hits_.push_back(hit2);
00074 
00075   points_.reserve(2);
00076   points_.push_back(point1);
00077   points_.push_back(point2);
00078 
00079   type_ = straightLine;
00080   inBarrel_ = true; // Not used for lines
00081   center_ = GlobalPoint(0,0,0);
00082   radius_ = 0;
00083   impactParameter_ = 0;
00084 
00085 }

RoadSearchCircleSeed::~RoadSearchCircleSeed (  ) 

Definition at line 87 of file RoadSearchCircleSeed.cc.

00087                                             {
00088 }


Member Function Documentation

std::vector<const TrackingRecHit*>::const_iterator RoadSearchCircleSeed::begin_hits (  )  const [inline]

Definition at line 59 of file RoadSearchCircleSeed.h.

References hits_.

Referenced by CompareDifferentHits().

00059 { return hits_.begin(); }

std::vector<GlobalPoint>::const_iterator RoadSearchCircleSeed::begin_points (  )  const [inline]

Definition at line 55 of file RoadSearchCircleSeed.h.

References points_.

00055 { return points_.begin(); }

double RoadSearchCircleSeed::calculateEta ( double  theta  )  const

Definition at line 118 of file RoadSearchCircleSeed.cc.

References funct::log(), and funct::tan().

Referenced by Eta().

00118                                                             {
00119   //
00120   // calculate eta from theta
00121   //
00122 
00123   return -1.*std::log(std::tan(theta/2.));
00124 }

double RoadSearchCircleSeed::calculateImpactParameter ( GlobalPoint center,
double  radius 
)

Definition at line 106 of file RoadSearchCircleSeed.cc.

References funct::abs(), d, funct::sqrt(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by RoadSearchCircleSeed().

00107                                                                      {
00108   //
00109   // calculate impact parameter to (0,0,0) from center and radius of circle
00110   //
00111 
00112   double d = std::sqrt( center.x() * center.x() +
00113                         center.y() * center.y() );
00114 
00115   return std::abs(d-radius);
00116 }

bool RoadSearchCircleSeed::calculateInBarrel (  ) 

Definition at line 90 of file RoadSearchCircleSeed.cc.

References hits_, and StripSubdetector::TEC.

Referenced by RoadSearchCircleSeed().

00090                                              {
00091   //
00092   // returns true if all hits are in the barrel,
00093   // otherwise returns false
00094   //
00095 
00096   for (std::vector<const TrackingRecHit*>::const_iterator hit = hits_.begin();
00097        hit != hits_.end(); ++hit) {
00098     if ((*hit)->geographicalId().subdetId() == StripSubdetector::TEC) {
00099       return false;
00100     }
00101   }
00102   
00103   return true;
00104 }

GlobalPoint RoadSearchCircleSeed::Center (  )  const [inline]

Definition at line 62 of file RoadSearchCircleSeed.h.

References center_.

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

00062 { return center_;}

bool RoadSearchCircleSeed::Compare ( const RoadSearchCircleSeed circle,
double  centerCut,
double  radiusCut,
unsigned int  differentHitsCut 
) const

Definition at line 316 of file RoadSearchCircleSeed.cc.

References CompareCenter(), CompareDifferentHits(), CompareRadius(), and HLT_VtxMuL3::result.

Referenced by RoadSearchSeedFinderAlgorithm::calculateCircleSeedsFromHits().

00319                                                                         {
00320   //
00321   // compare this circle with the input circle
00322   // compare: percentage of center difference of center average
00323   // compare: percentage of radius difference of radius average
00324   // compare: number of hits which don't overlap between the two circles
00325   //
00326 
00327   // return value
00328   bool result = false;
00329 
00330   result = CompareRadius(circle,radiusCut);
00331   if ( result ) {
00332     result = CompareCenter(circle,centerCut);
00333     if ( result ) {
00334       result = CompareDifferentHits(circle,differentHitsCut);
00335     }
00336   }
00337 
00338   return result;
00339 
00340 }

bool RoadSearchCircleSeed::CompareCenter ( const RoadSearchCircleSeed circle,
double  centerCut 
) const

Definition at line 342 of file RoadSearchCircleSeed.cc.

References Center(), center_, HLT_VtxMuL3::result, funct::sqrt(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by Compare().

00343                                                                  {
00344   //
00345   // compare this circle with the input circle
00346   // compare: percentage of center difference of center average
00347   //
00348 
00349   // return value
00350   bool result = false;
00351 
00352   double averageCenter = std::sqrt(((center_.x()+circle->Center().x())/2) *
00353                                    ((center_.x()+circle->Center().x())/2) +
00354                                    ((center_.y()+circle->Center().y())/2) *
00355                                    ((center_.y()+circle->Center().y())/2));
00356   double differenceCenter = std::sqrt((center_.x()-circle->Center().x()) *
00357                                       (center_.x()-circle->Center().x()) +
00358                                       (center_.y()-circle->Center().y()) *
00359                                       (center_.y()-circle->Center().y()));
00360 
00361   if ( differenceCenter/averageCenter <= centerCut ) {
00362     result = true;
00363   }
00364 
00365 //   edm::LogVerbatim("OLI") << "center difference: " << differenceCenter
00366 //                        << "center average: " << averageCenter
00367 //                        << "center percentage: " << differenceCenter/averageCenter
00368 //                        << " cut: " << centerCut
00369 //                        << " result: " << result;
00370 
00371   return result;
00372 
00373 }

bool RoadSearchCircleSeed::CompareDifferentHits ( const RoadSearchCircleSeed circle,
unsigned int  differentHitsCut 
) const

Definition at line 402 of file RoadSearchCircleSeed.cc.

References begin_hits(), counter(), end_hits(), hits_, and HLT_VtxMuL3::result.

Referenced by Compare().

00403                                                                                      {
00404   //
00405   // compare this circle with the input circle
00406   // compare: number of hits which don't overlap between the two circles
00407   //
00408 
00409   // return value
00410   bool result = false;
00411 
00412   // assume circles always have 3 hits
00413   unsigned int counter = 0;
00414   for ( std::vector<const TrackingRecHit*>::const_iterator hit1 = hits_.begin(),
00415           hit1End = hits_.end();
00416         hit1 != hit1End;
00417         ++hit1 ) {
00418     bool included = false;
00419     for ( std::vector<const TrackingRecHit*>::const_iterator hit2 = circle->begin_hits(),
00420             hit2End = circle->end_hits();
00421           hit2 != hit2End;
00422           ++hit2 ) {
00423       if ( *hit1 == *hit2 ) {
00424         included = true;
00425       }
00426     }
00427     if ( !included ) {
00428       ++counter;
00429     }
00430   }
00431 
00432   if ( counter <= differentHitsCut ) {
00433     result = true;
00434   }
00435 
00436 //   edm::LogVerbatim("OLI") << "hits: " << counter 
00437 //                        << " cut: " << differentHitsCut 
00438 //                        << " result: " << result;
00439 
00440   return result;
00441 
00442 }

bool RoadSearchCircleSeed::CompareRadius ( const RoadSearchCircleSeed circle,
double  radiusCut 
) const

Definition at line 375 of file RoadSearchCircleSeed.cc.

References funct::abs(), Radius(), radius_, and HLT_VtxMuL3::result.

Referenced by Compare().

00376                                                                  {
00377   //
00378   // compare: percentage of center difference of center average
00379   // compare: percentage of radius difference of radius average
00380   //
00381 
00382   // return value
00383   bool result = false;
00384 
00385   double averageRadius = (radius_ + circle->Radius() ) /2;
00386   double differenceRadius = std::abs(radius_ - circle->Radius());
00387   
00388   if ( differenceRadius/averageRadius <= radiusCut ) {
00389     result = true;
00390   }
00391 
00392 //   edm::LogVerbatim("OLI") << "radius difference: " << differenceRadius
00393 //                        << " radius average: " << averageRadius
00394 //                        << " radius percentage: " << differenceRadius/averageRadius
00395 //                        << " cut: " << radiusCut
00396 //                        << " result: " << result;
00397 
00398   return result;
00399 
00400 }

double RoadSearchCircleSeed::determinant ( double  array[][3],
unsigned int  bins 
)

std::vector<const TrackingRecHit*>::const_iterator RoadSearchCircleSeed::end_hits (  )  const [inline]

Definition at line 60 of file RoadSearchCircleSeed.h.

References hits_.

Referenced by CompareDifferentHits().

00060 { return hits_.end();   }

std::vector<GlobalPoint>::const_iterator RoadSearchCircleSeed::end_points (  )  const [inline]

Definition at line 56 of file RoadSearchCircleSeed.h.

References points_.

00056 { return points_.end();   }

double RoadSearchCircleSeed::Eta (  )  const [inline]

Definition at line 65 of file RoadSearchCircleSeed.h.

References calculateEta(), and Theta().

00065 { return calculateEta(Theta());}

const Roads::RoadSeed* RoadSearchCircleSeed::getSeed (  )  [inline]

Definition at line 70 of file RoadSearchCircleSeed.h.

References seed_.

00070 { return seed_;  }

const Roads::RoadSet* RoadSearchCircleSeed::getSet (  )  [inline]

Definition at line 73 of file RoadSearchCircleSeed.h.

References set_.

00073 { return set_;  }

std::vector<const TrackingRecHit*> RoadSearchCircleSeed::Hits (  )  const [inline]

Definition at line 58 of file RoadSearchCircleSeed.h.

References hits_.

00058 { return hits_; }

double RoadSearchCircleSeed::ImpactParameter (  )  const [inline]

Definition at line 64 of file RoadSearchCircleSeed.h.

References impactParameter_.

Referenced by RoadSearchSeedFinderAlgorithm::calculateCircleSeedsFromHits(), operator<<(), and Phi0().

00064 { return impactParameter_;}

bool RoadSearchCircleSeed::InBarrel (  )  const [inline]

Definition at line 67 of file RoadSearchCircleSeed.h.

References inBarrel_.

Referenced by RoadSearchSeedFinderAlgorithm::calculateCircleSeedsFromHits(), and operator<<().

00067 { return inBarrel_; }

double RoadSearchCircleSeed::Phi0 (  )  const

Definition at line 157 of file RoadSearchCircleSeed.cc.

References PV3DBase< T, PVType, FrameType >::barePhi(), center_, i, ImpactParameter(), PV3DBase< T, PVType, FrameType >::perp(), phi, Geom::pi(), points_, Radius(), straightLine, and type_.

00157                                         {
00158   //
00159   // calculate the angle in the x-y plane
00160   // of the momentum vector at the point of
00161   // closest approach to (0,0,0)
00162   //
00163   // Note:  A faster implementation would
00164   // calculate in the constructor, save, and
00165   // return the member here.  This implementation
00166   // minimizes the memory footprint.
00167   //
00168 
00169   // Calculate phi as the average phi of all
00170   // lines formed by combinations of hits if
00171   // this is a straight line
00172   if (type_ == straightLine) {
00173     std::vector<LineXY> lines;
00174     for (std::vector<GlobalPoint>::const_iterator point1 = points_.begin();
00175          point1 != points_.end(); ++point1) {
00176       for (std::vector<GlobalPoint>::const_iterator point2 = point1+1;
00177            point2 != points_.end(); ++point2) {
00178         lines.push_back(LineXY(*point1,*point2));
00179       }
00180     }
00181     double netPhi = 0.;
00182     for (std::vector<LineXY>::const_iterator line = lines.begin();
00183          line != lines.end(); ++line) {
00184       netPhi += line->Phi();
00185     }
00186     return netPhi/(double)lines.size();
00187   } // END calculation for linear seeds
00188 
00189   // This calculation is not valid for seeds which do not exit
00190   // the tracking detector (lines always exit)
00191   else if (2.*Radius()+ImpactParameter()<110) {
00192     return 100000.;
00193   }
00194 
00195   // circular seeds
00196   else {
00197     double phi = 100000.;
00198     double centerPhi = center_.barePhi();
00199 
00200     // Find the first hit in time, which determines the direction of
00201     // the momentum vector (tangent to the circle at the point of
00202     // closest approach, clockwise or counter-clockwise).
00203     // The first hit in time is always the hit with the smallest
00204     // value r as long as the track exits the tracking detector.
00205     GlobalPoint firstPoint = points_[0];
00206     for (unsigned int i=1; i<points_.size(); ++i) {
00207       if (firstPoint.perp() > points_[i].perp()) {
00208         firstPoint = points_[i];
00209       }
00210     }
00211     
00212     // Get the next hit, firstPoint is at the point of
00213     // closest approach and cannot be used to
00214     // determine the direction of the initial
00215     // momentum vector
00216     if (firstPoint.barePhi() == centerPhi) {
00217       GlobalPoint nextHit = points_[0];
00218       for (unsigned int i=1; i<points_.size(); ++i) {
00219         if (nextHit.perp()  == firstPoint.perp() || 
00220             (firstPoint.perp()!= points_[i].perp() &&
00221              nextHit.perp() >  points_[i].perp())) {
00222           nextHit = points_[i];
00223         }
00224       }
00225       firstPoint = nextHit;
00226     }
00227   
00228     // Find the direction of the momentum vector
00229     if (firstPoint.barePhi() > centerPhi) {
00230       // The momentum vector is tangent to
00231       // the track
00232       phi = centerPhi + Geom::pi()/2.;
00233       if (phi>Geom::pi()) {
00234         phi -= 2.*Geom::pi();
00235       }
00236     }
00237     // Other direction!
00238     else if (firstPoint.barePhi() < centerPhi) {
00239       // The momentum vector is tangent to
00240       // the track
00241       phi = centerPhi - Geom::pi()/2.;
00242       if (phi<-1.*Geom::pi()) {
00243         phi += 2.*Geom::pi();
00244       }
00245     }  
00246     return phi;
00247   } // END calculation for circular seeds
00248 }

std::vector<GlobalPoint> RoadSearchCircleSeed::Points (  )  const [inline]

Definition at line 54 of file RoadSearchCircleSeed.h.

References points_.

Referenced by operator<<().

00054 { return points_; }

std::string RoadSearchCircleSeed::print ( void   )  const

Definition at line 250 of file RoadSearchCircleSeed.cc.

References center_, counter(), impactParameter_, inBarrel_, points_, radius_, straightLine, type_, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

00250                                             {
00251   //
00252   // print function
00253   //
00254 
00255   std::ostringstream ost;
00256 
00257   if ( type_ == RoadSearchCircleSeed::straightLine ) {
00258     ost << "Straight Line: number of points: " << points_.size() << "\n";
00259     unsigned int counter = 0;
00260     for ( std::vector<GlobalPoint>::const_iterator point = points_.begin();
00261           point != points_.end();
00262           ++point ) {
00263       ++counter;
00264       ost << "    Point " << counter << ": " << point->x() << "," << point->y() << "\n";
00265     }
00266   } else {
00267     ost << "Circle: number of points: " << points_.size() << "\n";
00268     ost << "    Radius         : " << radius_  << "\n";
00269     ost << "    In the barrel  : " << inBarrel_ << "\n";
00270     ost << "    ImpactParameter: " << impactParameter_ << "\n";
00271     ost << "    Center         : " << center_.x() << "," << center_.y() << "\n";
00272     unsigned int counter = 0;
00273     for ( std::vector<GlobalPoint>::const_iterator point = points_.begin();
00274           point != points_.end();
00275           ++point ) {
00276       ++counter;
00277       ost << "    Point " << counter << "        : " << point->x() << "," << point->y() << "\n";
00278     }
00279   }
00280 
00281   return ost.str(); 
00282 }

double RoadSearchCircleSeed::Radius (  )  const [inline]

Definition at line 63 of file RoadSearchCircleSeed.h.

References radius_.

Referenced by RoadSearchSeedFinderAlgorithm::calculateCircleSeedsFromHits(), CompareRadius(), operator<<(), and Phi0().

00063 { return radius_;}

void RoadSearchCircleSeed::setSeed ( const Roads::RoadSeed input  )  [inline]

Definition at line 69 of file RoadSearchCircleSeed.h.

References seed_.

Referenced by RoadSearchSeedFinderAlgorithm::calculateCircleSeedsFromHits().

00069 { seed_ = input; }

void RoadSearchCircleSeed::setSet ( const Roads::RoadSet input  )  [inline]

Definition at line 72 of file RoadSearchCircleSeed.h.

References set_.

Referenced by RoadSearchSeedFinderAlgorithm::calculateCircleSeedsFromHits().

00072 { set_ = input; }

double RoadSearchCircleSeed::Theta (  )  const

Definition at line 126 of file RoadSearchCircleSeed.cc.

References points_.

Referenced by Eta().

00126                                          {
00127   //
00128   // calculate the theta of the seed
00129   // by taking the average theta of all
00130   // the lines formed by combinations of
00131   // hits in the seed
00132   // 
00133   // Note:  A faster implementation would
00134   // calculate in the constructor, save, and
00135   // return the member here.  This implementation
00136   // minimizes the memory footprint.
00137   //
00138 
00139   // Form all the possible lines
00140   std::vector<LineRZ> lines;
00141   for (std::vector<GlobalPoint>::const_iterator point1 = points_.begin();
00142        point1 != points_.end(); ++point1) {
00143     for (std::vector<GlobalPoint>::const_iterator point2 = point1+1;
00144          point2 != points_.end(); ++point2) {
00145       lines.push_back(LineRZ(*point1, *point2));
00146     }
00147   }
00148   
00149   double netTheta = 0.;
00150   for (std::vector<LineRZ>::const_iterator line = lines.begin();
00151        line != lines.end(); ++line){
00152     netTheta += line->Theta();
00153   }
00154   return netTheta/(double)lines.size();
00155 }

double RoadSearchCircleSeed::Type (  )  const [inline]

Definition at line 66 of file RoadSearchCircleSeed.h.

References type_.

Referenced by RoadSearchSeedFinderAlgorithm::calculateCircleSeedsFromHits(), and operator<<().

00066 { return type_; }


Member Data Documentation

GlobalPoint RoadSearchCircleSeed::center_ [private]

Definition at line 107 of file RoadSearchCircleSeed.h.

Referenced by Center(), CompareCenter(), Phi0(), print(), and RoadSearchCircleSeed().

std::vector<const TrackingRecHit*> RoadSearchCircleSeed::hits_ [private]

Definition at line 103 of file RoadSearchCircleSeed.h.

Referenced by begin_hits(), calculateInBarrel(), CompareDifferentHits(), end_hits(), Hits(), and RoadSearchCircleSeed().

double RoadSearchCircleSeed::impactParameter_ [private]

Definition at line 109 of file RoadSearchCircleSeed.h.

Referenced by ImpactParameter(), print(), and RoadSearchCircleSeed().

bool RoadSearchCircleSeed::inBarrel_ [private]

Definition at line 106 of file RoadSearchCircleSeed.h.

Referenced by InBarrel(), print(), and RoadSearchCircleSeed().

std::vector<GlobalPoint> RoadSearchCircleSeed::points_ [private]

Definition at line 101 of file RoadSearchCircleSeed.h.

Referenced by begin_points(), end_points(), Phi0(), Points(), print(), RoadSearchCircleSeed(), and Theta().

double RoadSearchCircleSeed::radius_ [private]

Definition at line 108 of file RoadSearchCircleSeed.h.

Referenced by CompareRadius(), print(), Radius(), and RoadSearchCircleSeed().

const Roads::RoadSeed* RoadSearchCircleSeed::seed_ [private]

Definition at line 111 of file RoadSearchCircleSeed.h.

Referenced by getSeed(), and setSeed().

const Roads::RoadSet* RoadSearchCircleSeed::set_ [private]

Definition at line 112 of file RoadSearchCircleSeed.h.

Referenced by getSet(), and setSet().

type RoadSearchCircleSeed::type_ [private]

Definition at line 105 of file RoadSearchCircleSeed.h.

Referenced by Phi0(), print(), RoadSearchCircleSeed(), and Type().


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