Mesh Oriented datABase  (version 5.4.1)
Array-based unstructured mesh datastructure
VerdictVector Class Reference

#include <VerdictVector.hpp>

Public Member Functions

 VerdictVector ()
 VerdictVector (const double x, const double y, const double z)
 VerdictVector (const double xyz[3])
 VerdictVector (const VerdictVector &tail, const VerdictVector &head)
 VerdictVector (const VerdictVector &copy_from)
void set (const double xv, const double yv, const double zv)
void set (const double xyz[3])
void set (const VerdictVector &tail, const VerdictVector &head)
void set (const VerdictVector &to_copy)
double x () const
double y () const
double z () const
void get_xyz (double &x, double &y, double &z)
void get_xyz (double xyz[3])
double & r ()
double & theta ()
void x (const double xv)
void y (const double yv)
void z (const double zv)
void r (const double xv)
void theta (const double yv)
void xy_to_rtheta ()
void rtheta_to_xy ()
void scale_angle (double gamma, double)
void blow_out (double gamma, double gamma2=0.0)
void rotate (double angle, double)
void reflect_about_xaxis (double dummy, double)
double normalize ()
VerdictVectorlength (const double new_length)
double length () const
double distance_between (const VerdictVector &test_vector)
double length_squared () const
double interior_angle (const VerdictVector &otherVector)
double vector_angle_quick (const VerdictVector &vec1, const VerdictVector &vec2)
double vector_angle (const VerdictVector &vector1, const VerdictVector &vector2) const
void perpendicular_z ()
void print_me ()
void orthogonal_vectors (VerdictVector &vector2, VerdictVector &vector3)
void next_point (const VerdictVector &direction, double distance, VerdictVector &out_point)
bool within_tolerance (const VerdictVector &vectorPtr2, double tolerance) const
VerdictVectoroperator+= (const VerdictVector &vec)
VerdictVectoroperator-= (const VerdictVector &vec)
VerdictVectoroperator*= (const VerdictVector &vec)
VerdictVectoroperator*= (const double scalar)
VerdictVectoroperator/= (const double scalar)
VerdictVector operator- () const
VerdictVectoroperator= (const VerdictVector &from)

Private Attributes

double xVal
double yVal
double zVal

Friends

VerdictVector operator~ (const VerdictVector &vec)
VerdictVector operator+ (const VerdictVector &v1, const VerdictVector &v2)
VerdictVector operator- (const VerdictVector &v1, const VerdictVector &v2)
VerdictVector operator* (const VerdictVector &v1, const VerdictVector &v2)
VerdictVector operator* (const VerdictVector &v1, const double sclr)
VerdictVector operator* (const double sclr, const VerdictVector &v1)
double operator% (const VerdictVector &v1, const VerdictVector &v2)
VerdictVector operator/ (const VerdictVector &v1, const double sclr)
int operator== (const VerdictVector &v1, const VerdictVector &v2)
int operator!= (const VerdictVector &v1, const VerdictVector &v2)
VerdictVector interpolate (const double param, const VerdictVector &v1, const VerdictVector &v2)

Detailed Description

Definition at line 35 of file VerdictVector.hpp.


Constructor & Destructor Documentation

Definition at line 325 of file VerdictVector.hpp.

Referenced by operator-().

: xVal( 0 ), yVal( 0 ), zVal( 0 ) {}
VerdictVector::VerdictVector ( const double  x,
const double  y,
const double  z 
) [inline]

Definition at line 332 of file VerdictVector.hpp.

    : xVal( xIn ), yVal( yIn ), zVal( zIn )
{
}
VerdictVector::VerdictVector ( const double  xyz[3])

Definition at line 438 of file VerdictVector.cpp.

: xVal( xyz[0] ), yVal( xyz[1] ), zVal( xyz[2] ) {}
VerdictVector::VerdictVector ( const VerdictVector tail,
const VerdictVector head 
) [inline]

Definition at line 327 of file VerdictVector.hpp.

    : xVal( head.xVal - tail.xVal ), yVal( head.yVal - tail.yVal ), zVal( head.zVal - tail.zVal )
{
}
VerdictVector::VerdictVector ( const VerdictVector copy_from) [inline]

Definition at line 320 of file VerdictVector.hpp.

    : xVal( copy_from.xVal ), yVal( copy_from.yVal ), zVal( copy_from.zVal )
{
}

Member Function Documentation

void VerdictVector::blow_out ( double  gamma,
double  gamma2 = 0.0 
)

Definition at line 133 of file VerdictVector.cpp.

References r(), rtheta_to_xy(), and xy_to_rtheta().

{
    // if gamma == 1, then
    // map on a circle : r'^2 = sqrt( 1 - (1-r)^2 )
    // if gamma ==0, then map back to itself
    // in between, linearly interpolate
    xy_to_rtheta();
    //  r() = sqrt( (2. - r()) * r() ) * gamma  + r() * (1-gamma);
    assert( gamma > 0.0 );
    // the following limits should really be roundoff-based
    if( r() > rmin * 1.001 && r() < 1.001 )
    {
        r() = rmin + pow( r(), gamma ) * ( 1.0 - rmin );
    }
    rtheta_to_xy();
}
double VerdictVector::distance_between ( const VerdictVector test_vector)

Definition at line 45 of file VerdictVector.cpp.

References x(), xVal, y(), yVal, z(), and zVal.

{
    double xv = xVal - test_vector.x();
    double yv = yVal - test_vector.y();
    double zv = zVal - test_vector.z();

    return ( sqrt( xv * xv + yv * yv + zv * zv ) );
}
void VerdictVector::get_xyz ( double &  x,
double &  y,
double &  z 
) [inline]

Definition at line 258 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

Referenced by orthogonal_vectors(), and v_tet_aspect_frobenius().

{
    xv = xVal;
    yv = yVal;
    zv = zVal;
}
void VerdictVector::get_xyz ( double  xyz[3]) [inline]

Definition at line 252 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

{
    xyz[0] = xVal;
    xyz[1] = yVal;
    xyz[2] = zVal;
}
double VerdictVector::interior_angle ( const VerdictVector otherVector)

Definition at line 64 of file VerdictVector.cpp.

References length(), and VERDICT_PI.

Referenced by v_tri_maximum_angle(), and v_tri_minimum_angle().

{
    double cosAngle = 0., angleRad = 0., len1, len2 = 0.;

    if( ( ( len1 = this->length() ) > 0 ) && ( ( len2 = otherVector.length() ) > 0 ) )
        cosAngle = ( *this % otherVector ) / ( len1 * len2 );
    else
    {
        assert( len1 > 0 );
        assert( len2 > 0 );
    }

    if( ( cosAngle > 1.0 ) && ( cosAngle < 1.0001 ) )
    {
        cosAngle = 1.0;
        angleRad = acos( cosAngle );
    }
    else if( cosAngle < -1.0 && cosAngle > -1.0001 )
    {
        cosAngle = -1.0;
        angleRad = acos( cosAngle );
    }
    else if( cosAngle >= -1.0 && cosAngle <= 1.0 )
        angleRad = acos( cosAngle );
    else
    {
        assert( cosAngle < 1.0001 && cosAngle > -1.0001 );
    }

    return ( ( angleRad * 180. ) / VERDICT_PI );
}
double VerdictVector::length ( ) const [inline]

Definition at line 477 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

Referenced by interior_angle(), length(), normalize(), and xy_to_rtheta().

{
    return ( sqrt( xVal * xVal + yVal * yVal + zVal * zVal ) );
}
void VerdictVector::next_point ( const VerdictVector direction,
double  distance,
VerdictVector out_point 
)

Definition at line 425 of file VerdictVector.cpp.

References normalize(), x(), xVal, y(), yVal, z(), and zVal.

{
    VerdictVector my_direction = direction;
    my_direction.normalize();

    // Determine next point in space
    out_point.x( xVal + ( distance * my_direction.x() ) );
    out_point.y( yVal + ( distance * my_direction.y() ) );
    out_point.z( zVal + ( distance * my_direction.z() ) );

    return;
}
double VerdictVector::normalize ( ) [inline]
VerdictVector & VerdictVector::operator*= ( const VerdictVector vec) [inline]

Definition at line 308 of file VerdictVector.hpp.

References x(), xVal, y(), yVal, z(), and zVal.

{
    double xcross, ycross, zcross;
    xcross = yVal * vector.z() - zVal * vector.y();
    ycross = zVal * vector.x() - xVal * vector.z();
    zcross = xVal * vector.y() - yVal * vector.x();
    xVal   = xcross;
    yVal   = ycross;
    zVal   = zcross;
    return *this;
}
VerdictVector & VerdictVector::operator*= ( const double  scalar) [inline]

Definition at line 382 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

{
    xVal *= scalar;
    yVal *= scalar;
    zVal *= scalar;
    return *this;
}
VerdictVector & VerdictVector::operator+= ( const VerdictVector vec) [inline]

Definition at line 292 of file VerdictVector.hpp.

References x(), xVal, y(), yVal, z(), and zVal.

{
    xVal += vector.x();
    yVal += vector.y();
    zVal += vector.z();
    return *this;
}
VerdictVector VerdictVector::operator- ( ) const [inline]

Definition at line 414 of file VerdictVector.hpp.

References VerdictVector(), xVal, yVal, and zVal.

{
    return VerdictVector( -xVal, -yVal, -zVal );
}
VerdictVector & VerdictVector::operator-= ( const VerdictVector vec) [inline]

Definition at line 300 of file VerdictVector.hpp.

References x(), xVal, y(), yVal, z(), and zVal.

{
    xVal -= vector.x();
    yVal -= vector.y();
    zVal -= vector.z();
    return *this;
}
VerdictVector & VerdictVector::operator/= ( const double  scalar) [inline]

Definition at line 391 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

{
    assert( scalar != 0 );
    xVal /= scalar;
    yVal /= scalar;
    zVal /= scalar;
    return *this;
}
VerdictVector & VerdictVector::operator= ( const VerdictVector from) [inline]

Definition at line 368 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

{
    xVal = from.xVal;
    yVal = from.yVal;
    zVal = from.zVal;
    return *this;
}
void VerdictVector::orthogonal_vectors ( VerdictVector vector2,
VerdictVector vector3 
)

Definition at line 354 of file VerdictVector.cpp.

References get_xyz(), normalize(), and set().

{
    double xv[3];
    unsigned short i    = 0;
    unsigned short imin = 0;
    double rmin         = 1.0E20;
    unsigned short iperm1[3];
    unsigned short iperm2[3];
    unsigned short cont_flag = 1;
    double vec1[3], vec2[3];
    double rmag;

    // Copy the input vector and normalize it
    VerdictVector vector1 = *this;
    vector1.normalize();

    // Initialize perm flags
    iperm1[0] = 1;
    iperm1[1] = 2;
    iperm1[2] = 0;
    iperm2[0] = 2;
    iperm2[1] = 0;
    iperm2[2] = 1;

    // Get into the array format we can work with
    vector1.get_xyz( vec1 );

    while( i < 3 && cont_flag )
    {
        if( fabs( vec1[i] ) < 1e-6 )
        {
            vec2[i]         = 1.0;
            vec2[iperm1[i]] = 0.0;
            vec2[iperm2[i]] = 0.0;
            cont_flag       = 0;
        }

        if( fabs( vec1[i] ) < rmin )
        {
            imin = i;
            rmin = fabs( vec1[i] );
        }
        ++i;
    }

    if( cont_flag )
    {
        xv[imin]         = 1.0;
        xv[iperm1[imin]] = 0.0;
        xv[iperm2[imin]] = 0.0;

        // Determine cross product
        vec2[0] = vec1[1] * xv[2] - vec1[2] * xv[1];
        vec2[1] = vec1[2] * xv[0] - vec1[0] * xv[2];
        vec2[2] = vec1[0] * xv[1] - vec1[1] * xv[0];

        // Unitize
        rmag = sqrt( vec2[0] * vec2[0] + vec2[1] * vec2[1] + vec2[2] * vec2[2] );
        vec2[0] /= rmag;
        vec2[1] /= rmag;
        vec2[2] /= rmag;
    }

    // Copy 1st orthogonal vector into VerdictVector vector2
    vector2.set( vec2 );

    // Cross vectors to determine last orthogonal vector
    vector3 = vector1 * vector2;
}
void VerdictVector::perpendicular_z ( ) [inline]

Definition at line 340 of file VerdictVector.hpp.

References x(), and y().

