libMesh::TensorValue< T > Class Template Reference
#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 &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 |
| 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) |
Detailed Description
template<typename T>
class libMesh::TensorValue< 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.
Constructor & Destructor Documentation
| 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.
| 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.
| 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.
| 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.
| 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.
| 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.
| libMesh::TensorValue< T >::TensorValue | ( | const TensorValue< T2 > & | p | ) | [inline] |
Copy-constructor.
Definition at line 205 of file tensor_value.h.
| libMesh::TensorValue< T >::TensorValue | ( | const TypeTensor< T2 > & | p | ) | [inline] |
Copy-constructor.
Definition at line 248 of file tensor_value.h.
| 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.
00258 : 00259 TypeTensor<T> (Complex (p_re(0,0), p_im(0,0)), 00260 Complex (p_re(0,1), p_im(0,1)), 00261 Complex (p_re(0,2), p_im(0,2)), 00262 Complex (p_re(1,0), p_im(1,0)), 00263 Complex (p_re(1,1), p_im(1,1)), 00264 Complex (p_re(1,2), p_im(1,2)), 00265 Complex (p_re(2,0), p_im(2,0)), 00266 Complex (p_re(2,1), p_im(2,1)), 00267 Complex (p_re(2,2), p_im(2,2))) 00268 { 00269 } #endif
Member Function Documentation
| 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.
Referenced by libMesh::TypeTensor< T >::operator+=().
00707 { 00708 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00709 _coords[i] += p._coords[i]; 00710 }
| 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().
00718 { 00719 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00720 _coords[i] += factor*p._coords[i]; 00721 00722 }
| void libMesh::TypeTensor< T >::assign | ( | const TypeTensor< T2 > & | p | ) | [inline, inherited] |
Assign to a tensor without creating a temporary.
Definition at line 572 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
00573 { 00574 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00575 _coords[i] = p._coords[i]; 00576 }
| 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 1015 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 1036 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
Referenced by libMesh::Sphere::Sphere().
01037 { 01038 #if LIBMESH_DIM == 1 01039 return _coords[0]; 01040 #endif 01041 01042 #if LIBMESH_DIM == 2 01043 return (_coords[0] * _coords[3] 01044 - _coords[1] * _coords[2]); 01045 #endif 01046 01047 #if LIBMESH_DIM == 3 01048 return (_coords[0] * _coords[4] * _coords[8] 01049 + _coords[1] * _coords[5] * _coords[6] 01050 + _coords[2] * _coords[3] * _coords[7] 01051 - _coords[0] * _coords[5] * _coords[7] 01052 - _coords[1] * _coords[3] * _coords[8] 01053 - _coords[2] * _coords[4] * _coords[6]); 01054 #endif 01055 }
| 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 601 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
00603 { 00604 #if LIBMESH_DIM < 3 00605 00606 if (i >= LIBMESH_DIM || j >= LIBMESH_DIM) 00607 { 00608 // libMesh::err << "ERROR: You are assigning to a tensor component" << std::endl 00609 // << "that is out of range for the compiled LIBMESH_DIM!" << std::endl 00610 // << " LIBMESH_DIM=" << LIBMESH_DIM << " , i=" << i << " , j=" << j << std::endl; 00611 libmesh_error(); 00612 } 00613 00614 #endif 00615 00616 libmesh_assert_less (i, LIBMESH_DIM); 00617 libmesh_assert_less (j, LIBMESH_DIM); 00618 00619 return _coords[i*LIBMESH_DIM+j]; 00620 }
| 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 582 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
00584 { 00585 libmesh_assert_less (i, 3); 00586 libmesh_assert_less (j, 3); 00587 00588 #if LIBMESH_DIM < 3 00589 const static T my_zero = 0; 00590 if (i >= LIBMESH_DIM || j >= LIBMESH_DIM) 00591 return my_zero; 00592 #endif 00593 00594 return _coords[i*LIBMESH_DIM+j]; 00595 }
| 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 979 of file type_tensor.h.
| 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 994 of file type_tensor.h.
| 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.
References libMesh::TypeTensor< T >::_coords.
00834 { 00835 typedef typename CompareTypes<T, Scalar>::supertype TS; 00836 00837 00838 #if LIBMESH_DIM == 1 00839 return TypeTensor<TS>(_coords[0]*factor); 00840 #endif 00841 00842 #if LIBMESH_DIM == 2 00843 return TypeTensor<TS>(_coords[0]*factor, 00844 _coords[1]*factor, 00845 _coords[2]*factor, 00846 _coords[3]*factor); 00847 #endif 00848 00849 #if LIBMESH_DIM == 3 00850 return TypeTensor<TS>(_coords[0]*factor, 00851 _coords[1]*factor, 00852 _coords[2]*factor, 00853 _coords[3]*factor, 00854 _coords[4]*factor, 00855 _coords[5]*factor, 00856 _coords[6]*factor, 00857 _coords[7]*factor, 00858 _coords[8]*factor); 00859 #endif 00860 }
| 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.
References libMesh::TypeTensor< T >::_coords.
00881 { 00882 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00883 _coords[i] *= factor; 00884 00885 return *this; 00886 }
| 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, and libMesh::TypeTensor< T >::TypeTensor().
00661 { 00662 00663 #if LIBMESH_DIM == 1 00664 return TypeTensor(_coords[0] + p._coords[0]); 00665 #endif 00666 00667 #if LIBMESH_DIM == 2 00668 return TypeTensor(_coords[0] + p._coords[0], 00669 _coords[1] + p._coords[1], 00670 0., 00671 _coords[2] + p._coords[2], 00672 _coords[3] + p._coords[3]); 00673 #endif 00674 00675 #if LIBMESH_DIM == 3 00676 return TypeTensor(_coords[0] + p._coords[0], 00677 _coords[1] + p._coords[1], 00678 _coords[2] + p._coords[2], 00679 _coords[3] + p._coords[3], 00680 _coords[4] + p._coords[4], 00681 _coords[5] + p._coords[5], 00682 _coords[6] + p._coords[6], 00683 _coords[7] + p._coords[7], 00684 _coords[8] + p._coords[8]); 00685 #endif 00686 00687 }
| const TypeTensor< T > & libMesh::TypeTensor< T >::operator+= | ( | const TypeTensor< T2 > & | p | ) | [inline, inherited] |
Add to this tensor.
Definition at line 694 of file type_tensor.h.
References libMesh::TypeTensor< T >::add().
00695 { 00696 this->add (p); 00697 00698 return *this; 00699 }
| TypeTensor< T > libMesh::TypeTensor< T >::operator- | ( | ) | const [inline, inherited] |
Return the opposite of a tensor
Definition at line 797 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords, and libMesh::TypeTensor< T >::TypeTensor().
00798 { 00799 00800 #if LIBMESH_DIM == 1 00801 return TypeTensor(-_coords[0]); 00802 #endif 00803 00804 #if LIBMESH_DIM == 2 00805 return TypeTensor(-_coords[0], 00806 -_coords[1], 00807 -_coords[2], 00808 -_coords[3]); 00809 #endif 00810 00811 #if LIBMESH_DIM == 3 00812 return TypeTensor(-_coords[0], 00813 -_coords[1], 00814 -_coords[2], 00815 -_coords[3], 00816 -_coords[4], 00817 -_coords[5], 00818 -_coords[6], 00819 -_coords[7], 00820 -_coords[8]); 00821 #endif 00822 00823 }
| 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, and libMesh::TypeTensor< T >::TypeTensor().
00731 { 00732 00733 #if LIBMESH_DIM == 1 00734 return TypeTensor(_coords[0] - p._coords[0]); 00735 #endif 00736 00737 #if LIBMESH_DIM == 2 00738 return TypeTensor(_coords[0] - p._coords[0], 00739 _coords[1] - p._coords[1], 00740 0., 00741 _coords[2] - p._coords[2], 00742 _coords[3] - p._coords[3]); 00743 #endif 00744 00745 #if LIBMESH_DIM == 3 00746 return TypeTensor(_coords[0] - p._coords[0], 00747 _coords[1] - p._coords[1], 00748 _coords[2] - p._coords[2], 00749 _coords[3] - p._coords[3], 00750 _coords[4] - p._coords[4], 00751 _coords[5] - p._coords[5], 00752 _coords[6] - p._coords[6], 00753 _coords[7] - p._coords[7], 00754 _coords[8] - p._coords[8]); 00755 #endif 00756 00757 }
| 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.
References libMesh::TypeTensor< T >::subtract().
00765 { 00766 this->subtract (p); 00767 00768 return *this; 00769 }
| 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.
References libMesh::TypeTensor< T >::_coords.
00898 { 00899 libmesh_assert_not_equal_to (factor, static_cast<T>(0.)); 00900 00901 typedef typename CompareTypes<T, Scalar>::supertype TS; 00902 00903 #if LIBMESH_DIM == 1 00904 return TypeTensor<TS>(_coords[0]/factor); 00905 #endif 00906 00907 #if LIBMESH_DIM == 2 00908 return TypeTensor<TS>(_coords[0]/factor, 00909 _coords[1]/factor, 00910 _coords[2]/factor, 00911 _coords[3]/factor); 00912 #endif 00913 00914 #if LIBMESH_DIM == 3 00915 return TypeTensor<TS>(_coords[0]/factor, 00916 _coords[1]/factor, 00917 _coords[2]/factor, 00918 _coords[3]/factor, 00919 _coords[4]/factor, 00920 _coords[5]/factor, 00921 _coords[6]/factor, 00922 _coords[7]/factor, 00923 _coords[8]/factor); 00924 #endif 00925 00926 }
| const TypeTensor< T > & libMesh::TypeTensor< T >::operator/= | ( | const T | factor | ) | [inline, inherited] |
Divide this tensor by a number, i.e. scale.
Definition at line 962 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
00963 { 00964 libmesh_assert_not_equal_to (factor, static_cast<T>(0.)); 00965 00966 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00967 _coords[i] /= factor; 00968 00969 return *this; 00970 }
| bool libMesh::TypeTensor< Complex >::operator< | ( | const TypeTensor< Complex > & | rhs | ) | const [inline, inherited] |
Definition at line 146 of file type_tensor.C.
00147 { 00148 for (unsigned int i=0; i<LIBMESH_DIM; i++) 00149 for (unsigned int j=0; j<LIBMESH_DIM; j++) 00150 { 00151 if ((*this)(i,j).real() < rhs(i,j).real()) 00152 return true; 00153 if ((*this)(i,j).real() > rhs(i,j).real()) 00154 return false; 00155 if ((*this)(i,j).imag() < rhs(i,j).imag()) 00156 return true; 00157 if ((*this)(i,j).imag() > rhs(i,j).imag()) 00158 return false; 00159 } 00160 return false; 00161 }
| bool libMesh::TypeTensor< Real >::operator< | ( | const TypeTensor< Real > & | rhs | ) | const [inline, inherited] |
Definition at line 113 of file type_tensor.C.
| bool libMesh::TypeTensor< T >::operator< | ( | const TypeTensor< T > & | rhs | ) | const [inherited] |
- Returns:
trueif this tensor is "less" than another. Useful for sorting.
| boostcopy::enable_if_c< ScalarTraits<Scalar>::value, TensorValue&>::type libMesh::TensorValue< T >::operator= | ( | const Scalar & | p | ) | [inline] |
Assignment-from-scalar operator. Used only to zero out tensors.
Reimplemented from libMesh::TypeTensor< T >.
Definition at line 131 of file tensor_value.h.
References libMesh::TypeTensor< T >::zero().
00132 { libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; }
| bool libMesh::TypeTensor< T >::operator== | ( | const TypeTensor< T > & | rhs | ) | const [inline, inherited] |
- Returns:
trueif two tensors are equal valued.
Definition at line 1098 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords, std::abs(), and libMesh::TOLERANCE.
01099 { 01100 #if LIBMESH_DIM == 1 01101 return (std::abs(_coords[0] - rhs._coords[0]) 01102 < TOLERANCE); 01103 #endif 01104 01105 #if LIBMESH_DIM == 2 01106 return ((std::abs(_coords[0] - rhs._coords[0]) + 01107 std::abs(_coords[1] - rhs._coords[1]) + 01108 std::abs(_coords[2] - rhs._coords[2]) + 01109 std::abs(_coords[3] - rhs._coords[3])) 01110 < 4.*TOLERANCE); 01111 #endif 01112 01113 #if LIBMESH_DIM == 3 01114 return ((std::abs(_coords[0] - rhs._coords[0]) + 01115 std::abs(_coords[1] - rhs._coords[1]) + 01116 std::abs(_coords[2] - rhs._coords[2]) + 01117 std::abs(_coords[3] - rhs._coords[3]) + 01118 std::abs(_coords[4] - rhs._coords[4]) + 01119 std::abs(_coords[5] - rhs._coords[5]) + 01120 std::abs(_coords[6] - rhs._coords[6]) + 01121 std::abs(_coords[7] - rhs._coords[7]) + 01122 std::abs(_coords[8] - rhs._coords[8])) 01123 < 9.*TOLERANCE); 01124 #endif 01125 01126 }
| bool libMesh::TypeTensor< Complex >::operator> | ( | const TypeTensor< Complex > & | rhs | ) | const [inline, inherited] |
Definition at line 166 of file type_tensor.C.
00167 { 00168 for (unsigned int i=0; i<LIBMESH_DIM; i++) 00169 for (unsigned int j=0; j<LIBMESH_DIM; j++) 00170 { 00171 if ((*this)(i,j).real() > rhs(i,j).real()) 00172 return true; 00173 if ((*this)(i,j).real() < rhs(i,j).real()) 00174 return false; 00175 if ((*this)(i,j).imag() > rhs(i,j).imag()) 00176 return true; 00177 if ((*this)(i,j).imag() < rhs(i,j).imag()) 00178 return false; 00179 } 00180 return false; 00181 }
| bool libMesh::TypeTensor< Real >::operator> | ( | const TypeTensor< Real > & | rhs | ) | const [inline, inherited] |
Definition at line 129 of file type_tensor.C.
| bool libMesh::TypeTensor< T >::operator> | ( | const TypeTensor< T > & | rhs | ) | const [inherited] |
- Returns:
trueif this tensor is "greater" than another.
| void libMesh::TypeTensor< T >::print | ( | std::ostream & | os = libMesh::out |
) | const [inline, inherited] |
Formatted print, by default to libMesh::out.
Definition at line 39 of file type_tensor.C.
00040 { 00041 #if LIBMESH_DIM == 1 00042 00043 os << "x=" << (*this)(0,0) << std::endl; 00044 00045 #endif 00046 #if LIBMESH_DIM == 2 00047 00048 os << "(xx,xy)=(" 00049 << std::setw(8) << (*this)(0,0) << ", " 00050 << std::setw(8) << (*this)(0,1) << ")" 00051 << std::endl; 00052 os << "(yx,yy)=(" 00053 << std::setw(8) << (*this)(1,0) << ", " 00054 << std::setw(8) << (*this)(1,1) << ")" 00055 << std::endl; 00056 00057 #endif 00058 #if LIBMESH_DIM == 3 00059 00060 os << "(xx,xy,xz)=(" 00061 << std::setw(8) << (*this)(0,0) << ", " 00062 << std::setw(8) << (*this)(0,1) << ", " 00063 << std::setw(8) << (*this)(0,2) << ")" 00064 << std::endl; 00065 os << "(yx,yy,yz)=(" 00066 << std::setw(8) << (*this)(1,0) << ", " 00067 << std::setw(8) << (*this)(1,1) << ", " 00068 << std::setw(8) << (*this)(1,2) << ")" 00069 << std::endl; 00070 os << "(zx,zy,zz)=(" 00071 << std::setw(8) << (*this)(2,0) << ", " 00072 << std::setw(8) << (*this)(2,1) << ", " 00073 << std::setw(8) << (*this)(2,2) << ")" 00074 << std::endl; 00075 #endif 00076 }
| 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::TypeTensor< T >::_coords, and libMesh::TypeVector< T >::_coords.
00647 { 00648 TypeVector<T> return_vector; 00649 00650 for(unsigned int j=0; j<LIBMESH_DIM; j++) 00651 return_vector._coords[j] = _coords[r*LIBMESH_DIM + j]; 00652 00653 return return_vector; 00654 }
| 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 1027 of file type_tensor.h.
References libMesh::TypeTensor< T >::size_sq().
Referenced by libMesh::System::calculate_norm().
01028 { 01029 return std::sqrt(this->size_sq()); 01030 }
| 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 1086 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords, libMesh::TensorTools::norm_sq(), libMesh::Real, and libMesh::Parallel::sum().
Referenced by libMesh::System::calculate_norm(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::HPCoarsenTest::select_refinement(), and libMesh::TypeTensor< T >::size().
01087 { 01088 Real sum = 0.; 01089 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 01090 sum += TensorTools::norm_sq(_coords[i]); 01091 return sum; 01092 }
| 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.
| 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.
| 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.
Referenced by libMesh::TypeTensor< T >::operator-=().
00777 { 00778 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00779 _coords[i] -= p._coords[i]; 00780 }
| 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().
00788 { 00789 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00790 _coords[i] -= factor*p._coords[i]; 00791 }
| T libMesh::TypeTensor< T >::tr | ( | ) | const [inline, inherited] |
Returns the trace of the tensor.
Definition at line 1059 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
| 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.
References libMesh::TypeTensor< T >::_coords, and libMesh::TypeTensor< T >::TypeTensor().
00933 { 00934 #if LIBMESH_DIM == 1 00935 return TypeTensor(_coords[0]); 00936 #endif 00937 00938 #if LIBMESH_DIM == 2 00939 return TypeTensor(_coords[0], 00940 _coords[2], 00941 _coords[1], 00942 _coords[3]); 00943 #endif 00944 00945 #if LIBMESH_DIM == 3 00946 return TypeTensor(_coords[0], 00947 _coords[3], 00948 _coords[6], 00949 _coords[1], 00950 _coords[4], 00951 _coords[7], 00952 _coords[2], 00953 _coords[5], 00954 _coords[8]); 00955 #endif 00956 }
| void libMesh::TypeTensor< T >::write_unformatted | ( | std::ostream & | out, | |
| const bool | newline = true | |||
| ) | const [inline, 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.
00085 { 00086 libmesh_assert (out_stream); 00087 00088 out_stream << std::setiosflags(std::ios::showpoint) 00089 << (*this)(0,0) << " " 00090 << (*this)(0,1) << " " 00091 << (*this)(0,2) << " "; 00092 if (newline) 00093 out_stream << '\n'; 00094 00095 out_stream << std::setiosflags(std::ios::showpoint) 00096 << (*this)(1,0) << " " 00097 << (*this)(1,1) << " " 00098 << (*this)(1,2) << " "; 00099 if (newline) 00100 out_stream << '\n'; 00101 00102 out_stream << std::setiosflags(std::ios::showpoint) 00103 << (*this)(2,0) << " " 00104 << (*this)(2,1) << " " 00105 << (*this)(2,2) << " "; 00106 if (newline) 00107 out_stream << '\n'; 00108 }
| void libMesh::TypeTensor< T >::zero | ( | ) | [inline, inherited] |
Zero the tensor in any dimension.
Definition at line 1076 of file type_tensor.h.
References libMesh::TypeTensor< T >::_coords.
Referenced by libMesh::TypeTensor< T >::operator=(), and libMesh::TensorValue< T >::operator=().
01077 { 01078 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 01079 _coords[i] = 0.; 01080 }
Friends And Related Function Documentation
| 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 336 of file type_tensor.h.
Member Data Documentation
T libMesh::TypeTensor< T >::_coords[LIBMESH_DIM *LIBMESH_DIM] [inherited] |
The coordinates of the TypeTensor
Definition at line 353 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 >::det(), libMesh::TypeTensor< T >::operator()(), libMesh::TypeTensor< T >::operator*(), libMesh::TypeTensor< T >::operator*=(), libMesh::TypeTensor< T >::operator+(), libMesh::TypeTensor< T >::operator-(), libMesh::TypeTensor< T >::operator/(), libMesh::TypeTensor< T >::operator/=(), libMesh::TypeTensor< T >::operator==(), libMesh::TypeTensor< T >::row(), libMesh::TypeTensor< T >::size_sq(), libMesh::TypeTensor< T >::subtract(), libMesh::TypeTensor< T >::subtract_scaled(), libMesh::TypeTensor< T >::tr(), libMesh::TypeTensor< T >::transpose(), libMesh::TypeTensor< T >::TypeTensor(), and libMesh::TypeTensor< T >::zero().
The documentation for this class was generated from the following file:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:40 UTC
Hosted By: