Mesh Oriented datABase
(version 5.4.1)
Array-based unstructured mesh datastructure
|
#include <Matrix3.hpp>
Public Member Functions | |
Matrix3 () | |
Matrix3 (Eigen::Matrix3d mat) | |
Matrix3 (double diagonal) | |
Matrix3 (const CartVect &diagonal) | |
Matrix3 (const std::vector< double > &diagonal) | |
Matrix3 (double v00, double v01, double v02, double v10, double v11, double v12, double v20, double v21, double v22) | |
Matrix3 (const Matrix3 &f) | |
template<typename Vector > | |
Matrix3 (const Vector &row0, const Vector &row1, const Vector &row2, const bool isRow) | |
Matrix3 (const double v[size]) | |
void | copyto (double v[Matrix3::size]) |
Matrix3 & | operator= (const Matrix3 &m) |
Matrix3 & | operator= (const double v[size]) |
double * | operator[] (unsigned i) |
const double * | operator[] (unsigned i) const |
double & | operator() (unsigned r, unsigned c) |
double | operator() (unsigned r, unsigned c) const |
double & | operator() (unsigned i) |
double | operator() (unsigned i) const |
double * | array () |
const double * | array () const |
Matrix3 & | operator+= (const Matrix3 &m) |
Matrix3 & | operator-= (const Matrix3 &m) |
Matrix3 & | operator*= (double s) |
Matrix3 & | operator/= (double s) |
Matrix3 & | operator*= (const Matrix3 &m) |
bool | is_symmetric () |
bool | is_positive_definite () |
template<typename Vector > | |
ErrorCode | eigen_decomposition (Vector &evals, Matrix3 &evecs) |
void | transpose_inplace () |
Matrix3 | transpose () const |
template<typename Vector > | |
void | copycol (int index, Vector &vol) |
void | swapcol (int srcindex, int destindex) |
template<typename Vector > | |
Vector | vcol (int index) const |
void | colscale (int index, double scale) |
void | rowscale (int index, double scale) |
CartVect | col (int index) const |
CartVect | row (int index) const |
double | determinant () const |
Matrix3 | inverse () const |
bool | invert () |
double | subdet (int r, int c) const |
void | print (std::ostream &s) const |
Static Public Attributes | |
static const int | size = 9 |
Private Attributes | |
Eigen::Matrix3d | _mat |
Friends | |
Matrix3 | operator+ (const Matrix3 &a, const Matrix3 &b) |
Matrix3 | operator- (const Matrix3 &a, const Matrix3 &b) |
Matrix3 | operator* (const Matrix3 &a, const Matrix3 &b) |
Definition at line 216 of file Matrix3.hpp.
moab::Matrix3::Matrix3 | ( | ) | [inline] |
Definition at line 231 of file Matrix3.hpp.
Referenced by inverse(), and transpose().
{ #ifndef MOAB_HAVE_LAPACK _mat.fill( 0.0 ); #else MOAB_DMEMZERO( _mat, Matrix3::size ); #endif }
moab::Matrix3::Matrix3 | ( | Eigen::Matrix3d | mat | ) | [inline] |
Definition at line 241 of file Matrix3.hpp.
: _mat( mat ) {}
moab::Matrix3::Matrix3 | ( | double | diagonal | ) | [inline] |
Definition at line 246 of file Matrix3.hpp.
moab::Matrix3::Matrix3 | ( | const CartVect & | diagonal | ) | [inline] |
Definition at line 256 of file Matrix3.hpp.
moab::Matrix3::Matrix3 | ( | const std::vector< double > & | diagonal | ) | [inline] |
Definition at line 273 of file Matrix3.hpp.
moab::Matrix3::Matrix3 | ( | double | v00, |
double | v01, | ||
double | v02, | ||
double | v10, | ||
double | v11, | ||
double | v12, | ||
double | v20, | ||
double | v21, | ||
double | v22 | ||
) | [inline] |
Definition at line 285 of file Matrix3.hpp.
moab::Matrix3::Matrix3 | ( | const Matrix3 & | f | ) | [inline] |
moab::Matrix3::Matrix3 | ( | const Vector & | row0, |
const Vector & | row1, | ||
const Vector & | row2, | ||
const bool | isRow | ||
) | [inline] |
Definition at line 323 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK if( isRow ) { _mat << row0[0], row0[1], row0[2], row1[0], row1[1], row1[2], row2[0], row2[1], row2[2]; } else { _mat << row0[0], row1[0], row2[0], row0[1], row1[1], row2[1], row0[2], row1[2], row2[2]; } #else MOAB_DMEMZERO( _mat, Matrix3::size ); if( isRow ) { _mat[0] = row0[0]; _mat[1] = row0[1]; _mat[2] = row0[2]; _mat[3] = row1[0]; _mat[4] = row1[1]; _mat[5] = row1[2]; _mat[6] = row2[0]; _mat[7] = row2[1]; _mat[8] = row2[2]; } else { _mat[0] = row0[0]; _mat[1] = row1[0]; _mat[2] = row2[0]; _mat[3] = row0[1]; _mat[4] = row1[1]; _mat[5] = row2[1]; _mat[6] = row0[2]; _mat[7] = row1[2]; _mat[8] = row2[2]; } #endif }
moab::Matrix3::Matrix3 | ( | const double | v[size] | ) | [inline] |
double* moab::Matrix3::array | ( | ) | [inline] |
Definition at line 469 of file Matrix3.hpp.
References _mat.
Referenced by moab::LinearTet::evaluate_reverse(), and moab::EvalSet::evaluate_reverse().
const double* moab::Matrix3::array | ( | ) | const [inline] |
Definition at line 478 of file Matrix3.hpp.
References _mat.
CartVect moab::Matrix3::col | ( | int | index | ) | const [inline] |
Definition at line 886 of file Matrix3.hpp.
Referenced by moab::OrientedBox::area(), moab::OrientedBox::axis(), moab::box_from_axes(), moab::OrientedBox::closest_location_in_box(), moab::OrientedBox::contained(), moab::OrientedBox::dimensions(), moab::OrientedBox::inner_radius(), moab::OrientedBox::inner_radius_squared(), moab::OrientedBox::make_hex(), moab::operator<<(), moab::OrientedBox::OrientedBox(), moab::OrientedBox::outer_radius(), moab::OrientedBox::outer_radius_squared(), moab::OrientedBox::scaled_axis(), and moab::OrientedBox::volume().
{ assert( index < Matrix3::size ); #ifndef MOAB_HAVE_LAPACK Eigen::Vector3d mvec = _mat.col( index ); return CartVect( mvec[0], mvec[1], mvec[2] ); #else switch( index ) { case 0: return CartVect( _mat[0], _mat[3], _mat[6] ); case 1: return CartVect( _mat[1], _mat[4], _mat[7] ); case 2: return CartVect( _mat[2], _mat[5], _mat[8] ); } return CartVect( 0.0 ); #endif }
void moab::Matrix3::colscale | ( | int | index, |
double | scale | ||
) | [inline] |
Definition at line 832 of file Matrix3.hpp.
Referenced by moab::box_from_axes(), and moab::OrientedBox::order_axes_by_length().
{ assert( index < Matrix3::size ); #ifndef MOAB_HAVE_LAPACK _mat.col( index ) *= scale; #else switch( index ) { case 0: _mat[0] *= scale; _mat[3] *= scale; _mat[6] *= scale; break; case 1: _mat[1] *= scale; _mat[4] *= scale; _mat[7] *= scale; break; case 2: _mat[2] *= scale; _mat[5] *= scale; _mat[8] *= scale; break; } #endif }
void moab::Matrix3::copycol | ( | int | index, |
Vector & | vol | ||
) | [inline] |
Definition at line 738 of file Matrix3.hpp.
References _mat.
void moab::Matrix3::copyto | ( | double | v[Matrix3::size] | ) | [inline] |
Definition at line 385 of file Matrix3.hpp.
Referenced by moab::LinearTet::initFcn(), and moab::LinearTri::initFcn().
double moab::Matrix3::determinant | ( | ) | const [inline] |
Definition at line 930 of file Matrix3.hpp.
References _mat.
Referenced by moab::Element::Map::det_ijacobian(), moab::Element::Map::det_jacobian(), moab::LinearTet::evaluate_reverse(), moab::LinearTri::evaluate_reverse(), moab::EvalSet::evaluate_reverse(), moab::Element::Map::ievaluate(), moab::LinearTet::initFcn(), moab::LinearTri::initFcn(), moab::element_utility::Spectral_hex_map< moab::Matrix3 >::integrate_scalar_field(), moab::Element::SpectralHex::integrate_scalar_field(), moab::ElemUtil::integrate_trilinear_hex(), moab::LinearQuad::integrateFcn(), moab::LinearHex::integrateFcn(), inverse(), invert(), moab::AffineXform::reflection(), moab::AffineXform::scale(), moab::Element::LinearTet::set_vertices(), moab::Element::LinearTri::set_vertices(), moab::ElemUtil::VolMap::solve_inverse(), moab::GeomUtil::VolMap::solve_inverse(), and test_spectral_hex().
ErrorCode moab::Matrix3::eigen_decomposition | ( | Vector & | evals, |
Matrix3 & | evecs | ||
) | [inline] |
Definition at line 589 of file Matrix3.hpp.
References _mat, is_symmetric(), MB_CHK_SET_ERR, MB_SUCCESS, and size.
Referenced by moab::OrientedBox::compute_from_covariance_data(), and moab::OrientedBox::compute_from_vertices().
{ const bool bisSymmetric = this->is_symmetric(); #ifndef MOAB_HAVE_LAPACK if( bisSymmetric ) { Eigen::SelfAdjointEigenSolver< Eigen::Matrix3d > eigensolver( this->_mat ); if( eigensolver.info() != Eigen::Success ) return MB_FAILURE; const Eigen::SelfAdjointEigenSolver< Eigen::Matrix3d >::RealVectorType& e3evals = eigensolver.eigenvalues(); evals[0] = e3evals( 0 ); evals[1] = e3evals( 1 ); evals[2] = e3evals( 2 ); evecs._mat = eigensolver.eigenvectors(); //.col(1) return MB_SUCCESS; } else { MB_CHK_SET_ERR( MB_FAILURE, "Unsymmetric matrix implementation with Eigen3 is currently not provided." ); // Eigen::EigenSolver<Eigen::Matrix3d> eigensolver(this->_mat, true); // if (eigensolver.info() != Eigen::Success) // return MB_FAILURE; // const Eigen::EigenSolver<Eigen::Matrix3d>::EigenvalueType& e3evals = // eigensolver.eigenvalues().real(); evals[0] = e3evals(0); evals[1] = e3evals(1); // evals[2] = e3evals(2); evecs._mat = eigensolver.eigenvectors().real(); //.col(1) // return MB_SUCCESS; } #else int info; /* Solve eigenproblem */ double devreal[3], drevecs[9]; if( !bisSymmetric ) { double devimag[3], dlevecs[9], dwork[102]; char dgeev_opts[2] = { 'N', 'V' }; int N = 3, LWORK = 102, NL = 1, NR = N; std::vector< double > devmat; devmat.assign( _mat, _mat + size ); MOAB_dgeev( &dgeev_opts[0], &dgeev_opts[1], &N, &devmat[0], &N, devreal, devimag, dlevecs, &NL, drevecs, &NR, dwork, &LWORK, &info ); // The result eigenvalues are ordered as high-->low evals[0] = devreal[2]; evals[1] = devreal[1]; evals[2] = devreal[0]; evecs._mat[0] = drevecs[6]; evecs._mat[1] = drevecs[3]; evecs._mat[2] = drevecs[0]; evecs._mat[3] = drevecs[7]; evecs._mat[4] = drevecs[4]; evecs._mat[5] = drevecs[1]; evecs._mat[6] = drevecs[8]; evecs._mat[7] = drevecs[5]; evecs._mat[8] = drevecs[2]; std::cout << "DGEEV: Optimal work vector: dsize = " << dwork[0] << ".\n"; } else { char dgeev_opts[2] = { 'V', 'L' }; const bool find_optimal = false; std::vector< int > iwork( 18 ); std::vector< double > devmat( 9, 0.0 ); std::vector< double > dwork( 38 ); int N = 3, lwork = 38, liwork = 18; devmat[0] = _mat[0]; devmat[1] = _mat[1]; devmat[2] = _mat[2]; devmat[4] = _mat[4]; devmat[5] = _mat[5]; devmat[8] = _mat[8]; if( find_optimal ) { int _lwork = -1; int _liwork = -1; double query_work_size = 0; int query_iwork_size = 0; // Make an empty call to find the optimal work vector size MOAB_dsyevd( &dgeev_opts[0], &dgeev_opts[1], &N, NULL, &N, NULL, &query_work_size, &_lwork, &query_iwork_size, &_liwork, &info ); lwork = (int)query_work_size; dwork.resize( lwork ); liwork = query_iwork_size; iwork.resize( liwork ); std::cout << "DSYEVD: Optimal work vector: dsize = " << lwork << ", and isize = " << liwork << ".\n"; } MOAB_dsyevd( &dgeev_opts[0], &dgeev_opts[1], &N, &devmat[0], &N, devreal, &dwork[0], &lwork, &iwork[0], &liwork, &info ); for( int i = 0; i < 9; ++i ) drevecs[i] = devmat[i]; // The result eigenvalues are ordered as low-->high, but vectors are in rows of A. evals[0] = devreal[0]; evals[1] = devreal[1]; evals[2] = devreal[2]; evecs._mat[0] = drevecs[0]; evecs._mat[3] = drevecs[1]; evecs._mat[6] = drevecs[2]; evecs._mat[1] = drevecs[3]; evecs._mat[4] = drevecs[4]; evecs._mat[7] = drevecs[5]; evecs._mat[2] = drevecs[6]; evecs._mat[5] = drevecs[7]; evecs._mat[8] = drevecs[8]; } if( !info ) { return MB_SUCCESS; } else { std::cout << "Failure in LAPACK_" << ( bisSymmetric ? "DSYEVD" : "DGEEV" ) << " call for eigen decomposition.\n"; std::cout << "Failed with error = " << info << ".\n"; return MB_FAILURE; } #endif }
Matrix3 moab::Matrix3::inverse | ( | ) | const [inline] |
Definition at line 940 of file Matrix3.hpp.
References _mat, determinant(), and Matrix3().
Referenced by moab::Element::Map::det_ijacobian(), moab::LinearTet::evaluate_reverse(), moab::LinearTri::evaluate_reverse(), moab::EvalSet::evaluate_reverse(), moab::Element::Map::ievaluate(), moab::Element::Map::ijacobian(), moab::LinearTet::initFcn(), moab::LinearTri::initFcn(), moab::AffineXform::inverse(), moab::Element::LinearTet::set_vertices(), moab::Element::LinearTri::set_vertices(), moab::ElemUtil::VolMap::solve_inverse(), and moab::GeomUtil::VolMap::solve_inverse().
{ #ifndef MOAB_HAVE_LAPACK return Matrix3( _mat.inverse() ); #else // return Matrix::compute_inverse( *this, this->determinant() ); Matrix3 m( 0.0 ); const double d_determinant = 1.0 / this->determinant(); m._mat[0] = d_determinant * ( _mat[4] * _mat[8] - _mat[5] * _mat[7] ); m._mat[1] = d_determinant * ( _mat[2] * _mat[7] - _mat[8] * _mat[1] ); m._mat[2] = d_determinant * ( _mat[1] * _mat[5] - _mat[4] * _mat[2] ); m._mat[3] = d_determinant * ( _mat[5] * _mat[6] - _mat[8] * _mat[3] ); m._mat[4] = d_determinant * ( _mat[0] * _mat[8] - _mat[6] * _mat[2] ); m._mat[5] = d_determinant * ( _mat[2] * _mat[3] - _mat[5] * _mat[0] ); m._mat[6] = d_determinant * ( _mat[3] * _mat[7] - _mat[6] * _mat[4] ); m._mat[7] = d_determinant * ( _mat[1] * _mat[6] - _mat[7] * _mat[0] ); m._mat[8] = d_determinant * ( _mat[0] * _mat[4] - _mat[3] * _mat[1] ); return m; #endif }
bool moab::Matrix3::invert | ( | ) | [inline] |
Definition at line 961 of file Matrix3.hpp.
References _mat, determinant(), moab::Util::is_finite(), and size.
{ bool invertible = false; double d_determinant; #ifndef MOAB_HAVE_LAPACK Eigen::Matrix3d invMat; _mat.computeInverseAndDetWithCheck( invMat, d_determinant, invertible ); if( !Util::is_finite( d_determinant ) ) return false; _mat = invMat; return invertible; #else d_determinant = this->determinant(); if( d_determinant > 1e-13 ) invertible = true; d_determinant = 1.0 / d_determinant; // invert the determinant std::vector< double > _m; _m.assign( _mat, _mat + size ); _mat[0] = d_determinant * ( _m[4] * _m[8] - _m[5] * _m[7] ); _mat[1] = d_determinant * ( _m[2] * _m[7] - _m[8] * _m[1] ); _mat[2] = d_determinant * ( _m[1] * _m[5] - _m[4] * _m[2] ); _mat[3] = d_determinant * ( _m[5] * _m[6] - _m[8] * _m[3] ); _mat[4] = d_determinant * ( _m[0] * _m[8] - _m[6] * _m[2] ); _mat[5] = d_determinant * ( _m[2] * _m[3] - _m[5] * _m[0] ); _mat[6] = d_determinant * ( _m[3] * _m[7] - _m[6] * _m[4] ); _mat[7] = d_determinant * ( _m[1] * _m[6] - _m[7] * _m[0] ); _mat[8] = d_determinant * ( _m[0] * _m[4] - _m[3] * _m[1] ); #endif return invertible; }
bool moab::Matrix3::is_positive_definite | ( | ) | [inline] |
Definition at line 569 of file Matrix3.hpp.
References _mat.
{ #ifndef MOAB_HAVE_LAPACK double subdet6 = _mat( 1 ) * _mat( 5 ) - _mat( 2 ) * _mat( 4 ); double subdet7 = _mat( 2 ) * _mat( 3 ) - _mat( 0 ) * _mat( 5 ); double subdet8 = _mat( 0 ) * _mat( 4 ) - _mat( 1 ) * _mat( 3 ); // Determinant:= d(6)*subdet6 + d(7)*subdet7 + d(8)*subdet8; const double det = _mat( 6 ) * subdet6 + _mat( 7 ) * subdet7 + _mat( 8 ) * subdet8; return _mat( 0 ) > 0 && subdet8 > 0 && det > 0; #else double subdet6 = _mat[1] * _mat[5] - _mat[2] * _mat[4]; double subdet7 = _mat[2] * _mat[3] - _mat[0] * _mat[5]; double subdet8 = _mat[0] * _mat[4] - _mat[1] * _mat[3]; // Determinant:= d(6)*subdet6 + d(7)*subdet7 + d(8)*subdet8; const double det = _mat[6] * subdet6 + _mat[7] * subdet7 + _mat[8] * subdet8; return _mat[0] > 0 && subdet8 > 0 && det > 0; #endif }
bool moab::Matrix3::is_symmetric | ( | ) | [inline] |
Definition at line 553 of file Matrix3.hpp.
Referenced by eigen_decomposition().
{ const double EPS = 1e-13; #ifndef MOAB_HAVE_LAPACK if( ( fabs( _mat( 1 ) - _mat( 3 ) ) < EPS ) && ( fabs( _mat( 2 ) - _mat( 6 ) ) < EPS ) && ( fabs( _mat( 5 ) - _mat( 7 ) ) < EPS ) ) return true; #else if( ( fabs( _mat[1] - _mat[3] ) < EPS ) && ( fabs( _mat[2] - _mat[6] ) < EPS ) && ( fabs( _mat[5] - _mat[7] ) < EPS ) ) return true; #endif else return false; }
double& moab::Matrix3::operator() | ( | unsigned | r, |
unsigned | c | ||
) | [inline] |
Definition at line 432 of file Matrix3.hpp.
References _mat.
double moab::Matrix3::operator() | ( | unsigned | r, |
unsigned | c | ||
) | const [inline] |
Definition at line 441 of file Matrix3.hpp.
References _mat.
double& moab::Matrix3::operator() | ( | unsigned | i | ) | [inline] |
Definition at line 450 of file Matrix3.hpp.
References _mat.
double moab::Matrix3::operator() | ( | unsigned | i | ) | const [inline] |
Definition at line 459 of file Matrix3.hpp.
References _mat.
Matrix3& moab::Matrix3::operator*= | ( | double | s | ) | [inline] |
Definition at line 509 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK _mat *= s; #else for( int i = 0; i < Matrix3::size; ++i ) _mat[i] *= s; #endif return *this; }
Definition at line 531 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK _mat *= m._mat; #else // Uncomment below if you want point-wise multiplication instead (.*) // for (int i=0; i < Matrix3::size; ++i) _mat[i] *= m._mat[i]; std::vector< double > dmat; dmat.assign( _mat, _mat + size ); _mat[0] = dmat[0] * m._mat[0] + dmat[1] * m._mat[3] + dmat[2] * m._mat[6]; _mat[1] = dmat[0] * m._mat[1] + dmat[1] * m._mat[4] + dmat[2] * m._mat[7]; _mat[2] = dmat[0] * m._mat[2] + dmat[1] * m._mat[5] + dmat[2] * m._mat[8]; _mat[3] = dmat[3] * m._mat[0] + dmat[4] * m._mat[3] + dmat[5] * m._mat[6]; _mat[4] = dmat[3] * m._mat[1] + dmat[4] * m._mat[4] + dmat[5] * m._mat[7]; _mat[5] = dmat[3] * m._mat[2] + dmat[4] * m._mat[5] + dmat[5] * m._mat[8]; _mat[6] = dmat[6] * m._mat[0] + dmat[7] * m._mat[3] + dmat[8] * m._mat[6]; _mat[7] = dmat[6] * m._mat[1] + dmat[7] * m._mat[4] + dmat[8] * m._mat[7]; _mat[8] = dmat[6] * m._mat[2] + dmat[7] * m._mat[5] + dmat[8] * m._mat[8]; #endif return *this; }
Definition at line 487 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK _mat += m._mat; #else for( int i = 0; i < Matrix3::size; ++i ) _mat[i] += m._mat[i]; #endif return *this; }
Definition at line 498 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK _mat -= m._mat; #else for( int i = 0; i < Matrix3::size; ++i ) _mat[i] -= m._mat[i]; #endif return *this; }
Matrix3& moab::Matrix3::operator/= | ( | double | s | ) | [inline] |
Definition at line 520 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK _mat /= s; #else for( int i = 0; i < Matrix3::size; ++i ) _mat[i] /= s; #endif return *this; }
Matrix3& moab::Matrix3::operator= | ( | const double | v[size] | ) | [inline] |
double* moab::Matrix3::operator[] | ( | unsigned | i | ) | [inline] |
Definition at line 414 of file Matrix3.hpp.
References _mat.
const double* moab::Matrix3::operator[] | ( | unsigned | i | ) | const [inline] |
Definition at line 423 of file Matrix3.hpp.
References _mat.
void moab::Matrix3::print | ( | std::ostream & | s | ) | const [inline] |
Definition at line 1007 of file Matrix3.hpp.
References _mat.
{ #ifndef MOAB_HAVE_LAPACK s << "| " << _mat( 0 ) << " " << _mat( 1 ) << " " << _mat( 2 ) << " | " << _mat( 3 ) << " " << _mat( 4 ) << " " << _mat( 5 ) << " | " << _mat( 6 ) << " " << _mat( 7 ) << " " << _mat( 8 ) << " |"; #else s << "| " << _mat[0] << " " << _mat[1] << " " << _mat[2] << " | " << _mat[3] << " " << _mat[4] << " " << _mat[5] << " | " << _mat[6] << " " << _mat[7] << " " << _mat[8] << " |"; #endif }
CartVect moab::Matrix3::row | ( | int | index | ) | const [inline] |
Definition at line 906 of file Matrix3.hpp.
{ assert( index < Matrix3::size ); #ifndef MOAB_HAVE_LAPACK Eigen::Vector3d mvec = _mat.row( index ); return CartVect( mvec[0], mvec[1], mvec[2] ); #else switch( index ) { case 0: return CartVect( _mat[0], _mat[1], _mat[2] ); case 1: return CartVect( _mat[3], _mat[4], _mat[5] ); case 2: return CartVect( _mat[6], _mat[7], _mat[8] ); } return CartVect( 0.0 ); #endif }
void moab::Matrix3::rowscale | ( | int | index, |
double | scale | ||
) | [inline] |
Definition at line 859 of file Matrix3.hpp.
{ assert( index < Matrix3::size ); #ifndef MOAB_HAVE_LAPACK _mat.row( index ) *= scale; #else switch( index ) { case 0: _mat[0] *= scale; _mat[1] *= scale; _mat[2] *= scale; break; case 1: _mat[3] *= scale; _mat[4] *= scale; _mat[5] *= scale; break; case 2: _mat[6] *= scale; _mat[7] *= scale; _mat[8] *= scale; break; } #endif }
double moab::Matrix3::subdet | ( | int | r, |
int | c | ||
) | const [inline] |
Definition at line 992 of file Matrix3.hpp.
{ assert( r >= 0 && c >= 0 ); if( r < 0 || c < 0 ) return DBL_MAX; #ifndef MOAB_HAVE_LAPACK const int r1 = ( r + 1 ) % 3, r2 = ( r + 2 ) % 3; const int c1 = ( c + 1 ) % 3, c2 = ( c + 2 ) % 3; return _mat( r1, c1 ) * _mat( r2, c2 ) - _mat( r1, c2 ) * _mat( r2, c1 ); #else const int r1 = Matrix3::size * ( ( r + 1 ) % 3 ), r2 = Matrix3::size * ( ( r + 2 ) % 3 ); const int c1 = ( c + 1 ) % 3, c2 = ( c + 2 ) % 3; return _mat[r1 + c1] * _mat[r2 + c2] - _mat[r1 + c2] * _mat[r2 + c1]; #endif }
void moab::Matrix3::swapcol | ( | int | srcindex, |
int | destindex | ||
) | [inline] |
Definition at line 764 of file Matrix3.hpp.
Referenced by moab::box_from_axes(), and moab::OrientedBox::order_axes_by_length().
{ assert( srcindex < Matrix3::size ); assert( destindex < Matrix3::size ); #ifndef MOAB_HAVE_LAPACK _mat.col( srcindex ).swap( _mat.col( destindex ) ); #else CartVect svol = this->vcol< CartVect >( srcindex ); CartVect dvol = this->vcol< CartVect >( destindex ); switch( srcindex ) { case 0: _mat[0] = dvol[0]; _mat[3] = dvol[1]; _mat[6] = dvol[2]; break; case 1: _mat[1] = dvol[0]; _mat[4] = dvol[1]; _mat[7] = dvol[2]; break; case 2: _mat[2] = dvol[0]; _mat[5] = dvol[1]; _mat[8] = dvol[2]; break; } switch( destindex ) { case 0: _mat[0] = svol[0]; _mat[3] = svol[1]; _mat[6] = svol[2]; break; case 1: _mat[1] = svol[0]; _mat[4] = svol[1]; _mat[7] = svol[2]; break; case 2: _mat[2] = svol[0]; _mat[5] = svol[1]; _mat[8] = svol[2]; break; } #endif }
Matrix3 moab::Matrix3::transpose | ( | ) | const [inline] |
Definition at line 721 of file Matrix3.hpp.
References _mat, and Matrix3().
Referenced by moab::OrientedBox::intersect_ray().
void moab::Matrix3::transpose_inplace | ( | ) | [inline] |
Vector moab::Matrix3::vcol | ( | int | index | ) | const [inline] |
Definition at line 813 of file Matrix3.hpp.
Definition at line 1051 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK return Matrix3( a._mat * b._mat ); #else return Matrix::mmult3( a, b ); #endif }
Definition at line 1027 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK return Matrix3( a._mat + b._mat ); #else Matrix3 s( a ); for( int i = 0; i < Matrix3::size; ++i ) s( i ) += b._mat[i]; return s; #endif }
Definition at line 1039 of file Matrix3.hpp.
{ #ifndef MOAB_HAVE_LAPACK return Matrix3( a._mat - b._mat ); #else Matrix3 s( a ); for( int i = 0; i < Matrix3::size; ++i ) s( i ) -= b._mat[i]; return s; #endif }
Eigen::Matrix3d moab::Matrix3::_mat [private] |
Definition at line 224 of file Matrix3.hpp.
Referenced by array(), col(), colscale(), copycol(), copyto(), determinant(), eigen_decomposition(), inverse(), invert(), is_positive_definite(), is_symmetric(), Matrix3(), operator()(), moab::operator*(), operator*=(), moab::operator+(), operator+=(), moab::operator-(), operator-=(), operator/=(), operator=(), operator[](), print(), row(), rowscale(), subdet(), swapcol(), transpose(), transpose_inplace(), and vcol().
const int moab::Matrix3::size = 9 [static] |
Definition at line 220 of file Matrix3.hpp.
Referenced by col(), colscale(), copyto(), eigen_decomposition(), moab::LinearTri::initFcn(), moab::LinearTet::initFcn(), invert(), Matrix3(), operator*=(), moab::operator+(), operator+=(), moab::operator-(), operator-=(), operator/=(), operator=(), row(), rowscale(), subdet(), swapcol(), and vcol().