{
    double temp = x();
    x( y() );
    y( -temp );
}
double & VerdictVector::r ( ) [inline]

Definition at line 264 of file VerdictVector.hpp.

References xVal.

Referenced by blow_out(), rtheta_to_xy(), scale_angle(), and xy_to_rtheta().

{
    return xVal;
}
void VerdictVector::r ( const double  xv) [inline]

Definition at line 284 of file VerdictVector.hpp.

References xVal.

{
    xVal = xv;
}
void VerdictVector::reflect_about_xaxis ( double  dummy,
double   
)

Definition at line 150 of file VerdictVector.cpp.

References yVal.

{
    yVal = -yVal;
}
void VerdictVector::rotate ( double  angle,
double   
)

Definition at line 126 of file VerdictVector.cpp.

References moab::angle(), rtheta_to_xy(), theta(), and xy_to_rtheta().

Definition at line 116 of file VerdictVector.cpp.

References r(), theta(), x(), and y().

Referenced by blow_out(), rotate(), and scale_angle().

{
    // careful about overwriting
    double x_ = r() * cos( theta() );
    double y_ = r() * sin( theta() );

    x( x_ );
    y( y_ );
}
void VerdictVector::scale_angle ( double  gamma,
double   
)

Definition at line 155 of file VerdictVector.cpp.

References r(), rtheta_to_xy(), theta(), TWO_VERDICT_PI, VERDICT_PI, and xy_to_rtheta().

{
    const double r_factor     = 0.3;
    const double theta_factor = 0.6;

    xy_to_rtheta();

    // if neary 2pi, treat as zero
    // some near zero stuff strays due to roundoff
    if( theta() > TWO_VERDICT_PI - 0.02 ) theta() = 0;
    // the above screws up on big sheets - need to overhaul at the sheet level

    if( gamma < 1 )
    {
        // squeeze together points of short radius so that
        // long chords won't cross them
        theta() += ( VERDICT_PI - theta() ) * ( 1 - gamma ) * theta_factor * ( 1 - r() );

        // push away from center of circle, again so long chords won't cross
        r( ( r_factor + r() ) / ( 1 + r_factor ) );

        // scale angle by gamma
        theta() *= gamma;
    }
    else
    {
        // scale angle by gamma, making sure points nearly 2pi are treated as zero
        double new_theta = theta() * gamma;
        if( new_theta < 2.5 * VERDICT_PI || r() < 0.2 ) theta( new_theta );
    }
    rtheta_to_xy();
}
void VerdictVector::set ( const double  xyz[3]) [inline]

Definition at line 354 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

{
    xVal = xyz[0];
    yVal = xyz[1];
    zVal = xyz[2];
}
void VerdictVector::set ( const VerdictVector tail,
const VerdictVector head 
) [inline]

Definition at line 361 of file VerdictVector.hpp.

References xVal, yVal, and zVal.

{
    xVal = head.xVal - tail.xVal;
    yVal = head.yVal - tail.yVal;
    zVal = head.zVal - tail.zVal;
}
void VerdictVector::set ( const VerdictVector to_copy) [inline]

Definition at line 376 of file VerdictVector.hpp.

{
    *this = to_copy;
}
double & VerdictVector::theta ( ) [inline]

Definition at line 268 of file VerdictVector.hpp.

References yVal.

Referenced by rotate(), rtheta_to_xy(), scale_angle(), and xy_to_rtheta().

{
    return yVal;
}
void VerdictVector::theta ( const double  yv) [inline]

Definition at line 288 of file VerdictVector.hpp.

References yVal.

{
    yVal = yv;
}
double VerdictVector::vector_angle ( const VerdictVector vector1,
const VerdictVector vector2 
) const

Definition at line 252 of file VerdictVector.cpp.

References moab::angle(), moab::dot(), length_squared(), normalize(), TWO_VERDICT_PI, and VERDICT_PI.

