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

Vector3D is the object that effeciently stores information about about three-deminsional vectors. It is also the parent class of MsqVertex. More...

#include <Vector3D.hpp>

+ Inheritance diagram for MBMesquite::Vector3D:

Public Member Functions

 Vector3D ()
 Vector3D (double xyz)
 Vector3D (double x, double y, double z)
 Vector3D (const double xyz[3])
 Vector3D (const Vector3D &to_copy)
double x () const
double y () const
double z () const
void get_coordinates (double &x, double &y, double &z) const
void get_coordinates (double xyz[3]) const
const double & operator[] (size_t index) const
void x (const double x)
void y (const double y)
void z (const double z)
void set (const double x, const double y, const double z)
void set (const double xyz[3])
void set (const Vector3D &to_copy)
double & operator[] (size_t index)
Vector3Doperator= (const Vector3D &to_copy)
Vector3Doperator= (const double &to_copy)
Vector3D operator- () const
Vector3Doperator*= (const double scalar)
Vector3Doperator/= (const double scalar)
 divides each Vector3D entry by the given scalar.
Vector3Doperator*= (const Vector3D &rhs)
Vector3Doperator+= (const Vector3D &rhs)
Vector3Doperator-= (const Vector3D &rhs)
int within_tolerance_box (const Vector3D &compare_to, double tolerance) const
double length_squared () const
double length () const
void set_length (const double new_length)
void normalize ()
Vector3D operator~ () const
const double * to_array () const
double * to_array ()

Static Public Member Functions

static double distance_between (const Vector3D &p1, const Vector3D &p2)
static double interior_angle (const Vector3D &a, const Vector3D &b, MsqError &err)
static Vector3D interpolate (const double param, const Vector3D &p1, const Vector3D &p2)

Protected Attributes

double mCoords [3]

Friends

const Vector3D operator+ (const Vector3D &lhs, const Vector3D &rhs)
const Vector3D operator- (const Vector3D &lhs, const Vector3D &rhs)
const Vector3D operator* (const Vector3D &lhs, const double scalar)
 lhs * scalar
const Vector3D operator* (const double scalar, const Vector3D &rhs)
 scalar * rhs
const Vector3D operator/ (const Vector3D &lhs, const double scalar)
double operator% (const Vector3D &v1, const Vector3D &v2)
 dot product
double inner (const Vector3D v1[], const Vector3D v2[], int n)
 dot product for array
double operator% (const double scalar, const Vector3D &v2)
 scalar * sum_i v2[i]
double operator% (const Vector3D &v1, const double scalar)
 scalar * sum_i v1[i]
const Vector3D operator* (const Vector3D &v1, const Vector3D &v2)
 cross product
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 plusEqTransAx (Vector3D &v, const Matrix3D &A, const Vector3D &x)
 \( v += A^T*x \)
void eqTransAx (Vector3D &v, const Matrix3D &A, const Vector3D &x)
bool operator== (const Vector3D &lhs, const Vector3D &rhs)
bool operator!= (const Vector3D &lhs, const Vector3D &rhs)
double length (const Vector3D *v, int n)
 L2 norm for an array of Vector3Ds.
double Linf (const Vector3D *v, int n)
 L inf norm for array of Vector3Ds.

Detailed Description

Vector3D is the object that effeciently stores information about about three-deminsional vectors. It is also the parent class of MsqVertex.

Definition at line 55 of file Vector3D.hpp.


Constructor & Destructor Documentation

Definition at line 168 of file Vector3D.hpp.

References mCoords.

Referenced by operator-().

{
    mCoords[0] = 0;
    mCoords[1] = 0;
    mCoords[2] = 0;
}
MBMesquite::Vector3D::Vector3D ( double  xyz) [inline]

Definition at line 174 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] = p_x;
    mCoords[1] = p_x;
    mCoords[2] = p_x;
}
MBMesquite::Vector3D::Vector3D ( double  x,
double  y,
double  z 
) [inline]

Definition at line 180 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] = p_x;
    mCoords[1] = p_y;
    mCoords[2] = p_z;
}
MBMesquite::Vector3D::Vector3D ( const double  xyz[3]) [inline]

Definition at line 186 of file Vector3D.hpp.

References mCoords.

{
    std::memcpy( mCoords, xyz, 3 * sizeof( double ) );
}
MBMesquite::Vector3D::Vector3D ( const Vector3D to_copy) [inline]

Definition at line 190 of file Vector3D.hpp.

References mCoords.

{
    std::memcpy( mCoords, to_copy.mCoords, 3 * sizeof( double ) );
}

