CMS 3D CMS Logo

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

ThirdHitPredictionFromInvParabola Class Reference

#include <ThirdHitPredictionFromInvParabola.h>

List of all members.

Public Types

typedef Basic2DVector< double > Point2D
typedef PixelRecoRange< float > Range
typedef PixelRecoRange< double > RangeD
typedef TkRotation2D< double > Rotation

Public Member Functions

void init (const GlobalPoint &P1, const GlobalPoint &P2, double ip, double curv)
Range operator() (double radius, int charge) const
Range rangeRPhi (double radius, int charge) const __attribute__((optimize(3
 ThirdHitPredictionFromInvParabola (const GlobalPoint &P1, const GlobalPoint &P2, double ip, double curv, double tolerance)

Public Attributes

Range fast math

Private Member Functions

double coeffA (double impactParameter, double charge) const
double coeffB (double impactParameter, double charge) const
void findPointAtCurve (double radius, double charge, double ip, double &u, double &v) const
double ipFromCurvature (double curvature, double charge) const
double predV (double u, double ip, double charge) const
Point2D transform (Point2D const &p) const
Point2D transformBack (Point2D const &p) const

Private Attributes

double dv
double overDu
double pv
double su
RangeD theIpRangeMinus
RangeD theIpRangePlus
Rotation theRotation
double theTolerance
double u1u2

Detailed Description

Definition at line 25 of file ThirdHitPredictionFromInvParabola.h.


Member Typedef Documentation

Definition at line 32 of file ThirdHitPredictionFromInvParabola.h.

Definition at line 30 of file ThirdHitPredictionFromInvParabola.h.

Definition at line 31 of file ThirdHitPredictionFromInvParabola.h.

Definition at line 29 of file ThirdHitPredictionFromInvParabola.h.


Constructor & Destructor Documentation

ThirdHitPredictionFromInvParabola::ThirdHitPredictionFromInvParabola ( const GlobalPoint P1,
const GlobalPoint P2,
double  ip,
double  curv,
double  tolerance 
)

Definition at line 24 of file ThirdHitPredictionFromInvParabola.cc.

References abs, and init().

  : theTolerance(torlerance)
{
  init(P1,P2,ip,std::abs(curv));
}

Member Function Documentation

double ThirdHitPredictionFromInvParabola::coeffA ( double  impactParameter,
double  charge 
) const [inline, private]

Definition at line 75 of file ThirdHitPredictionFromInvParabola.h.

References overDu, pv, and u1u2.

Referenced by findPointAtCurve(), and predV().

{
  return -charge*pv*overDu - u1u2*impactParameter;
}
double ThirdHitPredictionFromInvParabola::coeffB ( double  impactParameter,
double  charge 
) const [inline, private]

Definition at line 81 of file ThirdHitPredictionFromInvParabola.h.

References dv, overDu, and su.

Referenced by findPointAtCurve(), and predV().

{
  return charge*dv*overDu - su*impactParameter;
}
void ThirdHitPredictionFromInvParabola::findPointAtCurve ( double  radius,
double  charge,
double  ip,
double &  u,
double &  v 
) const [inline, private]

Definition at line 100 of file ThirdHitPredictionFromInvParabola.h.

References funct::A, alpha, coeffA(), coeffB(), delta, alignCSCRings::r, mathSSE::sqrt(), and v.

Referenced by rangeRPhi().

{
  //
  // assume u=(1-alpha^2/2)/r v=alpha/r
  // solve qudratic equation neglecting aplha^4 term
  //
  double A = coeffA(ip,c);
  double B = coeffB(ip,c);

  // double overR = 1./r;
  double ipOverR = ip/r; // *overR;

  double delta = 1-4*(0.5*B+ipOverR)*(-B+A*r-ipOverR);
  // double sqrtdelta = (delta > 0) ? std::sqrt(delta) : 0.;
  double sqrtdelta = std::sqrt(delta);
  double alpha = (c>0)?  (-c+sqrtdelta)/(B+2*ipOverR) :  (-c-sqrtdelta)/(B+2*ipOverR);

  v = alpha;  // *overR
  double d2 = 1. - v*v;  // overR*overR - v*v
  // u = (d2 > 0) ? std::sqrt(d2) : 0.;
  u = std::sqrt(d2);

  // u,v not rotated! not multiplied by 1/r
}
void ThirdHitPredictionFromInvParabola::init ( const GlobalPoint P1,
const GlobalPoint P2,
double  ip,
double  curv 
)

Definition at line 33 of file ThirdHitPredictionFromInvParabola.cc.

References PV3DBase< T, PVType, FrameType >::basicVector(), dv, PixelRecoRange< T >::intersection(), ipFromCurvature(), overDu, p1, p2, pv, PixelRecoRange< T >::sort(), su, theIpRangeMinus, theIpRangePlus, theRotation, transform(), u1u2, Basic2DVector< T >::x(), and Basic2DVector< T >::y().

Referenced by ThirdHitPredictionFromInvParabola().

{
//  GlobalVector aX = GlobalVector( P2.x()-P1.x(), P2.y()-P1.y(), 0.).unit();
 
  Point2D p1 = P1.basicVector().xy();
  Point2D p2 = P2.basicVector().xy();
  theRotation = Rotation(p1);
  p1 = transform(p1);  // (1./P1.xy().mag(),0); 
  p2 = transform(p2);

 
  u1u2 = p1.x()*p2.x();
  overDu = 1./(p2.x() - p1.x());
  pv = p1.y()*p2.x() - p2.y()*p1.x();
  dv = p2.y() - p1.y();
  su = p2.x() + p1.x();

  RangeD ipRange(-ip, ip); 
  ipRange.sort();
  
  double ipIntyPlus = ipFromCurvature(0.,1);
  double ipCurvPlus = ipFromCurvature(curv, 1);
  double ipCurvMinus = ipFromCurvature(curv, -1);

  
  RangeD ipRangePlus(ipIntyPlus, ipCurvPlus); ipRangePlus.sort();
  RangeD ipRangeMinus(-ipIntyPlus, ipCurvMinus); ipRangeMinus.sort();

  theIpRangePlus  = ipRangePlus.intersection(ipRange);
  theIpRangeMinus = ipRangeMinus.intersection(ipRange);
}
double ThirdHitPredictionFromInvParabola::ipFromCurvature ( double  curvature,
double  charge 
) const [inline, private]

Definition at line 87 of file ThirdHitPredictionFromInvParabola.h.

References overDu, pv, and u1u2.

Referenced by init().

{
  double overU1u2 = 1./u1u2;
  double inInf = -charge*pv*overDu*overU1u2;
  return inInf-curvature*overU1u2*0.5;
}
Range ThirdHitPredictionFromInvParabola::operator() ( double  radius,
int  charge 
) const [inline]

Definition at line 39 of file ThirdHitPredictionFromInvParabola.h.

References rangeRPhi().

{ return rangeRPhi(radius,charge); } 
double ThirdHitPredictionFromInvParabola::predV ( double  u,
double  ip,
double  charge 
) const [inline, private]

Definition at line 95 of file ThirdHitPredictionFromInvParabola.h.

References coeffA(), and coeffB().

{
  return -charge*( coeffA(ip,charge) - coeffB(ip,charge)*u - ip*u*u);
}
ThirdHitPredictionFromInvParabola::Range ThirdHitPredictionFromInvParabola::rangeRPhi ( double  radius,
int  charge 
) const

Definition at line 68 of file ThirdHitPredictionFromInvParabola.cc.

References DeDxDiscriminatorTools::charge(), PixelRecoRange< T >::empty(), findPointAtCurve(), i, PixelRecoRange< T >::intersection(), PixelRecoRange< T >::max(), PixelRecoRange< T >::min(), diffTwoXMLs::r1, diffTwoXMLs::r2, TkRotation2D< T >::rotateBack(), swap(), theIpRangeMinus, theIpRangePlus, theRotation, theTolerance, and v.

Referenced by operator()().

{
  double charge = icharge; // will (icharge>0) ? 1. : -1; be faster?

  RangeD ip = (charge > 0) ? theIpRangePlus : theIpRangeMinus;


  //  it will vectorize with gcc 4.7 (with -O3 -fno-math-errno)
  double ipv[2]={ip.min(),ip.max()};
  double u[2], v[2];
  for (int i=0; i!=2; ++i)
    findPointAtCurve(radius, charge, ipv[i],u[i],v[i]);

 
  double phi1 = theRotation.rotateBack(Point2D(u[0],v[0])).barePhi();
  double phi2 = phi1+(v[1]-v[0]); 
  
  if (ip.empty()) {
    Range r1(phi1*radius-theTolerance, phi1*radius+theTolerance); 
    Range r2(phi2*radius-theTolerance, phi2*radius+theTolerance); 
    return r1.intersection(r2);
  }

  if (phi2<phi1) std::swap(phi1, phi2); 
  return Range(radius*phi1-theTolerance, radius*phi2+theTolerance);
  
}
Point2D ThirdHitPredictionFromInvParabola::transform ( Point2D const &  p) const [inline, private]

Definition at line 52 of file ThirdHitPredictionFromInvParabola.h.

References Basic2DVector< T >::mag2(), TkRotation2D< T >::rotate(), and theRotation.

Referenced by init().

                                             {
    return theRotation.rotate(p)/p.mag2();
  }
Point2D ThirdHitPredictionFromInvParabola::transformBack ( Point2D const &  p) const [inline, private]

Member Data Documentation

Definition at line 63 of file ThirdHitPredictionFromInvParabola.h.

Referenced by coeffB(), and init().

Definition at line 41 of file ThirdHitPredictionFromInvParabola.h.

Definition at line 63 of file ThirdHitPredictionFromInvParabola.h.

Referenced by coeffA(), coeffB(), init(), and ipFromCurvature().

Definition at line 63 of file ThirdHitPredictionFromInvParabola.h.

Referenced by coeffA(), init(), and ipFromCurvature().

Definition at line 63 of file ThirdHitPredictionFromInvParabola.h.

Referenced by coeffB(), and init().

Definition at line 67 of file ThirdHitPredictionFromInvParabola.h.

Referenced by init(), and rangeRPhi().

Definition at line 67 of file ThirdHitPredictionFromInvParabola.h.

Referenced by init(), and rangeRPhi().

Definition at line 62 of file ThirdHitPredictionFromInvParabola.h.

Referenced by init(), rangeRPhi(), transform(), and transformBack().

Definition at line 68 of file ThirdHitPredictionFromInvParabola.h.

Referenced by rangeRPhi().

Definition at line 63 of file ThirdHitPredictionFromInvParabola.h.

Referenced by coeffA(), init(), and ipFromCurvature().