{
    // This routine does not assume that any of the input vectors are of unit
    // length. This routine does not normalize the input vectors.
    // Special cases:
    //     If the normal vector is zero length:
    //         If a new one can be computed from vectors 1 & 2:
    //             the normal is replaced with the vector cross product
    //         else the two vectors are colinear and zero or 2PI is returned.
    //     If the normal is colinear with either (or both) vectors
    //         a new one is computed with the cross products
    //         (and checked again).

    // Check for zero length normal vector
    VerdictVector normal = *this;
    double normal_lensq  = normal.length_squared();
    double len_tol       = 0.0000001;
    if( normal_lensq <= len_tol )
    {
        // null normal - make it the normal to the plane defined by vector1
        // and vector2. If still null, the vectors are colinear so check
        // for zero or 180 angle.
        normal       = vector1 * vector2;
        normal_lensq = normal.length_squared();
        if( normal_lensq <= len_tol )
        {
            double cosine = vector1 % vector2;
            if( cosine > 0.0 )
                return 0.0;
            else
                return VERDICT_PI;
        }
    }

    // Trap for normal vector colinear to one of the other vectors. If so,
    // use a normal defined by the two vectors.
    double dot_tol = 0.985;
    double dot     = vector1 % normal;
    if( dot * dot >= vector1.length_squared() * normal_lensq * dot_tol )
    {
        normal       = vector1 * vector2;
        normal_lensq = normal.length_squared();

        // Still problems if all three vectors were colinear
        if( normal_lensq <= len_tol )
        {
            double cosine = vector1 % vector2;
            if( cosine >= 0.0 )
                return 0.0;
            else
                return VERDICT_PI;
        }
    }
    else
    {
        // The normal and vector1 are not colinear, now check for vector2
        dot = vector2 % normal;
        if( dot * dot >= vector2.length_squared() * normal_lensq * dot_tol )
        {
            normal = vector1 * vector2;
        }
    }

    // Assume a plane such that the normal vector is the plane's normal.
    // Create yAxis perpendicular to both the normal and vector1. yAxis is
    // now in the plane. Create xAxis as the perpendicular to both yAxis and
    // the normal. xAxis is in the plane and is the projection of vector1
    // into the plane.

    normal.normalize();
    VerdictVector yAxis = normal;
    yAxis *= vector1;
    double yv = vector2 % yAxis;
    //  yAxis memory slot will now be used for xAxis
    yAxis *= normal;
    double xv = vector2 % yAxis;

    //  assert(x != 0.0 || y != 0.0);
    if( xv == 0.0 && yv == 0.0 )
    {
        return 0.0;
    }
    double angle = atan2( yv, xv );

    if( angle < 0.0 )
    {
        angle += TWO_VERDICT_PI;
    }
    return angle;
}
double VerdictVector::vector_angle_quick ( const VerdictVector vec1,
const VerdictVector vec2 
)

Definition at line 188 of file VerdictVector.cpp.

References moab::angle(), and TWO_VERDICT_PI.

{
    //- compute the angle between two vectors in the plane defined by this vector
    // build yAxis and xAxis such that xAxis is the projection of
    // vec1 onto the normal plane of this vector

    // NOTE: vec1 and vec2 are Vectors from the vertex of the angle along
    //       the two sides of the angle.
    //       The angle returned is the right-handed angle around this vector
    //       from vec1 to vec2.

    // NOTE: vector_angle_quick gives exactly the same answer as vector_angle below
    //       providing this vector is normalized.  It does so with two fewer
    //       cross-product evaluations and two fewer vector normalizations.
    //       This can be a substantial time savings if the function is called
    //       a significant number of times (e.g Hexer) ... (jrh 11/28/94)
    // NOTE: vector_angle() is much more robust. Do not use vector_angle_quick()
    //       unless you are very sure of the safety of your input vectors.

    VerdictVector ry = ( *this ) * vec1;
    VerdictVector rx = ry * ( *this );

    double xv = vec2 % rx;
    double yv = vec2 % ry;

    double angle;
    assert( xv != 0.0 || yv != 0.0 );

    angle = atan2( yv, xv );

    if( angle < 0.0 )
    {
        angle += TWO_VERDICT_PI;
    }
    return angle;
}
bool VerdictVector::within_tolerance ( const VerdictVector vectorPtr2,
double  tolerance 
) const

Definition at line 343 of file VerdictVector.cpp.

References x(), y(), and z().

