CMS 3D CMS Logo

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

SurveyPxbImageLocalFit Class Reference

Class to hold one picture of the BPix survey and the local fit. More...

#include <SurveyPxbImageLocalFit.h>

Inheritance diagram for SurveyPxbImageLocalFit:
SurveyPxbImage

List of all members.

Public Types

typedef unsigned int count_t
typedef std::vector< coord_tfidpoint_t
typedef std::vector< value_tlocalpars_t
typedef float pede_deriv_t
typedef int pede_label_t

Public Member Functions

void doFit (const fidpoint_t &fidpointvec)
 Invoke the fit.
void doFit (const fidpoint_t &fidpointvec, const pede_label_t &label1, const pede_label_t &label2)
void doFit (value_t x1, value_t y1, value_t g1, value_t x2, value_t y2, value_t g2)
value_t getChi2 ()
 returns the chi^2 of the fit
const pede_label_tgetGlobalDerivsLabelPtr (count_t i)
const pede_deriv_tgetGlobalDerivsPtr (count_t i)
pede_label_t getGlobalDerivsSize ()
const pede_deriv_tgetLocalDerivsPtr (count_t i)
pede_label_t getLocalDerivsSize ()
localpars_t getLocalParameters ()
 returns local parameters after fit
pede_deriv_t getResiduum (count_t i)
pede_deriv_t getSigma (count_t i)
bool isFitValid ()
 returns validity flag
void setGlobalDerivsToZero (count_t i)
void setLocalDerivsToZero (count_t i)
 SurveyPxbImageLocalFit ()
 SurveyPxbImageLocalFit (std::istringstream &iss)
 Constructor from istringstream.

Static Public Attributes

static const count_t nFidpoints = 4
static const count_t nGlD = 3
static const count_t nLcD = 4
static const count_t nLcPars = 4
static const count_t nMsrmts = 8

Private Member Functions

value_t dist (coord_t p1, coord_t p2)
 Distance.
void initFidPoints ()
 Initialise the fiducial points.

Private Attributes

localpars_t a_
 Local parameters.
value_t chi2_
 chi2 of the local fit
bool derivsValidFlag_
std::vector< coord_tfidpoints_
 True position of the fiducial points on a sensor wrt. local frame (u,v)
bool fitValidFlag_
 Validity Flag.
ROOT::Math::SMatrix
< pede_deriv_t, nMsrmts, nGlD
globalDerivsMatrix_
 Matrix with global derivs.
std::vector< pede_label_tlabelVec1_
 Vector with labels to global derivs.
std::vector< pede_label_tlabelVec2_
ROOT::Math::SMatrix
< pede_deriv_t, nMsrmts, nLcD
localDerivsMatrix_
 Matrix with local derivs.
ROOT::Math::SVector< value_t,
nMsrmts
r
 Vector of residuals.

Detailed Description

Class to hold one picture of the BPix survey and the local fit.

Definition at line 14 of file SurveyPxbImageLocalFit.h.


Member Typedef Documentation

typedef unsigned int SurveyPxbImageLocalFit::count_t

Reimplemented from SurveyPxbImage.

Definition at line 19 of file SurveyPxbImageLocalFit.h.

Definition at line 18 of file SurveyPxbImageLocalFit.h.

Definition at line 17 of file SurveyPxbImageLocalFit.h.

Definition at line 27 of file SurveyPxbImageLocalFit.h.

Definition at line 26 of file SurveyPxbImageLocalFit.h.


Constructor & Destructor Documentation

SurveyPxbImageLocalFit::SurveyPxbImageLocalFit ( ) [inline]

Definition at line 30 of file SurveyPxbImageLocalFit.h.

References initFidPoints().

SurveyPxbImageLocalFit::SurveyPxbImageLocalFit ( std::istringstream &  iss) [inline]

Constructor from istringstream.

Definition at line 37 of file SurveyPxbImageLocalFit.h.

References initFidPoints().


Member Function Documentation

value_t SurveyPxbImageLocalFit::dist ( coord_t  p1,
coord_t  p2 
) [inline, private]

Distance.

Definition at line 108 of file SurveyPxbImageLocalFit.h.

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

        {
            value_t dx = p1.x()-p2.x();
            value_t dy = p1.y()-p2.y();
            return sqrt(dx*dx+dy*dy);
        }
void SurveyPxbImageLocalFit::doFit ( const fidpoint_t fidpointvec)

Invoke the fit.

Definition at line 28 of file SurveyPxbImageLocalFit.cc.

References a, funct::A, a_, chi2_, funct::cos(), gather_cfg::cout, fidpoints_, fitValidFlag_, globalDerivsMatrix_, i, localDerivsMatrix_, SurveyPxbImage::measurementVec_, r, SurveyPxbImage::sigma_x_, SurveyPxbImage::sigma_y_, funct::sin(), mathSSE::sqrt(), and detailsBasic3DVector::y.

Referenced by doFit().

{
        fitValidFlag_ = false;

        // Calculate gamma of right module w.r.t left modules' fram
        const value_t dxr = fidpointvec[3].x()-fidpointvec[2].x();
        const value_t dyr = fidpointvec[3].y()-fidpointvec[2].y();
        const value_t gammar = atan(dyr/dxr);
        const value_t dxl = fidpointvec[1].x()-fidpointvec[0].x();
        const value_t dyl = fidpointvec[1].y()-fidpointvec[0].y();
        const value_t gammal = atan(dyl/dxl);
        const value_t gamma = gammar-gammal;
        //const value_t gamma = 0.; // Testhalber
        const value_t sing = sin(gamma);
        const value_t cosg = cos(gamma);
#ifdef DEBUG
        std::cout << "gamma: " << gamma << " gamma left: " << gammal << " gamma right: " << gammar << std::endl;
#endif


#ifdef DEBUG
        std::cout << "&fidpointvec: " << std::endl;
        for (count_t i=0; i!=fidpointvec.size(); i++)
            std::cout << i << ": " << fidpointvec[i] << std::endl;
        std::cout << "&fidpoints_: " << std::endl;
        for (count_t i=0; i!=fidpoints_.size(); i++)
            std::cout << i << ": " << fidpoints_[i] << std::endl;
#endif

        // Matrix of the local derivatives
        ROOT::Math::SMatrix<value_t,nMsrmts,nLcD> A; // 8x4
        A(0,0)=1.; A(0,1)=0; A(0,2)=+fidpointvec[0].x(); A(0,3)=+fidpointvec[0].y();
        A(1,0)=0.; A(1,1)=1; A(1,2)=+fidpointvec[0].y(); A(1,3)=-fidpointvec[0].x();
        A(2,0)=1.; A(2,1)=0; A(2,2)=+fidpointvec[1].x(); A(2,3)=+fidpointvec[1].y();
        A(3,0)=0.; A(3,1)=1; A(3,2)=+fidpointvec[1].y(); A(3,3)=-fidpointvec[1].x();
        A(4,0)=1.; A(4,1)=0; A(4,2)=+fidpointvec[2].x(); A(4,3)=+fidpointvec[2].y();
        A(5,0)=0.; A(5,1)=1; A(5,2)=+fidpointvec[2].y(); A(5,3)=-fidpointvec[2].x();
        A(6,0)=1.; A(6,1)=0; A(6,2)=+fidpointvec[3].x(); A(6,3)=+fidpointvec[3].y();
        A(7,0)=0.; A(7,1)=1; A(7,2)=+fidpointvec[3].y(); A(7,3)=-fidpointvec[3].x();
#ifdef DEBUG
        std::cout << "A: \n" << A << std::endl;
#endif

        // Covariance matrix
        ROOT::Math::SMatrix<value_t,nMsrmts,nMsrmts> W; // 8x8
        //ROOT::Math::MatRepSym<value_t,8> W;
        const value_t sigma_x2inv = 1./(sigma_x_*sigma_x_);
        const value_t sigma_y2inv = 1./(sigma_y_*sigma_y_);
        W(0,0) = sigma_x2inv;
        W(1,1) = sigma_y2inv;
        W(2,2) = sigma_x2inv;
        W(3,3) = sigma_y2inv;
        W(4,4) = sigma_x2inv;
        W(5,5) = sigma_y2inv;
        W(6,6) = sigma_x2inv;
        W(7,7) = sigma_y2inv;
#ifdef DEBUG
        std::cout << "W: \n" << W << std::endl;
#endif

        // Prepare for the fit
        ROOT::Math::SMatrix<value_t,nLcD,nLcD> ATWA; // 4x4
        ATWA = ROOT::Math::Transpose(A) * W * A;
        //ATWA = ROOT::Math::SimilarityT(A,W); // W muss symmterisch sein -> aendern.
        //std::cout << "ATWA: \n" << ATWA << std::endl;
        ROOT::Math::SMatrix<value_t,nLcD,nLcD> ATWAi; // 4x4
        int ifail = 0;
        ATWAi = ATWA.Inverse(ifail);
        if (ifail != 0)
        { // TODO: ifail Pruefung auf message logger ausgeben statt cout
                std::cout << "Problem singular - fit impossible." << std::endl;
                fitValidFlag_ = false;
                return;
        }
#ifdef DEBUG
        std::cout << "ATWA-1: \n" << ATWAi << ifail << std::endl;
#endif

        // Measurements
        ROOT::Math::SVector<value_t,nMsrmts> y; // 8
        y(0) = measurementVec_[0].x();
        y(1) = measurementVec_[0].y();
        y(2) = measurementVec_[1].x();
        y(3) = measurementVec_[1].y();
        y(4) = measurementVec_[2].x();
        y(5) = measurementVec_[2].y();
        y(6) = measurementVec_[3].x();
        y(7) = measurementVec_[3].y();
#ifdef DEBUG
        std::cout << "y: " << y << std::endl;
#endif

        // do the fit
        ROOT::Math::SVector<value_t,nLcD> a; // 4
        a = ATWAi * ROOT::Math::Transpose(A) * W * y;
        chi2_ = ROOT::Math::Dot(y,W*y)-ROOT::Math::Dot(a,ROOT::Math::Transpose(A)*W*y);
#ifdef DEBUG
        std::cout << "a: " << a 
                << " S= " << sqrt(a[2]*a[2]+a[3]*a[3]) 
                << " phi= " << atan(a[3]/a[2]) 
                << " chi2= " << chi2_ << std::endl;
        std::cout << "A*a: " << A*a << std::endl;
#endif
        a_.assign(a.begin(),a.end());

        // Calculate vector of residuals
        r = y - A*a;
#ifdef DEBUG
        std::cout << "r: " << r << std::endl;
#endif

        // Fill matrix for global fit with local derivatives
        localDerivsMatrix_(0,0)=1.; localDerivsMatrix_(0,1)=0; 
        localDerivsMatrix_(1,0)=0.; localDerivsMatrix_(1,1)=1; 
        localDerivsMatrix_(2,0)=1.; localDerivsMatrix_(2,1)=0; 
        localDerivsMatrix_(3,0)=0.; localDerivsMatrix_(3,1)=1; 
        localDerivsMatrix_(4,0)=1.; localDerivsMatrix_(4,1)=0; 
        localDerivsMatrix_(5,0)=0.; localDerivsMatrix_(5,1)=1; 
        localDerivsMatrix_(6,0)=1.; localDerivsMatrix_(6,1)=0; 
        localDerivsMatrix_(7,0)=0.; localDerivsMatrix_(7,1)=1; 
        localDerivsMatrix_(0,2)=+fidpointvec[0].x()+cosg*fidpoints_[0].x()-sing*fidpoints_[0].y();
        localDerivsMatrix_(0,3)=+fidpointvec[0].y()+cosg*fidpoints_[0].y()+sing*fidpoints_[0].x();
        localDerivsMatrix_(1,2)=+fidpointvec[0].y()+cosg*fidpoints_[0].y()+sing*fidpoints_[0].x();
        localDerivsMatrix_(1,3)=-fidpointvec[0].x()-cosg*fidpoints_[0].x()+sing*fidpoints_[0].y();
        localDerivsMatrix_(2,2)=+fidpointvec[1].x()+cosg*fidpoints_[1].x()-sing*fidpoints_[1].y();
        localDerivsMatrix_(2,3)=+fidpointvec[1].y()+cosg*fidpoints_[1].y()+sing*fidpoints_[1].x();
        localDerivsMatrix_(3,2)=+fidpointvec[1].y()+cosg*fidpoints_[1].y()+sing*fidpoints_[1].x();
        localDerivsMatrix_(3,3)=-fidpointvec[1].x()-cosg*fidpoints_[1].x()+sing*fidpoints_[1].y();
        localDerivsMatrix_(4,2)=+fidpointvec[2].x()+cosg*fidpoints_[2].x()-sing*fidpoints_[2].y();
        localDerivsMatrix_(4,3)=+fidpointvec[2].y()+cosg*fidpoints_[2].y()+sing*fidpoints_[2].x();
        localDerivsMatrix_(5,2)=+fidpointvec[2].y()+cosg*fidpoints_[2].y()+sing*fidpoints_[2].x();
        localDerivsMatrix_(5,3)=-fidpointvec[2].x()-cosg*fidpoints_[2].x()+sing*fidpoints_[2].y();
        localDerivsMatrix_(6,2)=+fidpointvec[3].x()+cosg*fidpoints_[3].x()-sing*fidpoints_[3].y();
        localDerivsMatrix_(6,3)=+fidpointvec[3].y()+cosg*fidpoints_[3].y()+sing*fidpoints_[3].x();
        localDerivsMatrix_(7,2)=+fidpointvec[3].y()+cosg*fidpoints_[3].y()+sing*fidpoints_[3].x();
        localDerivsMatrix_(7,3)=-fidpointvec[3].x()-cosg*fidpoints_[3].x()+sing*fidpoints_[3].y();

        // Fill vector with global derivatives and labels (8x3)
        globalDerivsMatrix_(0,0) = +a(2);
        globalDerivsMatrix_(0,1) = +a(3);
        globalDerivsMatrix_(0,2) = +cosg*(a(3)*fidpoints_[0].x()-a(2)*fidpoints_[0].y())
                                   -sing*(a(2)*fidpoints_[0].x()+a(3)*fidpoints_[0].y());
        globalDerivsMatrix_(1,0) = -a(3);
        globalDerivsMatrix_(1,1) = +a(2);
        globalDerivsMatrix_(1,2) = +cosg*(a(2)*fidpoints_[0].x()+a(3)*fidpoints_[0].y())
                                   -sing*(a(2)*fidpoints_[0].y()-a(3)*fidpoints_[0].x());
        globalDerivsMatrix_(2,0) = +a(2);
        globalDerivsMatrix_(2,1) = +a(3);
        globalDerivsMatrix_(2,2) = +cosg*(a(3)*fidpoints_[1].x()-a(2)*fidpoints_[1].y())
                                   -sing*(a(2)*fidpoints_[1].x()+a(3)*fidpoints_[1].y());
        globalDerivsMatrix_(3,0) = -a(3);
        globalDerivsMatrix_(3,1) = +a(2);
        globalDerivsMatrix_(3,2) = +cosg*(a(2)*fidpoints_[1].x()+a(3)*fidpoints_[1].y())
                                   -sing*(a(2)*fidpoints_[1].y()-a(3)*fidpoints_[1].x());

        globalDerivsMatrix_(4,0) = +a(2);
        globalDerivsMatrix_(4,1) = +a(3);
        globalDerivsMatrix_(4,2) = +cosg*(a(3)*fidpoints_[2].x()-a(2)*fidpoints_[2].y())
                                   -sing*(a(2)*fidpoints_[2].x()+a(3)*fidpoints_[2].y());
        globalDerivsMatrix_(5,0) = -a(3);
        globalDerivsMatrix_(5,1) = +a(2);
        globalDerivsMatrix_(5,2) = +cosg*(a(2)*fidpoints_[2].x()+a(3)*fidpoints_[2].y())
                                   -sing*(a(2)*fidpoints_[2].y()-a(3)*fidpoints_[2].x());
        globalDerivsMatrix_(6,0) = +a(2);
        globalDerivsMatrix_(6,1) = +a(3);
        globalDerivsMatrix_(6,2) = +cosg*(a(3)*fidpoints_[3].x()-a(2)*fidpoints_[3].y())
                                   -sing*(a(2)*fidpoints_[3].x()+a(3)*fidpoints_[3].y());
        globalDerivsMatrix_(7,0) = -a(3);
        globalDerivsMatrix_(7,1) = +a(2);
        globalDerivsMatrix_(7,2) = +cosg*(a(2)*fidpoints_[3].x()+a(3)*fidpoints_[3].y())
                                   -sing*(a(2)*fidpoints_[3].y()-a(3)*fidpoints_[3].x());

        fitValidFlag_ = true;
}
void SurveyPxbImageLocalFit::doFit ( const fidpoint_t fidpointvec,
const pede_label_t label1,
const pede_label_t label2 
)