Member Function Documentation

double MBMesquite::Vector3D::distance_between ( const Vector3D p1,
const Vector3D p2 
) [inline, static]

Definition at line 384 of file Vector3D.hpp.

References length().

Referenced by MBMesquite::MeshDomainAssoc::MeshDomainAssoc(), and MBMesquite::MeshUtil::meshes_are_different().

{
    Vector3D v = p2 - p1;
    return v.length();
}
void MBMesquite::Vector3D::get_coordinates ( double  xyz[3]) const [inline]

Definition at line 214 of file Vector3D.hpp.

References mCoords.

{
    std::memcpy( xyz, mCoords, 3 * sizeof( double ) );
}
double MBMesquite::Vector3D::interior_angle ( const Vector3D a,
const Vector3D b,
MsqError err 
) [static]

Definition at line 47 of file Vector3D.cpp.

References INTERNAL_ERROR, moab::Util::is_finite(), length(), and MSQ_SETERR.

Referenced by MBMesquite::MeshWriter::Projection::init().

{
    double len1      = lhs.length();
    double len2      = rhs.length();
    double angle_cos = ( lhs % rhs ) / ( len1 * len2 );
    if( !moab::Util::is_finite( angle_cos ) )
    {
        MSQ_SETERR( err )( MsqError::INTERNAL_ERROR );
        return 0.0;
    }

    // Adjust the cosine if slightly out of range
    if( ( angle_cos > 1.0 ) && ( angle_cos < 1.0001 ) )
    {
        angle_cos = 1.0;
    }
    else if( angle_cos < -1.0 && angle_cos > -1.0001 )
    {
        angle_cos = -1.0;
    }

    return std::acos( angle_cos );
}
Vector3D MBMesquite::Vector3D::interpolate ( const double  param,
const Vector3D p1,
const Vector3D p2 
) [inline, static]

Definition at line 482 of file Vector3D.hpp.

{
    return ( 1 - param ) * p1 + param * p2;
}
double MBMesquite::Vector3D::length ( ) const [inline]

Definition at line 401 of file Vector3D.hpp.

References mCoords.

Referenced by MBMesquite::MeshTransform::add_rotation(), SphericalDomainTest::check_closest_pt(), MBMesquite::MsqMeshEntity::check_element_orientation_corners(), SphericalDomainTest::check_normal(), check_results(), check_slaved_coords(), MBMesquite::MsqCircle::closest(), MBMesquite::MsqSphere::closest(), compare_node_coords(), compute_corner_area(), MBMesquite::MsqMeshEntity::compute_corner_normals(), MBMesquite::MsqMeshEntity::compute_signed_area(), distance_between(), distance_from_origin(), MBMesquite::AspectRatioGammaQualityMetric::evaluate(), MBMesquite::IdealWeightMeanRatio::evaluate(), MBMesquite::BoundedCylinderDomain::evaluate(), MBMesquite::IdealWeightInverseMeanRatio::evaluate(), MBMesquite::CylinderDomain::evaluate(), MBMesquite::ConicDomain::evaluate(), MBMesquite::EdgeLengthRangeQualityMetric::evaluate_common(), MBMesquite::IdealWeightMeanRatio::evaluate_with_gradient(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_gradient(), MBMesquite::IdealWeightMeanRatio::evaluate_with_Hessian(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_Hessian(), MBMesquite::IdealWeightMeanRatio::evaluate_with_Hessian_diagonal(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_Hessian_diagonal(), MBMesquite::NonSmoothDescent::find_plane_normal(), MBMesquite::MeshWriter::Projection::init(), interior_angle(), MBMesquite::MsqSphere::intersect(), main(), MBMesquite::MeshUtil::meshes_are_different(), MBMesquite::MsqPlane::MsqPlane(), MBMesquite::DomainUtil::non_coplanar_vertices(), MBMesquite::CircleDomain::position_from_length(), MBMesquite::MsqCircle::radial_vector(), set_length(), MBMesquite::SphericalDomain::snap_to(), MBMesquite::LineDomainTest::test_arc_length(), MBMesquite::CircleDomainTest::test_arc_length(), GeomPrimTest::test_circle_closest_to_point(), GeomPrimTest::test_circle_from_two_points(), ConicDomainTest::test_construct(), GeomPrimTest::test_line_basic(), ConicDomainTest::test_normal_at(), GeomPrimTest::test_plane_basic(), MBMesquite::CircleDomainTest::test_position_from_length(), MBMesquite::LineDomainTest::test_position_from_length(), VtkTest::test_read_quadratic(), MBMesquite::CircleDomainTest::test_snap_to(), CylinderDomainTest::test_x_normal_at(), CylinderDomainTest::test_z_normal_at(), untangle_function_2d(), and MBMesquite::SphericalDomain::vertex_normal_at().

{
    return std::sqrt( mCoords[0] * mCoords[0] + mCoords[1] * mCoords[1] + mCoords[2] * mCoords[2] );
}
Vector3D & MBMesquite::Vector3D::operator*= ( const double  scalar) [inline]

Definition at line 277 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] *= scalar;
    mCoords[1] *= scalar;
    mCoords[2] *= scalar;
    return *this;
}
Vector3D & MBMesquite::Vector3D::operator*= ( const Vector3D rhs) [inline]

Definition at line 292 of file Vector3D.hpp.

References mCoords.

{
    double new_coords[3] = { mCoords[1] * rhs.mCoords[2] - mCoords[2] * rhs.mCoords[1],
                             mCoords[2] * rhs.mCoords[0] - mCoords[0] * rhs.mCoords[2],
                             mCoords[0] * rhs.mCoords[1] - mCoords[1] * rhs.mCoords[0] };
    std::memcpy( mCoords, new_coords, 3 * sizeof( double ) );
    return *this;
}
Vector3D & MBMesquite::Vector3D::operator+= ( const Vector3D rhs) [inline]

Definition at line 300 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] += rhs.mCoords[0];
    mCoords[1] += rhs.mCoords[1];
    mCoords[2] += rhs.mCoords[2];
    return *this;
}
Vector3D MBMesquite::Vector3D::operator- ( ) const [inline]