{
    if( ( fabs( this->x() - vectorPtr2.x() ) < tolerance ) && ( fabs( this->y() - vectorPtr2.y() ) < tolerance ) &&
        ( fabs( this->z() - vectorPtr2.z() ) < tolerance ) )
    {
        return true;
    }

    return false;
}
void VerdictVector::x ( const double  xv) [inline]

Definition at line 272 of file VerdictVector.hpp.

References xVal.

{
    xVal = xv;
}

Definition at line 105 of file VerdictVector.cpp.

References length(), r(), theta(), TWO_VERDICT_PI, x(), and y().

Referenced by blow_out(), rotate(), and scale_angle().

{
    // careful about overwriting
    double r_     = length();
    double theta_ = atan2( y(), x() );
    if( theta_ < 0.0 ) theta_ += TWO_VERDICT_PI;

    r( r_ );
    theta( theta_ );
}
void VerdictVector::y ( const double  yv) [inline]

Definition at line 276 of file VerdictVector.hpp.

References yVal.

{
    yVal = yv;
}
void VerdictVector::z ( const double  zv) [inline]

Definition at line 280 of file VerdictVector.hpp.

References zVal.

{
    zVal = zv;
}

Friends And Related Function Documentation

VerdictVector interpolate ( const double  param,
const VerdictVector v1,
const VerdictVector v2 
) [friend]

Definition at line 98 of file VerdictVector.cpp.

{
    VerdictVector temp = ( 1.0 - param ) * v1;
    temp += param * v2;
    return temp;
}
int operator!= ( const VerdictVector v1,
const VerdictVector v2 
) [friend]

Definition at line 467 of file VerdictVector.hpp.

{
    return ( v1.xVal != v2.xVal || v1.yVal != v2.yVal || v1.zVal != v2.zVal );
}
double operator% ( const VerdictVector v1,
const VerdictVector v2 
) [friend]

Definition at line 495 of file VerdictVector.hpp.

{
    return ( vector1.x() * vector2.x() + vector1.y() * vector2.y() + vector1.z() * vector2.z() );
}
VerdictVector operator* ( const VerdictVector v1,
const VerdictVector v2 
) [friend]

Definition at line 439 of file VerdictVector.hpp.

{
    return VerdictVector( vector1 ) *= vector2;
}
VerdictVector operator* ( const VerdictVector v1,
const double  sclr 
) [friend]

Definition at line 445 of file VerdictVector.hpp.

{
    return VerdictVector( vector1 ) *= scalar;
}
VerdictVector operator* ( const double  sclr,
const VerdictVector v1 
) [friend]

Definition at line 451 of file VerdictVector.hpp.

{
    return VerdictVector( vector1 ) *= scalar;
}
VerdictVector operator+ ( const VerdictVector v1,
const VerdictVector v2 
) [friend]

Definition at line 419 of file VerdictVector.hpp.

{
    double xv = vector1.x() + vector2.x();
    double yv = vector1.y() + vector2.y();
    double zv = vector1.z() + vector2.z();
    return VerdictVector( xv, yv, zv );
    //  return VerdictVector(vector1) += vector2;
}
VerdictVector operator- ( const VerdictVector v1,
const VerdictVector v2 
) [friend]

Definition at line 428 of file VerdictVector.hpp.

{
    double xv = vector1.x() - vector2.x();
    double yv = vector1.y() - vector2.y();
    double zv = vector1.z() - vector2.z();
    return VerdictVector( xv, yv, zv );
    //  return VerdictVector(vector1) -= vector2;
}
VerdictVector operator/ ( const VerdictVector v1,
const double  sclr 
) [friend]

Definition at line 457 of file VerdictVector.hpp.

{
    return VerdictVector( vector1 ) /= scalar;
}
int operator== ( const VerdictVector v1,
const VerdictVector v2 
) [friend]

Definition at line 462 of file VerdictVector.hpp.

{
    return ( v1.xVal == v2.xVal && v1.yVal == v2.yVal && v1.zVal == v2.zVal );
}
VerdictVector operator~ ( const VerdictVector vec) [friend]

Definition at line 401 of file VerdictVector.hpp.

{
    double mag = sqrt( vec.xVal * vec.xVal + vec.yVal * vec.yVal + vec.zVal * vec.zVal );

    VerdictVector temp = vec;
    if( mag != 0.0 )
    {
        temp /= mag;
    }
    return temp;
}

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