Definition at line 15 of file SurveyPxbImageLocalFit.cc.

References doFit(), labelVec1_, and labelVec2_.

{
        labelVec1_.clear();
        labelVec1_.push_back(label1+0);
        labelVec1_.push_back(label1+1);
        labelVec1_.push_back(label1+5);
        labelVec2_.clear();
        labelVec2_.push_back(label2+0);
        labelVec2_.push_back(label2+1);
        labelVec2_.push_back(label2+5);
        doFit(fidpointvec);
}
void SurveyPxbImageLocalFit::doFit ( value_t  x1,
value_t  y1,
value_t  g1,
value_t  x2,
value_t  y2,
value_t  g2 
)

Definition at line 204 of file SurveyPxbImageLocalFit.cc.

References funct::cos(), doFit(), fidpoints_, and funct::sin().

{
        // Creating vectors with the global parameters of the two modules
        ROOT::Math::SVector<value_t,4> mod1, mod2;
        mod1(0)=u1;
        mod1(1)=v1;
        mod1(2)=cos(g1);
        mod1(3)=sin(g1);
        mod2(0)=u2;
        mod2(1)=v2;
        mod2(2)=cos(g2);
        mod2(3)=sin(g2);
        //std::cout << "mod1: " << mod1 << std::endl;
        //std::cout << "mod2: " << mod2 << std::endl;

        // Create a matrix for the transformed position of the fidpoints
        ROOT::Math::SMatrix<value_t,4,4> M1, M2;
        M1(0,0)=1.; M1(0,1)=0.; M1(0,2)=+fidpoints_[0].x(); M1(0,3)=-fidpoints_[0].y();
        M1(1,0)=0.; M1(1,1)=1.; M1(1,2)=+fidpoints_[0].y(); M1(1,3)=+fidpoints_[0].x();
        M1(2,0)=1.; M1(2,1)=0.; M1(2,2)=+fidpoints_[1].x(); M1(2,3)=-fidpoints_[1].y();
        M1(3,0)=0.; M1(3,1)=1.; M1(3,2)=+fidpoints_[1].y(); M1(3,3)=+fidpoints_[1].x();
        M2(0,0)=1.; M2(0,1)=0.; M2(0,2)=+fidpoints_[2].x(); M2(0,3)=-fidpoints_[2].y();
        M2(1,0)=0.; M2(1,1)=1.; M2(1,2)=+fidpoints_[2].y(); M2(1,3)=+fidpoints_[2].x();
        M2(2,0)=1.; M2(2,1)=0.; M2(2,2)=+fidpoints_[3].x(); M2(2,3)=-fidpoints_[3].y();
        M2(3,0)=0.; M2(3,1)=1.; M2(3,2)=+fidpoints_[3].y(); M2(3,3)=+fidpoints_[3].x();

        //std::cout << "M1:\n" << M1 << std::endl;
        //std::cout << "M2:\n" << M2 << std::endl;

        ROOT::Math::SVector<value_t,4> mod_tr1, mod_tr2;
        mod_tr1 = M1*mod2;
        mod_tr2 = M2*mod1;
        //std::cout << "mod_tr1: " << mod_tr1 << std::endl;
        //std::cout << "mod_tr2: " << mod_tr2 << std::endl;

        fidpoint_t fidpointvec;
        fidpointvec.push_back(coord_t(mod_tr1(0),mod_tr1(1)));
        fidpointvec.push_back(coord_t(mod_tr1(2),mod_tr1(3)));
        fidpointvec.push_back(coord_t(mod_tr2(0),mod_tr2(1)));
        fidpointvec.push_back(coord_t(mod_tr2(2),mod_tr2(3)));

        doFit(fidpointvec);
}
SurveyPxbImageLocalFit::value_t SurveyPxbImageLocalFit::getChi2 ( )

returns the chi^2 of the fit

Definition at line 254 of file SurveyPxbImageLocalFit.cc.

References chi2_, and fitValidFlag_.

{
        if (!fitValidFlag_) throw std::logic_error("SurveyPxbImageLocalFit::getChi2(): Fit is not valid. Call doFit(...) before calling this function.");
        return chi2_;
}
const pede_label_t* SurveyPxbImageLocalFit::getGlobalDerivsLabelPtr ( count_t  i) [inline]

Definition at line 61 of file SurveyPxbImageLocalFit.h.

References labelVec1_, and labelVec2_.

{return i<4 ? &labelVec1_[0] : &labelVec2_[0];};
const pede_deriv_t* SurveyPxbImageLocalFit::getGlobalDerivsPtr ( count_t  i) [inline]

Definition at line 60 of file SurveyPxbImageLocalFit.h.

References globalDerivsMatrix_, and nGlD.

{return globalDerivsMatrix_.Array()+i*nGlD;};
pede_label_t SurveyPxbImageLocalFit::getGlobalDerivsSize ( ) [inline]

Definition at line 58 of file SurveyPxbImageLocalFit.h.

References nGlD.

{return nGlD; };
const pede_deriv_t* SurveyPxbImageLocalFit::getLocalDerivsPtr ( count_t  i) [inline]

