CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/RecoPixelVertexing/PixelTriplets/plugins/ThirdHitPredictionFromInvLine.cc

Go to the documentation of this file.
00001 #include "ThirdHitPredictionFromInvLine.h"
00002 
00003 #include <cmath>
00004 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00005 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00006 
00007 #include "RecoTracker/TkMSParametrization/interface/PixelRecoUtilities.h"
00008 #include "RecoTracker/TkMSParametrization/interface/PixelRecoRange.h"
00009 
00010 template <class T> T sqr( T t) {return t*t;}
00011 
00012 typedef Basic3DVector<double> Point3D;
00013 typedef Basic2DVector<double> Point2D;
00014 typedef PixelRecoRange<double> Ranged;
00015 
00016 #include <iostream>
00017 
00018 using namespace std;
00019 
00020 ThirdHitPredictionFromInvLine::ThirdHitPredictionFromInvLine(
00021     const GlobalPoint & P1, const GlobalPoint & P2,
00022     double errorRPhiP1, double errorRPhiP2)
00023     : nPoints(0), theSum(0.),  theSumU(0.), theSumUU(0.), theSumV(0.), theSumUV(0.), theSumVV(0.),
00024       hasParameters(false),  theCurvatureValue(0.), theCurvatureError(0.), theChi2(0.)
00025 {
00026   GlobalVector aX = GlobalVector( P1.x(), P1.y(), 0.).unit();
00027   GlobalVector aY( -aX.y(), aX.x(), 0.);
00028   GlobalVector aZ( 0., 0., 1.);
00029   theRotation = Rotation(aX,aY,aZ);
00030 
00031   add(P1, errorRPhiP1);
00032   add(P2, errorRPhiP2);
00033 
00034 }
00035 
00036 GlobalPoint ThirdHitPredictionFromInvLine::crossing(double radius) const
00037 {
00038   double A = -(theSum*theSumUV-theSumU*theSumV)/(sqr(theSumU)-theSum*theSumUU);
00039   double B = (theSumU*theSumUV - theSumUU*theSumV)/(sqr(theSumU)-theSum*theSumUU);
00040   double delta =  sqr(2.*A*B) - 4*(1+sqr(A))*(sqr(B)-sqr(1/radius));
00041   double sqrtdelta = (delta > 0.) ? sqrt(delta) : 0.; 
00042   double u1 = (-2.*A*B + sqrtdelta)/2./(1+sqr(A)); 
00043   double v1 = A*u1+B;
00044   Point2D tmp = PointUV(u1,v1, &theRotation).unmap();
00045   return GlobalPoint(tmp.x(), tmp.y(), 0.);
00046 }
00047 
00048 
00049 void ThirdHitPredictionFromInvLine::add(const GlobalPoint &p, double errorRPhi)
00050 {
00051   double weigth = sqr(sqr(p.perp())/errorRPhi);
00052   add(PointUV(Point2D(p.x(),p.y()), &theRotation), weigth);
00053 }
00054 
00055 void ThirdHitPredictionFromInvLine::add(const ThirdHitPredictionFromInvLine::PointUV & point, double weigth)
00056 {
00057   hasParameters = false;
00058   nPoints++;
00059   theSum += weigth;
00060   theSumU += point.u()*weigth;
00061   theSumUU += sqr(point.u())*weigth;
00062   theSumV += point.v()*weigth;
00063   theSumUV +=  point.u()*point.v()*weigth; 
00064   theSumVV += sqr(point.v())*weigth;
00065 }
00066 
00067 void ThirdHitPredictionFromInvLine::remove(const GlobalPoint &p, double errorRPhi)
00068 {
00069   hasParameters = false;
00070   PointUV point(Point2D(p.x(),p.y()), &theRotation);
00071   double weigth = sqr(sqr(p.perp())/errorRPhi);
00072   nPoints--;
00073   theSum -= weigth;
00074   theSumU -= point.u()*weigth;
00075   theSumUU -= sqr(point.u())*weigth;
00076   theSumV -= point.v()*weigth;
00077   theSumUV -=  point.u()*point.v()*weigth;
00078   theSumVV -= sqr(point.v())*weigth;
00079 }
00080 
00081 void ThirdHitPredictionFromInvLine::print() const
00082 {
00083   std::cout <<" nPoints: " << nPoints <<" theSumU: "<< theSumU <<" theSumUU: "<<theSumUU
00084                                       <<" theSumV: "<< theSumV <<" theSumUV: "<<theSumUV<<std::endl;
00085 }
00086 
00087 void ThirdHitPredictionFromInvLine::check() const
00088 {
00089   if (hasParameters) return;
00090 
00091   long double D = theSumUU*theSum - theSumU*theSumU;
00092   long double A = (theSumUV*theSum-theSumU*theSumV)/D;
00093   long double B = (theSumUU*theSumV-theSumUV*theSumU)/D;
00094   double rho =  2.*fabs(B)/sqrt(1+sqr(A));
00095   double sigmaA2 = theSum/D;
00096   double sigmaB2 = theSumUU/D;
00097 
00098   hasParameters = true;
00099   theCurvatureError = sqrt( sqr(rho/B)*sigmaB2+ sqr(rho/(1+sqr(A)))*sigmaA2); 
00100   theCurvatureValue = 2.*fabs(B)/sqrt(1+sqr(A)); 
00101   theChi2 = theSumVV - 2*A*theSumUV - 2*B*theSumV + 2*A*B*theSumU + B*B*theSum + A*A*theSumUU;
00102 }
00103 
00104 /*
00105 GlobalPoint ThirdHitPredictionFromInvLine::center() const
00106 {
00107   long double den = theSumU*theSumUV - theSumUU*theSumV; 
00108   double a = (theSum*theSumUV-theSumU*theSumV)/2./den;
00109   double b = (sqr(theSumU)-theSum*theSumUU)/2./den;
00110   Point3D tmp = theRotation.multiplyInverse( Point2D(a,b) );
00111   return GlobalPoint(tmp.x(), tmp.y(), 0.);
00112 }
00113 */