Definition at line 273 of file Vector3D.hpp.

References mCoords, and Vector3D().

{
    return Vector3D( -mCoords[0], -mCoords[1], -mCoords[2] );
}
Vector3D & MBMesquite::Vector3D::operator-= ( const Vector3D rhs) [inline]

Definition at line 307 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] -= rhs.mCoords[0];
    mCoords[1] -= rhs.mCoords[1];
    mCoords[2] -= rhs.mCoords[2];
    return *this;
}
Vector3D & MBMesquite::Vector3D::operator/= ( const double  scalar) [inline]

divides each Vector3D entry by the given scalar.

Definition at line 285 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] /= scalar;
    mCoords[1] /= scalar;
    mCoords[2] /= scalar;
    return *this;
}
Vector3D & MBMesquite::Vector3D::operator= ( const Vector3D to_copy) [inline]

Reimplemented in MBMesquite::MsqVertex.

Definition at line 255 of file Vector3D.hpp.

References mCoords.

Referenced by MBMesquite::MsqVertex::operator=().

{
    mCoords[0] = to_copy.mCoords[0];
    mCoords[1] = to_copy.mCoords[1];
    mCoords[2] = to_copy.mCoords[2];
    //    memcpy(mCoords, to_copy.mCoords, 3*sizeof(double));
    return *this;
}
Vector3D & MBMesquite::Vector3D::operator= ( const double &  to_copy) [inline]

Definition at line 264 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] = to_copy;
    mCoords[1] = to_copy;
    mCoords[2] = to_copy;
    return *this;
}
const double & MBMesquite::Vector3D::operator[] ( size_t  index) const [inline]

Definition at line 218 of file Vector3D.hpp.

References mCoords.

{
    return mCoords[index];
}
double & MBMesquite::Vector3D::operator[] ( size_t  index) [inline]

Definition at line 250 of file Vector3D.hpp.

References mCoords.

{
    return mCoords[index];
}

Definition at line 142 of file Vector3D.hpp.

References MBMesquite::length().

    {
        return *this * ( 1.0 / length() );
    }
void MBMesquite::Vector3D::set ( const double  x,
const double  y,
const double  z 
) [inline]

Definition at line 236 of file Vector3D.hpp.

References mCoords.