Definition at line 59 of file SurveyPxbImageLocalFit.h.

References localDerivsMatrix_, and nLcD.

{return localDerivsMatrix_.Array()+i*nLcD; };
pede_label_t SurveyPxbImageLocalFit::getLocalDerivsSize ( ) [inline]

Definition at line 57 of file SurveyPxbImageLocalFit.h.

References nLcD.

{return nLcD; };
SurveyPxbImageLocalFit::localpars_t SurveyPxbImageLocalFit::getLocalParameters ( )

returns local parameters after fit

Definition at line 248 of file SurveyPxbImageLocalFit.cc.

References a_, and fitValidFlag_.

{
        if (!fitValidFlag_) throw std::logic_error("SurveyPxbImageLocalFit::getLocalParameters(): Fit is not valid. Call doFit(...) before calling this function.");
        return a_;
}
pede_deriv_t SurveyPxbImageLocalFit::getResiduum ( count_t  i) [inline]

Definition at line 62 of file SurveyPxbImageLocalFit.h.

References r.

{ return (pede_deriv_t) r(i); };
pede_deriv_t SurveyPxbImageLocalFit::getSigma ( count_t  i) [inline]

Definition at line 63 of file SurveyPxbImageLocalFit.h.

References SurveyPxbImage::sigma_x_, and SurveyPxbImage::sigma_y_.

