$extrastylesheet
#include <tensor_value.h>

Public Member Functions | |
| TensorValue () | |
| TensorValue (const T xx, const T xy=0, const T xz=0, const T yx=0, const T yy=0, const T yz=0, const T zx=0, const T zy=0, const T zz=0) | |
| template<typename Scalar > | |
| TensorValue (const Scalar xx, const Scalar xy=0, const Scalar xz=0, const Scalar yx=0, const Scalar yy=0, const Scalar yz=0, const Scalar zx=0, const Scalar zy=0, typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, const Scalar >::type zz=0) | |
| template<typename T2 > | |
| TensorValue (const TypeVector< T2 > &vx) | |
| template<typename T2 > | |
| TensorValue (const TypeVector< T2 > &vx, const TypeVector< T2 > &vy) | |
| template<typename T2 > | |
| TensorValue (const TypeVector< T2 > &vx, const TypeVector< T2 > &vy, const TypeVector< T2 > &vz) | |
| template<typename T2 > | |
| TensorValue (const TensorValue< T2 > &p) | |
| template<typename T2 > | |
| TensorValue (const TypeTensor< T2 > &p) | |
| TensorValue (const TypeTensor< Real > &p_re, const TypeTensor< Real > &p_im) | |
| template<typename Scalar > | |
| boostcopy::enable_if_c < ScalarTraits< Scalar > ::value, TensorValue & >::type | operator= (const Scalar &libmesh_dbg_var(p)) |
| template<typename T2 > | |
| void | assign (const TypeTensor< T2 > &) |
| const T & | operator() (const unsigned int i, const unsigned int j) const |
| T & | operator() (const unsigned int i, const unsigned int j) |
| ConstTypeTensorColumn< T > | slice (const unsigned int i) const |
| TypeTensorColumn< T > | slice (const unsigned int i) |
| TypeVector< T > | row (const unsigned int r) |
| template<typename T2 > | |
| TypeTensor< typename CompareTypes< T, T2 > ::supertype > | operator+ (const TypeTensor< T2 > &) const |
| template<typename T2 > | |
| const TypeTensor< T > & | operator+= (const TypeTensor< T2 > &) |
| template<typename T2 > | |
| void | add (const TypeTensor< T2 > &) |
| template<typename T2 > | |
| void | add_scaled (const TypeTensor< T2 > &, const T) |
| template<typename T2 > | |
| TypeTensor< typename CompareTypes< T, T2 > ::supertype > | operator- (const TypeTensor< T2 > &) const |
| TypeTensor< T > | operator- () const |
| template<typename T2 > | |
| const TypeTensor< T > & | operator-= (const TypeTensor< T2 > &) |
| template<typename T2 > | |
| void | subtract (const TypeTensor< T2 > &) |
| template<typename T2 > | |
| void | subtract_scaled (const TypeTensor< T2 > &, const T) |
| template<typename Scalar > | |
| boostcopy::enable_if_c < ScalarTraits< Scalar > ::value, TypeTensor< typename CompareTypes< T, Scalar > ::supertype > >::type | operator* (const Scalar) const |
| template<typename T2 > | |
| TypeTensor< T > | operator* (const TypeTensor< T2 > &) const |
| template<typename T2 > | |
| TypeVector< typename CompareTypes< T, T2 > ::supertype > | operator* (const TypeVector< T2 > &) const |
| template<typename Scalar > | |
| const TypeTensor< T > & | operator*= (const Scalar) |
| template<typename Scalar > | |
| boostcopy::enable_if_c < ScalarTraits< Scalar > ::value, TypeTensor< typename CompareTypes< T, Scalar > ::supertype > >::type | operator/ (const Scalar) const |
| const TypeTensor< T > & | operator/= (const T) |
| template<typename T2 > | |
| CompareTypes< T, T2 >::supertype | contract (const TypeTensor< T2 > &) const |
| TypeTensor< T > | transpose () const |
| TypeTensor< T > | inverse () const |
| Real | size () const |
| Real | size_sq () const |
| T | det () const |
| T | tr () const |
| void | zero () |
| bool | operator== (const TypeTensor< T > &rhs) const |
| bool | operator< (const TypeTensor< T > &rhs) const |
| template<> | |
| bool | operator< (const TypeTensor< Real > &rhs) const |
| template<> | |
| bool | operator< (const TypeTensor< Complex > &rhs) const |
| bool | operator> (const TypeTensor< T > &rhs) const |
| template<> | |
| bool | operator> (const TypeTensor< Real > &rhs) const |
| template<> | |
| bool | operator> (const TypeTensor< Complex > &rhs) const |
| void | print (std::ostream &os=libMesh::out) const |
| void | write_unformatted (std::ostream &out, const bool newline=true) const |
Public Attributes | |
| T | _coords [LIBMESH_DIM *LIBMESH_DIM] |
Friends | |
| std::ostream & | operator<< (std::ostream &os, const TypeTensor< T > &t) |
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space. The typedef RealTensorValue always defines a real-valued tensor, and NumberTensorValue defines a real or complex-valued tensor depending on how the library was configured.
Definition at line 40 of file tensor_value.h.
| libMesh::TensorValue< T >::TensorValue | ( | ) | [inline] |
Empty constructor. Gives the tensor 0 in LIBMESH_DIM dimensional T space.
Definition at line 157 of file tensor_value.h.
:
TypeTensor<T> ()
{
}
| libMesh::TensorValue< T >::TensorValue | ( | const T | xx, |
| const T | xy = 0, |
||
| const T | xz = 0, |
||
| const T | yx = 0, |
||
| const T | yy = 0, |
||
| const T | yz = 0, |
||
| const T | zx = 0, |
||
| const T | zy = 0, |
||
| const T | zz = 0 |
||
| ) | [inline, explicit] |
Constructor-from-T. By default sets higher dimensional entries to 0.
Definition at line 167 of file tensor_value.h.
:
TypeTensor<T> (xx,xy,xz,yx,yy,yz,zx,zy,zz)
{
}
| libMesh::TensorValue< T >::TensorValue | ( | const Scalar | xx, |
| const Scalar | xy = 0, |
||
| const Scalar | xz = 0, |
||
| const Scalar | yx = 0, |
||
| const Scalar | yy = 0, |
||
| const Scalar | yz = 0, |
||
| const Scalar | zx = 0, |
||
| const Scalar | zy = 0, |
||
| typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, const Scalar >::type | zz = 0 |
||
| ) | [inline, explicit] |
Constructor-from-scalars. By default sets higher dimensional entries to 0.
Definition at line 185 of file tensor_value.h.
:
TypeTensor<T> (xx,xy,xz,yx,yy,yz,zx,zy,zz)
{
}
| libMesh::TensorValue< T >::TensorValue | ( | const TypeVector< T2 > & | vx | ) | [inline] |
Constructor. Takes 1 row vector for LIBMESH_DIM=1
Definition at line 215 of file tensor_value.h.
:
TypeTensor<T> (vx)
{
}
| libMesh::TensorValue< T >::TensorValue | ( | const TypeVector< T2 > & | vx, |
| const TypeVector< T2 > & | vy | ||
| ) | [inline] |
Constructor. Takes 2 row vectors for LIBMESH_DIM=2
Definition at line 225 of file tensor_value.h.
:
TypeTensor<T> (vx, vy)
{
}
| libMesh::TensorValue< T >::TensorValue | ( | const TypeVector< T2 > & | vx, |
| const TypeVector< T2 > & | vy, | ||
| const TypeVector< T2 > & | vz | ||
| ) | [inline] |
Constructor. Takes 3 row vectors for LIBMESH_DIM=3
Definition at line 236 of file tensor_value.h.
:
TypeTensor<T> (vx, vy, vz)
{
}
| libMesh::TensorValue< T >::TensorValue | ( | const TensorValue< T2 > & | p | ) | [inline] |
| libMesh::TensorValue< T >::TensorValue | ( | const TypeTensor< T2 > & | p | ) | [inline] |
| libMesh::TensorValue< T >::TensorValue | ( | const TypeTensor< Real > & | p_re, |
| const TypeTensor< Real > & | p_im | ||
| ) | [inline] |
Constructor that takes two TypeTensor<Real> representing the real and imaginary part as arguments.
Definition at line 257 of file tensor_value.h.
: TypeTensor<T> (Complex (p_re(0,0), p_im(0,0)), Complex (p_re(0,1), p_im(0,1)), Complex (p_re(0,2), p_im(0,2)), Complex (p_re(1,0), p_im(1,0)), Complex (p_re(1,1), p_im(1,1)), Complex (p_re(1,2), p_im(1,2)), Complex (p_re(2,0), p_im(2,0)), Complex (p_re(2,1), p_im(2,1)), Complex (p_re(2,2), p_im(2,2))) { }
| void libMesh::TypeTensor< T >::add | ( | const TypeTensor< T2 > & | p | ) | [inline, inherited] |
Add to this tensor without creating a temporary.
Definition at line 706 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
{
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] += p._coords[i];
}
| void libMesh::TypeTensor< T >::add_scaled | ( | const TypeTensor< T2 > & | p, |
| const T | factor | ||
| ) | [inline, inherited] |
Add a scaled tensor to this tensor without creating a temporary.
Definition at line 717 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
Referenced by libMesh::HPCoarsenTest::add_projection(), libMesh::System::calculate_norm(), libMesh::MeshFunction::hessian(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::System::point_hessian(), and libMesh::HPCoarsenTest::select_refinement().
{
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] += factor*p._coords[i];
}
| void libMesh::TypeTensor< T >::assign | ( | const TypeTensor< T2 > & | p | ) | [inline, inherited] |
Assign to a tensor without creating a temporary.
Definition at line 577 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
{
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] = p._coords[i];
}
| CompareTypes< T, T2 >::supertype libMesh::TypeTensor< T >::contract | ( | const TypeTensor< T2 > & | t | ) | const [inline, inherited] |
Multiply 2 tensors together, i.e. dyadic product sum_ij Aij*Bij. The tensors may be of different types.
Multiply 2 tensors together, i.e. sum Aij*Bij. The tensors may be of different types.
Definition at line 1059 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords, and libMesh::Parallel::sum().
Referenced by libMesh::HPCoarsenTest::add_projection(), libMesh::TensorTools::inner_product(), and libMesh::HPCoarsenTest::select_refinement().
| T libMesh::TypeTensor< T >::det | ( | ) | const [inline, inherited] |
Returns the determinant of the tensor. Because these are 3x3 tensors at most, we don't do an LU decomposition like DenseMatrix does.
Definition at line 1080 of file type_tensor.h.
Referenced by libMesh::Sphere::Sphere().
{
#if LIBMESH_DIM == 1
return _coords[0];
#endif
#if LIBMESH_DIM == 2
return (_coords[0] * _coords[3]
- _coords[1] * _coords[2]);
#endif
#if LIBMESH_DIM == 3
return (_coords[0] * _coords[4] * _coords[8]
+ _coords[1] * _coords[5] * _coords[6]
+ _coords[2] * _coords[3] * _coords[7]
- _coords[0] * _coords[5] * _coords[7]
- _coords[1] * _coords[3] * _coords[8]
- _coords[2] * _coords[4] * _coords[6]);
#endif
}
| TypeTensor< T > libMesh::TypeTensor< T >::inverse | ( | ) | const [inline, inherited] |
Returns the inverse of the tensor.
Definition at line 962 of file type_tensor.h.
{
#if LIBMESH_DIM == 1
libmesh_assert_not_equal_to(_coords[0], static_cast<T>(0.));
return TypeTensor(1. / _coords[0]);
#endif
#if LIBMESH_DIM == 2
// Get convenient reference to this.
const TypeTensor<T> & A = *this;
// Use temporary variables, avoid multiple accesses
T a = A(0,0), b = A(0,1),
c = A(1,0), d = A(1,1);
// Make sure det = ad - bc is not zero
T my_det = a*d - b*c;
libmesh_assert_not_equal_to(my_det, static_cast<T>(0.));
return TypeTensor(d/my_det, -b/my_det, -c/my_det, a/my_det);
#endif
#if LIBMESH_DIM == 3
// Get convenient reference to this.
const TypeTensor<T> & A = *this;
T a11 = A(0,0), a12 = A(0,1), a13 = A(0,2),
a21 = A(1,0), a22 = A(1,1), a23 = A(1,2),
a31 = A(2,0), a32 = A(2,1), a33 = A(2,2);
T my_det = a11*(a33*a22-a32*a23) - a21*(a33*a12-a32*a13) + a31*(a23*a12-a22*a13);
libmesh_assert_not_equal_to(my_det, static_cast<T>(0.));
// Inline comment characters are for lining up columns.
return TypeTensor( (a33*a22-a32*a23)/my_det, -(a33*a12-a32*a13)/my_det, (a23*a12-a22*a13)/my_det,
-(a33*a21-a31*a23)/my_det, (a33*a11-a31*a13)/my_det, -(a23*a11-a21*a13)/my_det,
(a32*a21-a31*a22)/my_det, -(a32*a11-a31*a12)/my_det, (a22*a11-a21*a12)/my_det);
#endif
}
| const T & libMesh::TypeTensor< T >::operator() | ( | const unsigned int | i, |
| const unsigned int | j | ||
| ) | const [inline, inherited] |
Return the
element of the tensor.
Definition at line 587 of file type_tensor.h.
{
libmesh_assert_less (i, 3);
libmesh_assert_less (j, 3);
#if LIBMESH_DIM < 3
const static T my_zero = 0;
if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
return my_zero;
#endif
return _coords[i*LIBMESH_DIM+j];
}
| T & libMesh::TypeTensor< T >::operator() | ( | const unsigned int | i, |
| const unsigned int | j | ||
| ) | [inline, inherited] |
Return a writeable reference to the
element of the tensor.
Definition at line 606 of file type_tensor.h.
{
#if LIBMESH_DIM < 3
if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
libmesh_error_msg("ERROR: You are assigning to a tensor component that is out of range for the compiled LIBMESH_DIM!");
#endif
libmesh_assert_less (i, LIBMESH_DIM);
libmesh_assert_less (j, LIBMESH_DIM);
return _coords[i*LIBMESH_DIM+j];
}
| boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< typename CompareTypes< T, Scalar >::supertype > >::type libMesh::TypeTensor< T >::operator* | ( | const Scalar | factor | ) | const [inline, inherited] |
Multiply a tensor by a number, i.e. scale.
Definition at line 833 of file type_tensor.h.
{
typedef typename CompareTypes<T, Scalar>::supertype TS;
#if LIBMESH_DIM == 1
return TypeTensor<TS>(_coords[0]*factor);
#endif
#if LIBMESH_DIM == 2
return TypeTensor<TS>(_coords[0]*factor,
_coords[1]*factor,
_coords[2]*factor,
_coords[3]*factor);
#endif
#if LIBMESH_DIM == 3
return TypeTensor<TS>(_coords[0]*factor,
_coords[1]*factor,
_coords[2]*factor,
_coords[3]*factor,
_coords[4]*factor,
_coords[5]*factor,
_coords[6]*factor,
_coords[7]*factor,
_coords[8]*factor);
#endif
}
| TypeTensor< T > libMesh::TypeTensor< T >::operator* | ( | const TypeTensor< T2 > & | p | ) | const [inline, inherited] |
Multiply 2 tensors together, i.e. matrix product. The tensors may be of different types.
Definition at line 1038 of file type_tensor.h.
{
TypeTensor<T> returnval;
for (unsigned int i=0; i<LIBMESH_DIM; i++)
for (unsigned int j=0; j<LIBMESH_DIM; j++)
for (unsigned int k=0; k<LIBMESH_DIM; k++)
returnval(i,j) += (*this)(i,k)*p(k,j);
return returnval;
}
| TypeVector< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::operator* | ( | const TypeVector< T2 > & | p | ) | const [inline, inherited] |
Multiply a tensor and vector together, i.e. matrix-vector product. The tensor and vector may be of different types.
Definition at line 1023 of file type_tensor.h.
{
TypeVector<typename CompareTypes<T,T2>::supertype> returnval;
for (unsigned int i=0; i<LIBMESH_DIM; i++)
for (unsigned int j=0; j<LIBMESH_DIM; j++)
returnval(i) += (*this)(i,j)*p(j);
return returnval;
}
| const TypeTensor< T > & libMesh::TypeTensor< T >::operator*= | ( | const Scalar | factor | ) | [inline, inherited] |
Multiply this tensor by a number, i.e. scale.
Definition at line 880 of file type_tensor.h.
{
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] *= factor;
return *this;
}
| TypeTensor< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::operator+ | ( | const TypeTensor< T2 > & | p | ) | const [inline, inherited] |
Add two tensors.
Definition at line 660 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
{
#if LIBMESH_DIM == 1
return TypeTensor(_coords[0] + p._coords[0]);
#endif
#if LIBMESH_DIM == 2
return TypeTensor(_coords[0] + p._coords[0],
_coords[1] + p._coords[1],
0.,
_coords[2] + p._coords[2],
_coords[3] + p._coords[3]);
#endif
#if LIBMESH_DIM == 3
return TypeTensor(_coords[0] + p._coords[0],
_coords[1] + p._coords[1],
_coords[2] + p._coords[2],
_coords[3] + p._coords[3],
_coords[4] + p._coords[4],
_coords[5] + p._coords[5],
_coords[6] + p._coords[6],
_coords[7] + p._coords[7],
_coords[8] + p._coords[8]);
#endif
}
| const TypeTensor< T > & libMesh::TypeTensor< T >::operator+= | ( | const TypeTensor< T2 > & | p | ) | [inline, inherited] |
| TypeTensor< typename CompareTypes< T, T2 >::supertype > libMesh::TypeTensor< T >::operator- | ( | const TypeTensor< T2 > & | p | ) | const [inline, inherited] |
Subtract two tensors.
Definition at line 730 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
{
#if LIBMESH_DIM == 1
return TypeTensor(_coords[0] - p._coords[0]);
#endif
#if LIBMESH_DIM == 2
return TypeTensor(_coords[0] - p._coords[0],
_coords[1] - p._coords[1],
0.,
_coords[2] - p._coords[2],
_coords[3] - p._coords[3]);
#endif
#if LIBMESH_DIM == 3
return TypeTensor(_coords[0] - p._coords[0],
_coords[1] - p._coords[1],
_coords[2] - p._coords[2],
_coords[3] - p._coords[3],
_coords[4] - p._coords[4],
_coords[5] - p._coords[5],
_coords[6] - p._coords[6],
_coords[7] - p._coords[7],
_coords[8] - p._coords[8]);
#endif
}
| TypeTensor< T > libMesh::TypeTensor< T >::operator- | ( | ) | const [inline, inherited] |
Return the opposite of a tensor
Definition at line 797 of file type_tensor.h.
{
#if LIBMESH_DIM == 1
return TypeTensor(-_coords[0]);
#endif
#if LIBMESH_DIM == 2
return TypeTensor(-_coords[0],
-_coords[1],
-_coords[2],
-_coords[3]);
#endif
#if LIBMESH_DIM == 3
return TypeTensor(-_coords[0],
-_coords[1],
-_coords[2],
-_coords[3],
-_coords[4],
-_coords[5],
-_coords[6],
-_coords[7],
-_coords[8]);
#endif
}
| const TypeTensor< T > & libMesh::TypeTensor< T >::operator-= | ( | const TypeTensor< T2 > & | p | ) | [inline, inherited] |
Subtract from this tensor.
Definition at line 764 of file type_tensor.h.
{
this->subtract (p);
return *this;
}
| boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< typename CompareTypes< T, Scalar >::supertype > >::type libMesh::TypeTensor< T >::operator/ | ( | const Scalar | factor | ) | const [inline, inherited] |
Divide a tensor by a number, i.e. scale.
Definition at line 897 of file type_tensor.h.
{
libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
typedef typename CompareTypes<T, Scalar>::supertype TS;
#if LIBMESH_DIM == 1
return TypeTensor<TS>(_coords[0]/factor);
#endif
#if LIBMESH_DIM == 2
return TypeTensor<TS>(_coords[0]/factor,
_coords[1]/factor,
_coords[2]/factor,
_coords[3]/factor);
#endif
#if LIBMESH_DIM == 3
return TypeTensor<TS>(_coords[0]/factor,
_coords[1]/factor,
_coords[2]/factor,
_coords[3]/factor,
_coords[4]/factor,
_coords[5]/factor,
_coords[6]/factor,
_coords[7]/factor,
_coords[8]/factor);
#endif
}
| const TypeTensor< T > & libMesh::TypeTensor< T >::operator/= | ( | const T | factor | ) | [inline, inherited] |
Divide this tensor by a number, i.e. scale.
Definition at line 1006 of file type_tensor.h.
{
libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] /= factor;
return *this;
}
| bool libMesh::TypeTensor< Real >::operator< | ( | const TypeTensor< Real > & | rhs | ) | const [inherited] |
Definition at line 113 of file type_tensor.C.
{
for (unsigned int i=0; i<LIBMESH_DIM; i++)
for (unsigned int j=0; j<LIBMESH_DIM; j++)
{
if ((*this)(i,j) < rhs(i,j))
return true;
if ((*this)(i,j) > rhs(i,j))
return false;
}
return false;
}
| bool libMesh::TypeTensor< Complex >::operator< | ( | const TypeTensor< Complex > & | rhs | ) | const [inherited] |
Definition at line 146 of file type_tensor.C.
{
for (unsigned int i=0; i<LIBMESH_DIM; i++)
for (unsigned int j=0; j<LIBMESH_DIM; j++)
{
if ((*this)(i,j).real() < rhs(i,j).real())
return true;
if ((*this)(i,j).real() > rhs(i,j).real())
return false;
if ((*this)(i,j).imag() < rhs(i,j).imag())
return true;
if ((*this)(i,j).imag() > rhs(i,j).imag())
return false;
}
return false;
}
| bool libMesh::TypeTensor< T >::operator< | ( | const TypeTensor< T > & | rhs | ) | const [inherited] |
true if this tensor is "less" than another. Useful for sorting. | boostcopy::enable_if_c< ScalarTraits<Scalar>::value, TensorValue&>::type libMesh::TensorValue< T >::operator= | ( | const Scalar & | libmesh_dbg_varp | ) | [inline] |
Assignment-from-scalar operator. Used only to zero out tensors.
Definition at line 131 of file tensor_value.h.
References libMesh::TypeTensor< T >::zero().
{ libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; }
| bool libMesh::TypeTensor< T >::operator== | ( | const TypeTensor< T > & | rhs | ) | const [inline, inherited] |
true if two tensors are equal valued. Definition at line 1142 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords, std::abs(), and libMesh::TOLERANCE.
{
#if LIBMESH_DIM == 1
return (std::abs(_coords[0] - rhs._coords[0])
< TOLERANCE);
#endif
#if LIBMESH_DIM == 2
return ((std::abs(_coords[0] - rhs._coords[0]) +
std::abs(_coords[1] - rhs._coords[1]) +
std::abs(_coords[2] - rhs._coords[2]) +
std::abs(_coords[3] - rhs._coords[3]))
< 4.*TOLERANCE);
#endif
#if LIBMESH_DIM == 3
return ((std::abs(_coords[0] - rhs._coords[0]) +
std::abs(_coords[1] - rhs._coords[1]) +
std::abs(_coords[2] - rhs._coords[2]) +
std::abs(_coords[3] - rhs._coords[3]) +
std::abs(_coords[4] - rhs._coords[4]) +
std::abs(_coords[5] - rhs._coords[5]) +
std::abs(_coords[6] - rhs._coords[6]) +
std::abs(_coords[7] - rhs._coords[7]) +
std::abs(_coords[8] - rhs._coords[8]))
< 9.*TOLERANCE);
#endif
}
| bool libMesh::TypeTensor< Real >::operator> | ( | const TypeTensor< Real > & | rhs | ) | const [inherited] |
Definition at line 129 of file type_tensor.C.
{
for (unsigned int i=0; i<LIBMESH_DIM; i++)
for (unsigned int j=0; j<LIBMESH_DIM; j++)
{
if ((*this)(i,j) > rhs(i,j))
return true;
if ((*this)(i,j) < rhs(i,j))
return false;
}
return false;
}
| bool libMesh::TypeTensor< Complex >::operator> | ( | const TypeTensor< Complex > & | rhs | ) | const [inherited] |
Definition at line 166 of file type_tensor.C.
{
for (unsigned int i=0; i<LIBMESH_DIM; i++)
for (unsigned int j=0; j<LIBMESH_DIM; j++)
{
if ((*this)(i,j).real() > rhs(i,j).real())
return true;
if ((*this)(i,j).real() < rhs(i,j).real())
return false;
if ((*this)(i,j).imag() > rhs(i,j).imag())
return true;
if ((*this)(i,j).imag() < rhs(i,j).imag())
return false;
}
return false;
}
| bool libMesh::TypeTensor< T >::operator> | ( | const TypeTensor< T > & | rhs | ) | const [inherited] |
true if this tensor is "greater" than another. | void libMesh::TypeTensor< T >::print | ( | std::ostream & | os = libMesh::out | ) | const [inherited] |
Formatted print, by default to libMesh::out.
Definition at line 39 of file type_tensor.C.
{
#if LIBMESH_DIM == 1
os << "x=" << (*this)(0,0) << std::endl;
#endif
#if LIBMESH_DIM == 2
os << "(xx,xy)=("
<< std::setw(8) << (*this)(0,0) << ", "
<< std::setw(8) << (*this)(0,1) << ")"
<< std::endl;
os << "(yx,yy)=("
<< std::setw(8) << (*this)(1,0) << ", "
<< std::setw(8) << (*this)(1,1) << ")"
<< std::endl;
#endif
#if LIBMESH_DIM == 3
os << "(xx,xy,xz)=("
<< std::setw(8) << (*this)(0,0) << ", "
<< std::setw(8) << (*this)(0,1) << ", "
<< std::setw(8) << (*this)(0,2) << ")"
<< std::endl;
os << "(yx,yy,yz)=("
<< std::setw(8) << (*this)(1,0) << ", "
<< std::setw(8) << (*this)(1,1) << ", "
<< std::setw(8) << (*this)(1,2) << ")"
<< std::endl;
os << "(zx,zy,zz)=("
<< std::setw(8) << (*this)(2,0) << ", "
<< std::setw(8) << (*this)(2,1) << ", "
<< std::setw(8) << (*this)(2,2) << ")"
<< std::endl;
#endif
}
| TypeVector< T > libMesh::TypeTensor< T >::row | ( | const unsigned int | r | ) | [inline, inherited] |
Return one row of the tensor as a TypeVector.
Definition at line 646 of file type_tensor.h.
References libMesh::TypeVector< T >::_coords.
{
TypeVector<T> return_vector;
for(unsigned int j=0; j<LIBMESH_DIM; j++)
return_vector._coords[j] = _coords[r*LIBMESH_DIM + j];
return return_vector;
}
| Real libMesh::TypeTensor< T >::size | ( | ) | const [inline, inherited] |
Returns the Frobenius norm of the tensor, i.e. the square-root of the sum of the elements squared.
Definition at line 1071 of file type_tensor.h.
Referenced by libMesh::System::calculate_norm().
{
return std::sqrt(this->size_sq());
}
| Real libMesh::TypeTensor< T >::size_sq | ( | ) | const [inline, inherited] |
Returns the Frobenius norm of the tensor squared, i.e. sum of the element magnitudes squared.
Definition at line 1130 of file type_tensor.h.
References libMesh::TensorTools::norm_sq(), libMesh::Real, and libMesh::Parallel::sum().
Referenced by libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::System::calculate_norm(), libMesh::ExactErrorEstimator::find_squared_element_error(), and libMesh::HPCoarsenTest::select_refinement().
{
Real sum = 0.;
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
sum += TensorTools::norm_sq(_coords[i]);
return sum;
}
| ConstTypeTensorColumn< T > libMesh::TypeTensor< T >::slice | ( | const unsigned int | i | ) | const [inline, inherited] |
Return a proxy for the
column of the tensor.
Definition at line 626 of file type_tensor.h.
{
libmesh_assert_less (i, LIBMESH_DIM);
return ConstTypeTensorColumn<T>(*this, i);
}
| TypeTensorColumn< T > libMesh::TypeTensor< T >::slice | ( | const unsigned int | i | ) | [inline, inherited] |
Return a writeable proxy for the
column of the tensor.
Definition at line 636 of file type_tensor.h.
{
libmesh_assert_less (i, LIBMESH_DIM);
return TypeTensorColumn<T>(*this, i);
}
| void libMesh::TypeTensor< T >::subtract | ( | const TypeTensor< T2 > & | p | ) | [inline, inherited] |
Subtract from this tensor without creating a temporary.
Definition at line 776 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
{
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] -= p._coords[i];
}
| void libMesh::TypeTensor< T >::subtract_scaled | ( | const TypeTensor< T2 > & | p, |
| const T | factor | ||
| ) | [inline, inherited] |
Subtract a scaled value from this tensor without creating a temporary.
Definition at line 787 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
Referenced by libMesh::HPCoarsenTest::select_refinement().
{
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] -= factor*p._coords[i];
}
| T libMesh::TypeTensor< T >::tr | ( | ) | const [inline, inherited] |
Returns the trace of the tensor.
Definition at line 1103 of file type_tensor.h.
| TypeTensor< T > libMesh::TypeTensor< T >::transpose | ( | ) | const [inline, inherited] |
The transpose (with complex numbers not conjugated) of the tensor.
Definition at line 932 of file type_tensor.h.
{
#if LIBMESH_DIM == 1
return TypeTensor(_coords[0]);
#endif
#if LIBMESH_DIM == 2
return TypeTensor(_coords[0],
_coords[2],
_coords[1],
_coords[3]);
#endif
#if LIBMESH_DIM == 3
return TypeTensor(_coords[0],
_coords[3],
_coords[6],
_coords[1],
_coords[4],
_coords[7],
_coords[2],
_coords[5],
_coords[8]);
#endif
}
| void libMesh::TypeTensor< T >::write_unformatted | ( | std::ostream & | out, |
| const bool | newline = true |
||
| ) | const [inherited] |
Unformatted print to the stream out. Simply prints the elements of the tensor separated by spaces and newlines.
Definition at line 83 of file type_tensor.C.
References libMesh::libmesh_assert().
{
libmesh_assert (out_stream);
out_stream << std::setiosflags(std::ios::showpoint)
<< (*this)(0,0) << " "
<< (*this)(0,1) << " "
<< (*this)(0,2) << " ";
if (newline)
out_stream << '\n';
out_stream << std::setiosflags(std::ios::showpoint)
<< (*this)(1,0) << " "
<< (*this)(1,1) << " "
<< (*this)(1,2) << " ";
if (newline)
out_stream << '\n';
out_stream << std::setiosflags(std::ios::showpoint)
<< (*this)(2,0) << " "
<< (*this)(2,1) << " "
<< (*this)(2,2) << " ";
if (newline)
out_stream << '\n';
}
| void libMesh::TypeTensor< T >::zero | ( | ) | [inline, inherited] |
Zero the tensor in any dimension.
Definition at line 1120 of file type_tensor.h.
Referenced by libMesh::TensorValue< T >::operator=(), and libMesh::TypeTensor< T >::operator=().
{
for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
_coords[i] = 0.;
}
| std::ostream& operator<< | ( | std::ostream & | os, |
| const TypeTensor< T > & | t | ||
| ) | [friend, inherited] |
Formatted print as above but allows you to do std::cout << t << std::endl;
Definition at line 341 of file type_tensor.h.
{
t.print(os);
return os;
}
T libMesh::TypeTensor< T >::_coords[LIBMESH_DIM *LIBMESH_DIM] [inherited] |
The coordinates of the TypeTensor
Definition at line 358 of file type_tensor.h.
Referenced by libMesh::TypeTensor< T >::add(), libMesh::TypeTensor< T >::add_scaled(), libMesh::TypeTensor< T >::assign(), libMesh::TypeTensor< T >::contract(), libMesh::TypeTensor< T >::operator+(), libMesh::TypeTensor< T >::operator-(), libMesh::TypeTensor< T >::operator==(), libMesh::TypeTensor< T >::subtract(), libMesh::TypeTensor< T >::subtract_scaled(), and libMesh::TypeTensor< T >::TypeTensor().