Referenced by MBMesquite::average_corner_diagonals(), HigherOrderTest::basic_hex_test(), HigherOrderTest::basic_quad_test(), HigherOrderTest::basic_tet_test(), MBMesquite::PointDomain::closest_point(), MBMesquite::SphericalDomain::closest_point(), MBMesquite::ParallelHelperImpl::comm_smoothed_vtx_b(), MBMesquite::ParallelHelperImpl::comm_smoothed_vtx_b_no_all(), MBMesquite::ParallelHelperImpl::comm_smoothed_vtx_nb(), MBMesquite::ParallelHelperImpl::comm_smoothed_vtx_nb_no_all(), MBMesquite::ParallelHelperImpl::comm_smoothed_vtx_tnb(), MBMesquite::ParallelHelperImpl::comm_smoothed_vtx_tnb_no_all(), MBMesquite::PointDomain::element_normal_at(), MBMesquite::XYRectangle::element_normal_at(), MBMesquite::CylinderDomain::evaluate(), get_ideal_quad(), get_ideal_tet(), main(), MBMesquite::MsqCommonIGeom::move_to(), MBMesquite::MsqCommonIGeom::normal(), MBMesquite::SphericalDomain::snap_to(), MsqHessianTest::test_axpy(), MsqHessianTest::test_cg_solver(), MsqVertexTest::test_hex_vertices(), MsqMeshEntityTest::test_hex_vertices(), PatchDataTest::test_move_free_vertices_constrained(), PatchDataTest::test_movement_function(), HigherOrderTest::test_quad_basic_mid_spin(), Matrix3DTest::test_times_vector(), Matrix3DTest::test_vector_times(), MeshInterfaceTest::test_vertices(), CylinderDomainTest::test_x_normal_at(), BoundedCylinderDomainTest::test_x_snap_to(), CylinderDomainTest::test_x_snap_to(), CylinderDomainTest::test_z_normal_at(), CylinderDomainTest::test_z_snap_to(), BoundedCylinderDomainTest::test_z_snap_to(), SphereDomainArg::value(), CylinderDomainArg::value(), PlanarDomainArg::value(), LineDomainArg::value(), CircleDomainArg::value(), MBMesquite::PointDomain::vertex_normal_at(), MBMesquite::SphericalDomain::vertex_normal_at(), MBMesquite::XYRectangle::vertex_normal_at(), MBMesquite::TagVertexMesh::vertices_get_coordinates(), MBMesquite::MsqIMesh::vertices_get_coordinates(), MBMesquite::MsqMOAB::vertices_get_coordinates(), and MBMesquite::ArrayMesh::vertices_get_coordinates().

{
    mCoords[0] = p_x;
    mCoords[1] = p_y;
    mCoords[2] = p_z;
}
void MBMesquite::Vector3D::set ( const double  xyz[3]) [inline]

Definition at line 242 of file Vector3D.hpp.

References mCoords.

{
    std::memcpy( mCoords, xyz, 3 * sizeof( double ) );
}
void MBMesquite::Vector3D::set ( const Vector3D to_copy) [inline]

Definition at line 246 of file Vector3D.hpp.

References mCoords.

{
    std::memcpy( mCoords, to_copy.mCoords, 3 * sizeof( double ) );
}
void MBMesquite::Vector3D::set_length ( const double  new_length) [inline]

Definition at line 471 of file Vector3D.hpp.

References length().

Referenced by normalize().

{
    double factor = new_length / length();
    *this *= factor;
}
double* MBMesquite::Vector3D::to_array ( ) [inline]

Definition at line 158 of file Vector3D.hpp.

    {
        return mCoords;
    }
int MBMesquite::Vector3D::within_tolerance_box ( const Vector3D compare_to,
double  tolerance 
) const [inline]

Definition at line 389 of file Vector3D.hpp.

References mCoords.

Referenced by VtkTest::check_4quad_block(), VtkTest::check_8hex_block(), and smooth_mesh().

{
    return ( ( std::fabs( this->mCoords[0] - compare_to.mCoords[0] ) < tolerance ) &&
             ( std::fabs( this->mCoords[1] - compare_to.mCoords[1] ) < tolerance ) &&
             ( std::fabs( this->mCoords[2] - compare_to.mCoords[2] ) < tolerance ) );
}
void MBMesquite::Vector3D::x ( const double  x) [inline]

Definition at line 224 of file Vector3D.hpp.

References mCoords.

{
    mCoords[0] = p_x;
}
void MBMesquite::Vector3D::y ( const double  y) [inline]

Definition at line 228 of file Vector3D.hpp.

References mCoords.

{
    mCoords[1] = p_y;
}
void MBMesquite::Vector3D::z ( const double  z) [inline]

Definition at line 232 of file Vector3D.hpp.

References mCoords.

{
    mCoords[2] = p_z;
}

Friends And Related Function Documentation

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 inner ( const Vector3D  v1[],
const Vector3D  v2[],
int  n 
) [friend]

dot product for array

Dot product for arrays of Vector3Ds. see also operator% .

Definition at line 344 of file Vector3D.hpp.

{
    int i;
    double dot = 0;
    for( i = 0; i < n; ++i )
        dot += lhs[i].mCoords[0] * rhs[i].mCoords[0] + lhs[i].mCoords[1] * rhs[i].mCoords[1] +
               lhs[i].mCoords[2] * rhs[i].mCoords[2];
    return dot;
}
double length ( const Vector3D v,
int  n 
) [friend]

L2 norm for an array of Vector3Ds.

Definition at line 434 of file Vector3D.hpp.

{
    return std::sqrt( length_squared( v, n ) );
}
double Linf ( const Vector3D v,
int  n 
) [friend]

L inf norm for array of Vector3Ds.

Definition at line 443 of file Vector3D.hpp.

{
    double max = 0;
    // loop over the length of the array
    for( int i = 0; i < n; ++i )
    {
        if( max < std::fabs( v[i][0] ) ) max = std::fabs( v[i][0] );
        if( max < std::fabs( v[i][1] ) ) max = std::fabs( v[i][1] );
        if( max < std::fabs( v[i][2] ) ) max = std::fabs( v[i][2] );
    }
    // return the value of the largest entry in the array
    return max;
}
bool operator!= ( const Vector3D lhs,
const Vector3D rhs 
) [friend]

Definition at line 492 of file Vector3D.hpp.

{
    return v1.mCoords[0] != v2.mCoords[0] || v1.mCoords[1] != v2.mCoords[1] || v1.mCoords[2] != v2.mCoords[2];
}
double operator% ( const Vector3D v1,
const Vector3D v2 
) [friend]

dot product

Definition at line 337 of file Vector3D.hpp.

{
    return ( lhs.mCoords[0] * rhs.mCoords[0] + lhs.mCoords[1] * rhs.mCoords[1] + lhs.mCoords[2] * rhs.mCoords[2] );
}
double operator% ( const double  scalar,
const Vector3D v2 
) [friend]

scalar * sum_i v2[i]

Definition at line 363 of file Vector3D.hpp.

{
    return ( scalar * ( rhs.mCoords[0] + rhs.mCoords[1] + rhs.mCoords[2] ) );
}
double operator% ( const Vector3D v1,
const double  scalar 
) [friend]

scalar * sum_i v1[i]

Definition at line 368 of file Vector3D.hpp.

{
    return ( scalar * ( lhs.mCoords[0] + lhs.mCoords[1] + lhs.mCoords[2] ) );
}
const Vector3D operator* ( const Vector3D lhs,
const double  scalar 
) [friend]

lhs * scalar

Definition at line 324 of file Vector3D.hpp.

{
    return Vector3D( lhs.x() * scalar, lhs.y() * scalar, lhs.z() * scalar );
}
const Vector3D operator* ( const double  scalar,
const Vector3D rhs 
) [friend]

scalar * rhs

Definition at line 328 of file Vector3D.hpp.

{
    return Vector3D( rhs.x() * scalar, rhs.y() * scalar, rhs.z() * scalar );
}
const Vector3D operator* ( const Vector3D v1,
const Vector3D v2 
) [friend]

cross product

Definition at line 373 of file Vector3D.hpp.

{
    return Vector3D( lhs.mCoords[1] * rhs.mCoords[2] - lhs.mCoords[2] * rhs.mCoords[1],
                     lhs.mCoords[2] * rhs.mCoords[0] - lhs.mCoords[0] * rhs.mCoords[2],
                     lhs.mCoords[0] * rhs.mCoords[1] - lhs.mCoords[1] * rhs.mCoords[0] );
}
const Vector3D operator+ ( const Vector3D lhs,
const Vector3D rhs 
) [friend]

Definition at line 316 of file Vector3D.hpp.

{
    return Vector3D( lhs.x() + rhs.x(), lhs.y() + rhs.y(), lhs.z() + rhs.z() );
}
const Vector3D operator- ( const Vector3D lhs,
const Vector3D rhs 
) [friend]

Definition at line 320 of file Vector3D.hpp.

{
    return Vector3D( lhs.x() - rhs.x(), lhs.y() - rhs.y(), lhs.z() - rhs.z() );
}
const Vector3D operator/ ( const Vector3D lhs,
const double  scalar 
) [friend]

Definition at line 332 of file Vector3D.hpp.

{
    assert( scalar != 0 );
    return Vector3D( lhs.x() / scalar, lhs.y() / scalar, lhs.z() / scalar );
}
bool operator== ( const Vector3D lhs,
const Vector3D rhs 
) [friend]

Definition at line 487 of file Vector3D.hpp.

{
    return v1.mCoords[0] == v2.mCoords[0] && v1.mCoords[1] == v2.mCoords[1] && v1.mCoords[2] == v2.mCoords[2];
}
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];
}

Member Data Documentation

List of all members.


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