MOAB: Mesh Oriented datABase  (version 5.4.1)
MBMesquite::Matrix3D Class Reference

3*3 Matric class, row-oriented, 0-based [i][j] indexing. More...

#include <Matrix3D.hpp>

Public Member Functions

 Matrix3D ()
 Default constructor sets all entries to 0.
 Matrix3D (const Matrix3D &A)
 Matrix3D (double value)
 sets all entries of the matrix to value.
 Matrix3D (double a00, double a01, double a02, double a10, double a11, double a12, double a20, double a21, double a22)
 Matrix3D (const Vector3D &col1, const Vector3D &col2, const Vector3D &col3)
 Matrix3D (double radians, const Vector3D &axis)
 Matrix3D (const double *v)
 Matrix3D (const char *s)
 Matrix3D (const SymMatrix3D &m)
 ~Matrix3D ()
Matrix3Doperator= (const Matrix3D &A)
Matrix3Doperator= (const SymMatrix3D &m)
Matrix3Doperator= (double scalar)
Matrix3Doperator= (const char *s)
void zero ()
 Sets all entries to zero (more efficient than assignement).
void identity ()
void set_column (int j, const Vector3D &c)
 Sets column j (0, 1 or 2) to Vector3D c.
double column_length (int i) const
 returns the column length -- i is 0-based.
double sub_det (int r, int c) const
Matrix3Dtranspose ()
Matrix3Dequal_mult_elem (const Matrix3D &A)
Matrix3Dassign_product (const Matrix3D &A, const Matrix3D &B)
const Matrix3D operator* (double s) const
 multiplies each entry by the scalar s
void operator+= (const Matrix3D &rhs)
void operator+= (const SymMatrix3D &rhs)
void operator-= (const Matrix3D &rhs)
void operator-= (const SymMatrix3D &rhs)
void operator*= (double s)
 multiplies each entry by the scalar s
Matrix3Dplus_transpose_equal (const Matrix3D &B)
 \( += B^T \)
Matrix3Douter_product (const Vector3D &v1, const Vector3D &v2)
 Computes \( A = v_1 v_2^T \).
void fill_lower_triangle ()
size_t num_rows () const
size_t num_cols () const
double * operator[] (unsigned i)
 returns a pointer to a row.
const double * operator[] (unsigned i) const
 returns a pointer to a row.
double & operator() (unsigned short r, unsigned short c)
double operator() (unsigned short r, unsigned short c) const
Vector3D row (unsigned r) const
Vector3D column (unsigned c) const
bool positive_definite () const
SymMatrix3D upper () const
SymMatrix3D lower () const

Protected Member Functions

void copy (const double *v)
void set (double val)
void set_values (const char *s)

Protected Attributes

double v_ [9]

Friends