{ return i%2 ? sigma_x_ : sigma_y_ ; };
void SurveyPxbImageLocalFit::initFidPoints ( ) [inline, private]

Initialise the fiducial points.

Definition at line 99 of file SurveyPxbImageLocalFit.h.

References fidpoints_.

Referenced by SurveyPxbImageLocalFit().

        {
                fidpoints_[0] = coord_t(-0.91,-3.30);
                fidpoints_[1] = coord_t(+0.91,-3.30);
                fidpoints_[2] = coord_t(-0.91,+3.30);
                fidpoints_[3] = coord_t(+0.91,+3.30);
        }
bool SurveyPxbImageLocalFit::isFitValid ( ) [inline]

returns validity flag

Definition at line 49 of file SurveyPxbImageLocalFit.h.

References fitValidFlag_.

{ return fitValidFlag_; };
void SurveyPxbImageLocalFit::setGlobalDerivsToZero ( count_t  i)

Definition at line 267 of file SurveyPxbImageLocalFit.cc.

References globalDerivsMatrix_, i, nGlD, and nMsrmts.

{
    if (!(j < nGlD)) throw std::range_error("SurveyPxbImageLocalFit::setLocalDerivsToZero(j): j out of range.");
    for(count_t i=0; i!=nMsrmts; i++)
        globalDerivsMatrix_(i,j)=0;
}
void SurveyPxbImageLocalFit::setLocalDerivsToZero ( count_t  i)

Definition at line 260 of file SurveyPxbImageLocalFit.cc.

References i, localDerivsMatrix_, nLcD, and nMsrmts.

{
    if (!(j < nLcD)) throw std::range_error("SurveyPxbImageLocalFit::setLocalDerivsToZero(j): j out of range.");
    for(count_t i=0; i!=nMsrmts; i++)
        localDerivsMatrix_(i,j)=0;
}

Member Data Documentation

Local parameters.

Definition at line 70 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), and getLocalParameters().

chi2 of the local fit

Definition at line 90 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), and getChi2().

Definition at line 96 of file SurveyPxbImageLocalFit.h.

True position of the fiducial points on a sensor wrt. local frame (u,v)

Definition at line 87 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), and initFidPoints().

Validity Flag.

Definition at line 93 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), getChi2(), getLocalParameters(), and isFitValid().

Matrix with global derivs.

Definition at line 77 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), getGlobalDerivsPtr(), and setGlobalDerivsToZero().

Vector with labels to global derivs.

Definition at line 84 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), and getGlobalDerivsLabelPtr().

Definition at line 84 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), and getGlobalDerivsLabelPtr().

Matrix with local derivs.

Definition at line 81 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), getLocalDerivsPtr(), and setLocalDerivsToZero().

Definition at line 24 of file SurveyPxbImageLocalFit.h.

Definition at line 23 of file SurveyPxbImageLocalFit.h.

ROOT::Math::SVector<value_t, nMsrmts> SurveyPxbImageLocalFit::r [private]

Vector of residuals.

Definition at line 73 of file SurveyPxbImageLocalFit.h.

Referenced by doFit(), and getResiduum().