bool operator== (const Matrix3D &lhs, const Matrix3D &rhs)
bool operator!= (const Matrix3D &lhs, const Matrix3D &rhs)
Matrix3D operator- (const Matrix3D &A)
double Frobenius_2 (const Matrix3D &A)
 Return the square of the Frobenius norm of A, i.e. sum (diag (A' * A))
Matrix3D transpose (const Matrix3D &A)
const Matrix3D operator+ (const Matrix3D &A, const Matrix3D &B)
const Matrix3D operator- (const Matrix3D &A, const Matrix3D &B)
const Matrix3D operator* (const Matrix3D &A, const Matrix3D &B)
const Matrix3D mult_element (const Matrix3D &A, const Matrix3D &B)
 Multiplies entry by entry. This is NOT a matrix multiplication.
void matmult (Matrix3D &C, const Matrix3D &A, const Matrix3D &B)
 \( C = A \times B \)
const Vector3D operator* (const Matrix3D &A, const Vector3D &x)
 Computes \( A v \) .
const Vector3D operator* (const Vector3D &x, const Matrix3D &A)
 Computes \( v^T A \) .
const Matrix3D operator* (double s, const Matrix3D &A)
Matrix3D plus_transpose (const Matrix3D &A, const Matrix3D &B)
 \( + B^T \)
void eqAx (Vector3D &v, const Matrix3D &A, const Vector3D &x)
 \( v = A*x \)
void plusEqAx (Vector3D &v, const Matrix3D &A, const Vector3D &x)
 \( v += A*x \)
void eqTransAx (Vector3D &v, const Matrix3D &A, const Vector3D &x)
void plusEqTransAx (Vector3D &v, const Matrix3D &A, const Vector3D &x)
 \( v += A^T*x \)
void plusEqaA (Matrix3D &B, const double a, const Matrix3D &A)
 \( B += a*A \)
double det (const Matrix3D &A)
 determinant of matrix A, det(A).
void inv (Matrix3D &B, const Matrix3D &A)
 \( B = A^{-1} \)
void timesInvA (Matrix3D &B, const Matrix3D &A)
 \( B *= A^{-1} \)
void QR (Matrix3D &Q, Matrix3D &R, const Matrix3D &A)
 \( Q*R = A \)

Detailed Description

3*3 Matric class, row-oriented, 0-based [i][j] indexing.

Since the size of the object is fixed at compile time, the Matrix3D object is as fast as a double[9] array.

Definition at line 67 of file Matrix3D.hpp.


Constructor & Destructor Documentation

Default constructor sets all entries to 0.

Definition at line 103 of file Matrix3D.hpp.

    {
        zero();
    }
MBMesquite::Matrix3D::Matrix3D ( const Matrix3D A) [inline]

Definition at line 108 of file Matrix3D.hpp.

References v_.

    {
        copy( A.v_ );
    }
MBMesquite::Matrix3D::Matrix3D ( double  value) [inline]

sets all entries of the matrix to value.

Definition at line 114 of file Matrix3D.hpp.

References value().

    {
        set( value );
    }
MBMesquite::Matrix3D::Matrix3D ( double  a00,
double  a01,
double  a02,
double  a10,
double  a11,
double  a12,
double  a20,
double  a21,
double  a22 
) [inline]

Definition at line 119 of file Matrix3D.hpp.

    {
        v_[0] = a00;
        v_[1] = a01;
        v_[2] = a02;
        v_[3] = a10;
        v_[4] = a11;
        v_[5] = a12;
        v_[6] = a20;
        v_[7] = a21;
        v_[8] = a22;
    }
MBMesquite::Matrix3D::Matrix3D ( const Vector3D col1,
const Vector3D col2,
const Vector3D col3 
) [inline]

Definition at line 140 of file Matrix3D.hpp.

    {
        set_column( 0, col1 );
        set_column( 1, col2 );
        set_column( 2, col3 );
    }
MBMesquite::Matrix3D::Matrix3D ( double  radians,
const Vector3D axis 
) [inline]

Definition at line 147 of file Matrix3D.hpp.

References MBMesquite::Vector3D::normalize().

    {
        Vector3D v( axis );
        v.normalize();
        const double c = std::cos( radians );
        const double s = std::sin( radians );
        v_[0]          = c + ( 1.0 - c ) * v[0] * v[0];
        v_[1]          = -v[2] * s + ( 1.0 - c ) * v[0] * v[1];
        v_[2]          = v[1] * s + ( 1.0 - c ) * v[0] * v[2];
        v_[3]          = v[2] * s + ( 1.0 - c ) * v[0] * v[1];
        v_[4]          = c + ( 1.0 - c ) * v[1] * v[1];
        v_[5]          = -v[0] * s + ( 1.0 - c ) * v[1] * v[2];
        v_[6]          = -v[1] * s + ( 1.0 - c ) * v[0] * v[2];
        v_[7]          = v[0] * s + ( 1.0 - c ) * v[1] * v[2];
        v_[8]          = c + ( 1.0 - c ) * v[2] * v[2];
    }
MBMesquite::Matrix3D::Matrix3D ( const double *  v) [inline]

sets matrix entries to values in array.

Parameters:
vis an array of 9 doubles.

Definition at line 166 of file Matrix3D.hpp.

    {
        copy( v );
    }
MBMesquite::Matrix3D::Matrix3D ( const char *  s) [inline]

for test purposes, matrices can be instantiated as

 Matrix3D A("3 2 1  4 5 6  9 8 7"); 

Definition at line 173 of file Matrix3D.hpp.

    {
        set_values( s );
    }

Definition at line 178 of file Matrix3D.hpp.

    {
        *this = m;
    }

Definition at line 184 of file Matrix3D.hpp.

{}

Member Function Documentation

Matrix3D & MBMesquite::Matrix3D::assign_product ( const Matrix3D A,
const Matrix3D B 
) [inline]

Definition at line 724 of file Matrix3D.hpp.

References v_.

Referenced by MBMesquite::matmult(), and MBMesquite::operator*().

{
    v_[0] = A.v_[0] * B.v_[0] + A.v_[1] * B.v_[3] + A.v_[2] * B.v_[6];
    v_[1] = A.v_[0] * B.v_[1] + A.v_[1] * B.v_[4] + A.v_[2] * B.v_[7];
    v_[2] = A.v_[0] * B.v_[2] + A.v_[1] * B.v_[5] + A.v_[2] * B.v_[8];
    v_[3] = A.v_[3] * B.v_[0] + A.v_[4] * B.v_[3] + A.v_[5] * B.v_[6];
    v_[4] = A.v_[3] * B.v_[1] + A.v_[4] * B.v_[4] + A.v_[5] * B.v_[7];
    v_[5] = A.v_[3] * B.v_[2] + A.v_[4] * B.v_[5] + A.v_[5] * B.v_[8];
    v_[6] = A.v_[6] * B.v_[0] + A.v_[7] * B.v_[3] + A.v_[8] * B.v_[6];
    v_[7] = A.v_[6] * B.v_[1] + A.v_[7] * B.v_[4] + A.v_[8] * B.v_[7];
    v_[8] = A.v_[6] * B.v_[2] + A.v_[7] * B.v_[5] + A.v_[8] * B.v_[8];
    return *this;
}
Vector3D MBMesquite::Matrix3D::column ( unsigned  c) const [inline]

Definition at line 354 of file Matrix3D.hpp.

Referenced by TMPDerivsTest::test_set_scaled_sum_outer_product_2D(), and TMPDerivsTest::test_set_scaled_sum_outer_product_3D().

    {
        return Vector3D( v_[c], v_[c + 3], v_[c + 6] );
    }
double MBMesquite::Matrix3D::column_length ( int  i) const [inline]

returns the column length -- i is 0-based.

Definition at line 254 of file Matrix3D.hpp.

    {
        return sqrt( v_[0 + i] * v_[0 + i] + v_[3 + i] * v_[3 + i] + v_[6 + i] * v_[6 + i] );
    }
void MBMesquite::Matrix3D::copy ( const double *  v) [inline, protected]

Definition at line 72 of file Matrix3D.hpp.

    {
        v_[0] = v[0];
        v_[1] = v[1];
        v_[2] = v[2];
        v_[3] = v[3];
        v_[4] = v[4];
        v_[5] = v[5];
        v_[6] = v[6];
        v_[7] = v[7];
        v_[8] = v[8];
    }

Definition at line 462 of file Matrix3D.hpp.

References v_.

Referenced by MBMesquite::mult_element().

{
    v_[0] *= A.v_[0];
    v_[1] *= A.v_[1];
    v_[2] *= A.v_[2];
    v_[3] *= A.v_[3];
    v_[4] *= A.v_[4];
    v_[5] *= A.v_[5];
    v_[6] *= A.v_[6];
    v_[7] *= A.v_[7];
    v_[8] *= A.v_[8];
    return *this;
}
void MBMesquite::Matrix3D::identity ( ) [inline]

Definition at line 232 of file Matrix3D.hpp.

    {
        v_[0] = 1.;
        v_[1] = 0.;
        v_[2] = 0.;
        v_[3] = 0.;
        v_[4] = 1.;
        v_[5] = 0.;
        v_[6] = 0.;
        v_[7] = 0.;
        v_[8] = 1.;
    }

Definition at line 365 of file Matrix3D.hpp.

    {
        return SymMatrix3D( v_[0], v_[3], v_[6], v_[4], v_[7], v_[8] );
    }
size_t MBMesquite::Matrix3D::num_cols ( ) const [inline]

Definition at line 323 of file Matrix3D.hpp.

    {
        return 3;
    }
size_t MBMesquite::Matrix3D::num_rows ( ) const [inline]

Definition at line 319 of file Matrix3D.hpp.

    {
        return 3;
    }
double& MBMesquite::Matrix3D::operator() ( unsigned short  r,
unsigned short  c 
) [inline]

Definition at line 340 of file Matrix3D.hpp.

    {
        return v_[3 * r + c];
    }
double MBMesquite::Matrix3D::operator() ( unsigned short  r,
unsigned short  c 
) const [inline]

Definition at line 344 of file Matrix3D.hpp.

    {
        return v_[3 * r + c];
    }
const Matrix3D MBMesquite::Matrix3D::operator* ( double  s) const [inline]

multiplies each entry by the scalar s

Definition at line 711 of file Matrix3D.hpp.

{
    Matrix3D temp( *this );
    temp *= s;
    return temp;
}
void MBMesquite::Matrix3D::operator*= ( double  s) [inline]

multiplies each entry by the scalar s

Definition at line 581 of file Matrix3D.hpp.

References v_.

{
    v_[0] *= s;
    v_[1] *= s;
    v_[2] *= s;
    v_[3] *= s;
    v_[4] *= s;
    v_[5] *= s;
    v_[6] *= s;
    v_[7] *= s;
    v_[8] *= s;
}
void MBMesquite::Matrix3D::operator+= ( const Matrix3D rhs) [inline]

Definition at line 528 of file Matrix3D.hpp.

References v_.

{
    v_[0] += rhs.v_[0];
    v_[1] += rhs.v_[1];
    v_[2] += rhs.v_[2];
    v_[3] += rhs.v_[3];
    v_[4] += rhs.v_[4];
    v_[5] += rhs.v_[5];
    v_[6] += rhs.v_[6];
    v_[7] += rhs.v_[7];
    v_[8] += rhs.v_[8];
}
void MBMesquite::Matrix3D::operator+= ( const SymMatrix3D rhs) [inline]

Definition at line 541 of file Matrix3D.hpp.

References v_.

{
    v_[0] += rhs[0];
    v_[1] += rhs[1];
    v_[2] += rhs[2];
    v_[3] += rhs[1];
    v_[4] += rhs[3];
    v_[5] += rhs[4];
    v_[6] += rhs[2];
    v_[7] += rhs[4];
    v_[8] += rhs[5];
}
void MBMesquite::Matrix3D::operator-= ( const Matrix3D rhs) [inline]

Definition at line 554 of file Matrix3D.hpp.

References v_.

{
    v_[0] -= rhs.v_[0];
    v_[1] -= rhs.v_[1];
    v_[2] -= rhs.v_[2];
    v_[3] -= rhs.v_[3];
    v_[4] -= rhs.v_[4];
    v_[5] -= rhs.v_[5];
    v_[6] -= rhs.v_[6];
    v_[7] -= rhs.v_[7];
    v_[8] -= rhs.v_[8];
}
void MBMesquite::Matrix3D::operator-= ( const SymMatrix3D rhs) [inline]

Definition at line 567 of file Matrix3D.hpp.

References v_.

{
    v_[0] -= rhs[0];
    v_[1] -= rhs[1];
    v_[2] -= rhs[2];
    v_[3] -= rhs[1];
    v_[4] -= rhs[3];
    v_[5] -= rhs[4];
    v_[6] -= rhs[2];
    v_[7] -= rhs[4];
    v_[8] -= rhs[5];
}
Matrix3D& MBMesquite::Matrix3D::operator= ( const Matrix3D A) [inline]

Definition at line 187 of file Matrix3D.hpp.

References v_.

    {
        copy( A.v_ );
        return *this;
    }
Matrix3D& MBMesquite::Matrix3D::operator= ( const SymMatrix3D m) [inline]

Definition at line 193 of file Matrix3D.hpp.

    {
        v_[0] = m[0];
        v_[1] = v_[3] = m[1];
        v_[2] = v_[6] = m[2];
        v_[4]         = m[3];
        v_[5] = v_[7] = m[4];
        v_[8]         = m[5];
        return *this;
    }
Matrix3D& MBMesquite::Matrix3D::operator= ( double  scalar) [inline]

Definition at line 204 of file Matrix3D.hpp.

    {
        set( scalar );
        return *this;
    }
Matrix3D& MBMesquite::Matrix3D::operator= ( const char *  s) [inline]

for test purposes, matrices can be assigned as follows

 A = "3 2 1  4 5 6  9 8 7"; 

Definition at line 212 of file Matrix3D.hpp.

    {
        set_values( s );
        return *this;
    }
double* MBMesquite::Matrix3D::operator[] ( unsigned  i) [inline]

returns a pointer to a row.

Definition at line 329 of file Matrix3D.hpp.

    {
        return v_ + 3 * i;
    }
const double* MBMesquite::Matrix3D::operator[] ( unsigned  i) const [inline]

returns a pointer to a row.

Definition at line 335 of file Matrix3D.hpp.

    {
        return v_ + 3 * i;
    }
Matrix3D & MBMesquite::Matrix3D::outer_product ( const Vector3D v1,
const Vector3D v2 
) [inline]

Computes \( A = v_1 v_2^T \).

Definition at line 635 of file Matrix3D.hpp.

References v_.

Referenced by MBMesquite::MeshTransform::add_rotation(), MBMesquite::PMeanPMetric::average_with_Hessian(), PMeanPTemplateTest::check_result(), MBMesquite::LPtoPTemplate::evaluate_with_Hessian(), MBMesquite::PowerQualityMetric::evaluate_with_Hessian(), MBMesquite::MultiplyQualityMetric::evaluate_with_Hessian(), MBMesquite::pmean_corner_hessians(), pmean_of_triangle_corner_hessians(), MBMesquite::sum_sqr_corner_hessians(), PMeanPMetricTest::test_hessian(), CompositeOFTest::test_multiply_hess_diagonal(), and Matrix3DTest::test_outer_product().

{
    // remember, matrix entries are v_[0] to v_[8].

    // diagonal
    v_[0] = v1[0] * v2[0];
    v_[4] = v1[1] * v2[1];
    v_[8] = v1[2] * v2[2];

    // upper triangular part
    v_[1] = v1[0] * v2[1];
    v_[2] = v1[0] * v2[2];
    v_[5] = v1[1] * v2[2];

    // lower triangular part
    v_[3] = v2[0] * v1[1];
    v_[6] = v2[0] * v1[2];
    v_[7] = v2[1] * v1[2];

    return *this;
}

\( += B^T \)

Definition at line 595 of file Matrix3D.hpp.

References v_.

Referenced by MBMesquite::MsqHessian::add(), MBMesquite::MultiplyQualityMetric::evaluate_with_Hessian(), MBMesquite::plus_transpose(), and positive_definite().

{
    if( &b == this )
    {
        v_[0] *= 2.0;
        v_[1] += v_[3];
        v_[2] += v_[6];
        v_[3] = v_[1];
        v_[4] *= 2.0;
        v_[5] += v_[7];
        v_[6] = v_[2];
        v_[7] = v_[5];
        v_[8] *= 2.0;
    }
    else
    {
        v_[0] += b.v_[0];
        v_[1] += b.v_[3];
        v_[2] += b.v_[6];

        v_[3] += b.v_[1];
        v_[4] += b.v_[4];
        v_[5] += b.v_[7];

        v_[6] += b.v_[2];
        v_[7] += b.v_[5];
        v_[8] += b.v_[8];
    }
    return *this;
}
bool MBMesquite::Matrix3D::positive_definite ( ) const [inline]

Definition at line 894 of file Matrix3D.hpp.

References det, plus_transpose_equal(), and sub_det().

{
    // A = B + C
    // where
    // B = (A + transpose(A))/2
    // C = (A - transpose(A))/2
    // B is always a symmetric matrix and
    // A is positive definite iff B is positive definite.
    Matrix3D B( *this );
    B.plus_transpose_equal( *this );
    B *= 0.5;

    // Sylvester's Criterion for positive definite symmetric matrix
    return ( B[0][0] > 0.0 ) && ( B.sub_det( 2, 2 ) > 0.0 ) && ( det( B ) > 0.0 );
}
void MBMesquite::Matrix3D::set ( double  val) [inline, protected]

Definition at line 85 of file Matrix3D.hpp.

    {
        v_[0] = val;
        v_[1] = val;
        v_[2] = val;
        v_[3] = val;
        v_[4] = val;
        v_[5] = val;
        v_[6] = val;
        v_[7] = val;
        v_[8] = val;
    }
void MBMesquite::Matrix3D::set_column ( int  j,
const Vector3D c 
) [inline]

Sets column j (0, 1 or 2) to Vector3D c.

Definition at line 246 of file Matrix3D.hpp.

    {
        v_[0 + j] = c[0];
        v_[3 + j] = c[1];
        v_[6 + j] = c[2];
    }
void MBMesquite::Matrix3D::set_values ( const char *  s) [inline, protected]

Definition at line 394 of file Matrix3D.hpp.

{
    std::istringstream ins( s );
    ins >> *this;
}
double MBMesquite::Matrix3D::sub_det ( int  r,
int  c 
) const [inline]

Definition at line 259 of file Matrix3D.hpp.

Referenced by positive_definite().

    {
        int r1 = 3 * ( ( r + 1 ) % 3 );
        int r2 = 3 * ( ( r + 2 ) % 3 );
        int c1 = ( ( c + 1 ) % 3 );
        int c2 = ( ( c + 2 ) % 3 );
        return v_[r1 + c1] * v_[r2 + c2] - v_[r2 + c1] * v_[r1 + c2];
    }

Definition at line 491 of file Matrix3D.hpp.

References t, and v_.

{
    double t;
    t     = v_[1];
    v_[1] = v_[3];
    v_[3] = t;
    t     = v_[2];
    v_[2] = v_[6];
    v_[6] = t;
    t     = v_[5];
    v_[5] = v_[7];
    v_[7] = t;
    return *this;
}
void MBMesquite::Matrix3D::zero ( ) [inline]

Sets all entries to zero (more efficient than assignement).

Definition at line 219 of file Matrix3D.hpp.

Referenced by MBMesquite::MsqHessian::zero_out().

    {
        v_[0] = 0.;
        v_[1] = 0.;
        v_[2] = 0.;
        v_[3] = 0.;
        v_[4] = 0.;
        v_[5] = 0.;
        v_[6] = 0.;
        v_[7] = 0.;
        v_[8] = 0.;
    }

Friends And Related Function Documentation

double det ( const Matrix3D A) [friend]

determinant of matrix A, det(A).

Definition at line 805 of file Matrix3D.hpp.

Referenced by positive_definite().

{
    return ( A.v_[0] * ( A.v_[4] * A.v_[8] - A.v_[7] * A.v_[5] ) - A.v_[1] * ( A.v_[3] * A.v_[8] - A.v_[6] * A.v_[5] ) +
             A.v_[2] * ( A.v_[3] * A.v_[7] - A.v_[6] * A.v_[4] ) );
}
void eqAx ( Vector3D v,
const Matrix3D A,
const Vector3D x 
) [friend]

\( v = A*x \)

Definition at line 764 of file Matrix3D.hpp.

{
    v.mCoords[0] = A.v_[0] * x[0] + A.v_[1] * x.mCoords[1] + A.v_[2] * x.mCoords[2];
    v.mCoords[1] = A.v_[3] * x[0] + A.v_[4] * x.mCoords[1] + A.v_[5] * x.mCoords[2];
    v.mCoords[2] = A.v_[6] * x[0] + A.v_[7] * x.mCoords[1] + A.v_[8] * x.mCoords[2];
}
void eqTransAx ( Vector3D v,
const Matrix3D A,
const Vector3D x 
) [friend]

Definition at line 778 of file Matrix3D.hpp.

{
    v.mCoords[0] = A.v_[0] * x.mCoords[0] + A.v_[3] * x.mCoords[1] + A.v_[6] * x.mCoords[2];
    v.mCoords[1] = A.v_[1] * x.mCoords[0] + A.v_[4] * x.mCoords[1] + A.v_[7] * x.mCoords[2];
    v.mCoords[2] = A.v_[2] * x.mCoords[0] + A.v_[5] * x.mCoords[1] + A.v_[8] * x.mCoords[2];
}
double Frobenius_2 ( const Matrix3D A) [friend]

Return the square of the Frobenius norm of A, i.e. sum (diag (A' * A))

Definition at line 485 of file Matrix3D.hpp.

{
    return A.v_[0] * A.v_[0] + A.v_[1] * A.v_[1] + A.v_[2] * A.v_[2] + A.v_[3] * A.v_[3] + A.v_[4] * A.v_[4] +
           A.v_[5] * A.v_[5] + A.v_[6] * A.v_[6] + A.v_[7] * A.v_[7] + A.v_[8] * A.v_[8];
}
void inv ( Matrix3D B,
const Matrix3D A 
) [friend]

\( B = A^{-1} \)

Definition at line 811 of file Matrix3D.hpp.

{
    double inv_detA = 1.0 / ( det( A ) );
    // First row of Ainv
    Ainv.v_[0] = inv_detA * ( A.v_[4] * A.v_[8] - A.v_[5] * A.v_[7] );
    Ainv.v_[1] = inv_detA * ( A.v_[2] * A.v_[7] - A.v_[8] * A.v_[1] );
    Ainv.v_[2] = inv_detA * ( A.v_[1] * A.v_[5] - A.v_[4] * A.v_[2] );
    // Second row of Ainv
    Ainv.v_[3] = inv_detA * ( A.v_[5] * A.v_[6] - A.v_[8] * A.v_[3] );
    Ainv.v_[4] = inv_detA * ( A.v_[0] * A.v_[8] - A.v_[6] * A.v_[2] );
    Ainv.v_[5] = inv_detA * ( A.v_[2] * A.v_[3] - A.v_[5] * A.v_[0] );
    // Third row of Ainv
    Ainv.v_[6] = inv_detA * ( A.v_[3] * A.v_[7] - A.v_[6] * A.v_[4] );
    Ainv.v_[7] = inv_detA * ( A.v_[1] * A.v_[6] - A.v_[7] * A.v_[0] );
    Ainv.v_[8] = inv_detA * ( A.v_[0] * A.v_[4] - A.v_[3] * A.v_[1] );
}
void matmult ( Matrix3D C,
const Matrix3D A,
const Matrix3D B 
) [friend]

\( C = A \times B \)

Definition at line 739 of file Matrix3D.hpp.

{
    C.assign_product( A, B );
}
const Matrix3D mult_element ( const Matrix3D A,
const Matrix3D B 
) [friend]

Multiplies entry by entry. This is NOT a matrix multiplication.

Definition at line 477 of file Matrix3D.hpp.

{
    Matrix3D tmp( A );
    tmp.equal_mult_elem( B );
    return tmp;
}
bool operator!= ( const Matrix3D lhs,
const Matrix3D rhs 
) [friend]

Definition at line 409 of file Matrix3D.hpp.

{
    return !( lhs == rhs );
}
const Matrix3D operator* ( const Matrix3D A,
const Matrix3D B 
) [friend]
Returns:
A*B

Definition at line 665 of file Matrix3D.hpp.

{
    Matrix3D tmp;
    tmp.assign_product( A, B );
    return tmp;
}
const Vector3D operator* ( const Matrix3D A,
const Vector3D x 
) [friend]

Computes \( A v \) .

Definition at line 745 of file Matrix3D.hpp.

{
    Vector3D tmp;
    eqAx( tmp, A, x );
    return tmp;
}
const Vector3D operator* ( const Vector3D x,
const Matrix3D A 
) [friend]

Computes \( v^T A \) .

This function implicitly considers the transpose of vector x times the matrix A and it is implicit that the returned vector must be transposed.

Definition at line 757 of file Matrix3D.hpp.

{
    Vector3D tmp;
    eqTransAx( tmp, A, x );
    return tmp;
}
const Matrix3D operator* ( double  s,
const Matrix3D A 
) [friend]

friend function to allow for commutatative property of scalar mulitplication.

Definition at line 719 of file Matrix3D.hpp.

{
    return ( A.operator*( s ) );
}
const Matrix3D operator+ ( const Matrix3D A,
const Matrix3D B 
) [friend]
Returns:
A+B

Definition at line 420 of file Matrix3D.hpp.

{
    Matrix3D tmp( A );
    tmp += B;
    return tmp;
}
Matrix3D operator- ( const Matrix3D A) [friend]

Definition at line 414 of file Matrix3D.hpp.

{
    return Matrix3D( -A.v_[0], -A.v_[1], -A.v_[2], -A.v_[3], -A.v_[4], -A.v_[5], -A.v_[6], -A.v_[7], -A.v_[8] );
}
const Matrix3D operator- ( const Matrix3D A,
const Matrix3D B 
) [friend]
Returns:
A-B

Definition at line 440 of file Matrix3D.hpp.

{
    Matrix3D tmp( A );
    tmp -= B;
    return tmp;
}
bool operator== ( const Matrix3D lhs,
const Matrix3D rhs 
) [friend]

Definition at line 403 of file Matrix3D.hpp.

{
    return lhs.v_[0] == rhs.v_[0] && lhs.v_[1] == rhs.v_[1] && lhs.v_[2] == rhs.v_[2] && lhs.v_[3] == rhs.v_[3] &&
           lhs.v_[4] == rhs.v_[4] && lhs.v_[5] == rhs.v_[5] && lhs.v_[6] == rhs.v_[6] && lhs.v_[7] == rhs.v_[7] &&
           lhs.v_[8] == rhs.v_[8];
}
Matrix3D plus_transpose ( const Matrix3D A,
const Matrix3D B 
) [friend]

\( + B^T \)

Definition at line 627 of file Matrix3D.hpp.

{
    Matrix3D tmp( A );
    tmp.plus_transpose_equal( B );
    return tmp;
}
void plusEqaA ( Matrix3D B,
const double  a,
const Matrix3D A 
) [friend]

\( B += a*A \)

Definition at line 792 of file Matrix3D.hpp.

{
    B.v_[0] += a * A.v_[0];
    B.v_[1] += a * A.v_[1];
    B.v_[2] += a * A.v_[2];
    B.v_[3] += a * A.v_[3];
    B.v_[4] += a * A.v_[4];
    B.v_[5] += a * A.v_[5];
    B.v_[6] += a * A.v_[6];
    B.v_[7] += a * A.v_[7];
    B.v_[8] += a * A.v_[8];
}
void plusEqAx ( Vector3D v,
const Matrix3D A,
const Vector3D x 
) [friend]

\( v += A*x \)

Definition at line 771 of file Matrix3D.hpp.

{
    v.mCoords[0] += A.v_[0] * x[0] + A.v_[1] * x.mCoords[1] + A.v_[2] * x.mCoords[2];
    v.mCoords[1] += A.v_[3] * x[0] + A.v_[4] * x.mCoords[1] + A.v_[5] * x.mCoords[2];
    v.mCoords[2] += A.v_[6] * x[0] + A.v_[7] * x.mCoords[1] + A.v_[8] * x.mCoords[2];
}
void plusEqTransAx ( Vector3D v,
const Matrix3D A,
const Vector3D x 
) [friend]

\( v += A^T*x \)

Definition at line 785 of file Matrix3D.hpp.

{
    v.mCoords[0] += A.v_[0] * x.mCoords[0] + A.v_[3] * x.mCoords[1] + A.v_[6] * x.mCoords[2];
    v.mCoords[1] += A.v_[1] * x.mCoords[0] + A.v_[4] * x.mCoords[1] + A.v_[7] * x.mCoords[2];
    v.mCoords[2] += A.v_[2] * x.mCoords[0] + A.v_[5] * x.mCoords[1] + A.v_[8] * x.mCoords[2];
}
void QR ( Matrix3D Q,
Matrix3D R,
const Matrix3D A 
) [friend]

\( Q*R = A \)

Definition at line 836 of file Matrix3D.hpp.

{
    // Compute the QR factorization of A.  This code uses the
    // Modified Gram-Schmidt method for computing the factorization.
    // The Householder version is more stable, but costs twice as many
    // floating point operations.

    Q = A;

    R[0][0]         = sqrt( Q[0][0] * Q[0][0] + Q[1][0] * Q[1][0] + Q[2][0] * Q[2][0] );
    double temp_dbl = 1.0 / R[0][0];
    R[1][0]         = 0.0L;
    R[2][0]         = 0.0L;
    // Q[0][0] /= R[0][0];
    // Q[1][0] /= R[0][0];
    // Q[2][0] /= R[0][0];
    Q[0][0] *= temp_dbl;
    Q[1][0] *= temp_dbl;
    Q[2][0] *= temp_dbl;

    R[0][1] = Q[0][0] * Q[0][1] + Q[1][0] * Q[1][1] + Q[2][0] * Q[2][1];
    Q[0][1] -= Q[0][0] * R[0][1];
    Q[1][1] -= Q[1][0] * R[0][1];
    Q[2][1] -= Q[2][0] * R[0][1];

    R[0][2] = Q[0][0] * Q[0][2] + Q[1][0] * Q[1][2] + Q[2][0] * Q[2][2];
    Q[0][2] -= Q[0][0] * R[0][2];
    Q[1][2] -= Q[1][0] * R[0][2];
    Q[2][2] -= Q[2][0] * R[0][2];

    R[1][1]  = sqrt( Q[0][1] * Q[0][1] + Q[1][1] * Q[1][1] + Q[2][1] * Q[2][1] );
    temp_dbl = 1.0 / R[1][1];
    R[2][1]  = 0.0L;
    //     Q[0][1] /= R[1][1];
    //     Q[1][1] /= R[1][1];
    //     Q[2][1] /= R[1][1];
    Q[0][1] *= temp_dbl;
    Q[1][1] *= temp_dbl;
    Q[2][1] *= temp_dbl;

    R[1][2] = Q[0][1] * Q[0][2] + Q[1][1] * Q[1][2] + Q[2][1] * Q[2][2];
    Q[0][2] -= Q[0][1] * R[1][2];
    Q[1][2] -= Q[1][1] * R[1][2];
    Q[2][2] -= Q[2][1] * R[1][2];

    R[2][2]  = sqrt( Q[0][2] * Q[0][2] + Q[1][2] * Q[1][2] + Q[2][2] * Q[2][2] );
    temp_dbl = 1.0 / R[2][2];

    //     Q[0][2] /= R[2][2];
    //     Q[1][2] /= R[2][2];
    //     Q[2][2] /= R[2][2];
    Q[0][2] *= temp_dbl;
    Q[1][2] *= temp_dbl;
    Q[2][2] *= temp_dbl;

    return;
}
void timesInvA ( Matrix3D B,
const Matrix3D A 
) [friend]

\( B *= A^{-1} \)

Definition at line 828 of file Matrix3D.hpp.

{

    Matrix3D Ainv;
    inv( Ainv, A );
    B = B * Ainv;
}
Matrix3D transpose ( const Matrix3D A) [friend]

Definition at line 506 of file Matrix3D.hpp.

{
    Matrix3D S;
    //     size_t i;
    //     for (i=0; i<3; ++i) {
    //         S[size_t(0)][i] = A[i][0];
    //         S[size_t(1)][i] = A[i][1];
    //         S[size_t(2)][i] = A[i][2];
    //     }
    S.v_[0] = A.v_[0];
    S.v_[1] = A.v_[3];
    S.v_[2] = A.v_[6];
    S.v_[3] = A.v_[1];
    S.v_[4] = A.v_[4];
    S.v_[5] = A.v_[7];
    S.v_[6] = A.v_[2];
    S.v_[7] = A.v_[5];
    S.v_[8] = A.v_[8];

    return S;
}

Member Data Documentation

List of all members.


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines