20 #ifndef LIBMESH_TYPE_TENSOR_H 21 #define LIBMESH_TYPE_TENSOR_H 38 template <
unsigned int N,
typename T>
class TypeNTensor;
50 template <
typename T2>
81 template <
typename Scalar>
92 const Scalar>::type zz=0);
102 template<
typename T2>
106 template<
typename T2>
126 template<
typename T2>
137 template<
typename T2>
145 template <
typename Scalar>
150 { libmesh_assert_equal_to (p, Scalar(0)); this->
zero();
return *
this; }
155 const T &
operator () (
const unsigned int i,
const unsigned int j)
const;
160 T &
operator () (
const unsigned int i,
const unsigned int j);
182 template<
typename T2>
191 template<
typename T2>
197 template<
typename T2>
203 template <
typename T2>
211 template<
typename T2>
220 template<
typename T2>
226 template<
typename T2>
233 template <
typename T2>
246 template <
typename Scalar>
261 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
272 template <
typename Scalar>
291 template <
typename T2>
299 template <
typename T2>
311 template <
typename T2>
321 template <
typename T2>
349 #ifdef LIBMESH_ENABLE_DEPRECATED 365 #ifdef LIBMESH_ENABLE_DEPRECATED 402 bool operator < (const TypeTensor<T> & rhs)
const;
421 friend std::ostream & operator << (std::ostream & os, const TypeTensor<T> & t)
442 template <
typename T>
466 for (
unsigned int i=0; i != LIBMESH_DIM; ++i)
473 const unsigned int _j;
477 template <
typename T>
492 const T &
slice (
const unsigned int i)
const 497 const unsigned int _j;
502 template <
typename T>
525 template <
typename T>
543 #elif LIBMESH_DIM == 1 544 libmesh_assert_equal_to (xy, 0);
545 libmesh_assert_equal_to (yx, 0);
546 libmesh_assert_equal_to (yy, 0);
559 libmesh_assert_equal_to (xz, 0);
560 libmesh_assert_equal_to (yz, 0);
561 libmesh_assert_equal_to (zx, 0);
562 libmesh_assert_equal_to (zy, 0);
563 libmesh_assert_equal_to (zz, 0);
568 template <
typename T>
569 template <
typename Scalar>
581 const Scalar>::type zz)
605 template <
typename T>
606 template<
typename T2>
611 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
616 template <
typename T>
617 template <
typename T2>
620 libmesh_assert_equal_to (LIBMESH_DIM, 1);
624 template <
typename T>
625 template <
typename T2>
629 libmesh_assert_equal_to (LIBMESH_DIM, 2);
636 template <
typename T>
637 template <
typename T2>
642 libmesh_assert_equal_to (LIBMESH_DIM, 3);
657 template <
typename T>
665 template <
typename T>
666 template<
typename T2>
670 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
676 template <
typename T>
679 const unsigned int j)
const 681 libmesh_assert_less (i, 3);
682 libmesh_assert_less (j, 3);
685 const static T my_zero = 0;
686 if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
690 return _coords[i*LIBMESH_DIM+j];
695 template <
typename T>
698 const unsigned int j)
702 if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
703 libmesh_error_msg(
"ERROR: You are assigning to a tensor component that is out of range for the compiled LIBMESH_DIM!");
707 libmesh_assert_less (i, LIBMESH_DIM);
708 libmesh_assert_less (j, LIBMESH_DIM);
710 return _coords[i*LIBMESH_DIM+j];
714 template <
typename T>
719 libmesh_assert_less (i, LIBMESH_DIM);
724 template <
typename T>
729 libmesh_assert_less (i, LIBMESH_DIM);
734 template <
typename T>
741 for (
unsigned int j=0; j<LIBMESH_DIM; j++)
742 return_vector.
_coords[j] = _coords[r*LIBMESH_DIM + j];
744 return return_vector;
747 template <
typename T>
748 template<
typename T2>
782 template <
typename T>
783 template<
typename T2>
794 template <
typename T>
795 template<
typename T2>
799 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
805 template <
typename T>
806 template <
typename T2>
810 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
811 _coords[i] += factor*p.
_coords[i];
817 template <
typename T>
818 template<
typename T2>
852 template <
typename T>
853 template <
typename T2>
864 template <
typename T>
865 template <
typename T2>
869 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
875 template <
typename T>
876 template <
typename T2>
880 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
881 _coords[i] -= factor*p.
_coords[i];
886 template <
typename T>
918 template <
typename T>
919 template <
typename Scalar>
926 typedef decltype((*
this)(0, 0) * factor) TS;
955 template <
typename T,
typename Scalar>
966 template <
typename T>
967 template <
typename Scalar>
969 typename boostcopy::enable_if_c<
971 TypeTensor<typename CompareTypes<T, Scalar>::supertype>>::type
974 libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1005 template <
typename T>
1009 #if LIBMESH_DIM == 1 1013 #if LIBMESH_DIM == 2 1020 #if LIBMESH_DIM == 3 1035 template <
typename T>
1039 #if LIBMESH_DIM == 1 1040 if (_coords[0] == static_cast<T>(0.))
1041 libmesh_convergence_failure();
1045 #if LIBMESH_DIM == 2 1050 T a =
A(0,0), b =
A(0,1),
1051 c =
A(1,0), d =
A(1,1);
1054 T my_det = a*d - b*c;
1056 if (my_det == static_cast<T>(0.))
1057 libmesh_convergence_failure();
1059 return TypeTensor(d/my_det, -b/my_det, -c/my_det, a/my_det);
1062 #if LIBMESH_DIM == 3 1066 T a11 =
A(0,0), a12 =
A(0,1), a13 =
A(0,2),
1067 a21 =
A(1,0), a22 =
A(1,1), a23 =
A(1,2),
1068 a31 =
A(2,0), a32 =
A(2,1), a33 =
A(2,2);
1070 T my_det = a11*(a33*a22-a32*a23) - a21*(a33*a12-a32*a13) + a31*(a23*a12-a22*a13);
1072 if (my_det == static_cast<T>(0.))
1073 libmesh_convergence_failure();
1076 return TypeTensor( (a33*a22-a32*a23)/my_det, -(a33*a12-a32*a13)/my_det, (a23*a12-a22*a13)/my_det,
1077 -(a33*a21-a31*a23)/my_det, (a33*a11-a31*a13)/my_det, -(a23*a11-a21*a13)/my_det,
1078 (a32*a21-a31*a22)/my_det, -(a32*a11-a31*a12)/my_det, (a22*a11-a21*a12)/my_det);
1084 template <
typename T>
1088 #if LIBMESH_DIM == 1 1089 if (_coords[0] == static_cast<T>(0.))
1090 libmesh_convergence_failure();
1091 x(0) = b(0) / _coords[0];
1094 #if LIBMESH_DIM == 2 1095 T my_det = _coords[0]*_coords[3] - _coords[1]*_coords[2];
1097 if (my_det == static_cast<T>(0.))
1098 libmesh_convergence_failure();
1100 T my_det_inv = 1./my_det;
1102 x(0) = my_det_inv*( _coords[3]*b(0) - _coords[1]*b(1));
1103 x(1) = my_det_inv*(-_coords[2]*b(0) + _coords[0]*b(1));
1106 #if LIBMESH_DIM == 3 1109 _coords[0]*(_coords[8]*_coords[4] - _coords[7]*_coords[5])
1111 -_coords[3]*(_coords[8]*_coords[1] - _coords[7]*_coords[2]) +
1113 +_coords[6]*(_coords[5]*_coords[1] - _coords[4]*_coords[2]);
1115 if (my_det == static_cast<T>(0.))
1116 libmesh_convergence_failure();
1118 T my_det_inv = 1./my_det;
1119 x(0) = my_det_inv*((_coords[8]*_coords[4] - _coords[7]*_coords[5])*b(0) -
1120 (_coords[8]*_coords[1] - _coords[7]*_coords[2])*b(1) +
1121 (_coords[5]*_coords[1] - _coords[4]*_coords[2])*b(2));
1123 x(1) = my_det_inv*((_coords[6]*_coords[5] - _coords[8]*_coords[3])*b(0) +
1124 (_coords[8]*_coords[0] - _coords[6]*_coords[2])*b(1) -
1125 (_coords[5]*_coords[0] - _coords[3]*_coords[2])*b(2));
1127 x(2) = my_det_inv*((_coords[7]*_coords[3] - _coords[6]*_coords[4])*b(0) -
1128 (_coords[7]*_coords[0] - _coords[6]*_coords[1])*b(1) +
1129 (_coords[4]*_coords[0] - _coords[3]*_coords[1])*b(2));
1135 template <
typename T>
1139 libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1141 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1142 _coords[i] /= factor;
1150 template <
typename T>
1151 template <
typename T2>
1157 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
1158 for (
unsigned int j=0; j<LIBMESH_DIM; j++)
1159 returnval(i) += (*this)(i,j)*p(j);
1166 template <
typename T>
1167 template <
typename T2>
1172 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
1173 for (
unsigned int j=0; j<LIBMESH_DIM; j++)
1174 for (
unsigned int k=0; k<LIBMESH_DIM; k++)
1175 returnval(i,j) += (*this)(i,k)*p(k,j);
1180 template <
typename T>
1181 template <
typename T2>
1186 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
1187 for (
unsigned int j=0; j<LIBMESH_DIM; j++)
1188 for (
unsigned int k=0; k<LIBMESH_DIM; k++)
1189 temp(i,j) += (*this)(i,k)*p(k,j);
1200 template <
typename T>
1201 template <
typename T2>
1207 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1208 sum += _coords[i]*t.
_coords[i];
1214 #ifdef LIBMESH_ENABLE_DEPRECATED 1215 template <
typename T>
1219 libmesh_deprecated();
1220 return this->norm();
1226 template <
typename T>
1230 return std::sqrt(this->
norm_sq());
1235 template <
typename T>
1239 #if LIBMESH_DIM == 1 1243 #if LIBMESH_DIM == 2 1244 return (_coords[0] * _coords[3]
1245 - _coords[1] * _coords[2]);
1248 #if LIBMESH_DIM == 3 1249 return (_coords[0] * _coords[4] * _coords[8]
1250 + _coords[1] * _coords[5] * _coords[6]
1251 + _coords[2] * _coords[3] * _coords[7]
1252 - _coords[0] * _coords[5] * _coords[7]
1253 - _coords[1] * _coords[3] * _coords[8]
1254 - _coords[2] * _coords[4] * _coords[6]);
1258 template <
typename T>
1262 #if LIBMESH_DIM == 1 1266 #if LIBMESH_DIM == 2 1267 return _coords[0] + _coords[3];
1270 #if LIBMESH_DIM == 3 1271 return _coords[0] + _coords[4] + _coords[8];
1275 template <
typename T>
1279 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1285 #ifdef LIBMESH_ENABLE_DEPRECATED 1286 template <
typename T>
1290 libmesh_deprecated();
1297 template <
typename T>
1302 for (
unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1309 template <
typename T>
1313 #if LIBMESH_DIM == 1 1318 #if LIBMESH_DIM == 2 1326 #if LIBMESH_DIM == 3 1343 #endif // LIBMESH_TYPE_TENSOR_H void write_unformatted(std::ostream &out, const bool newline=true) const
const T & operator()(const unsigned int i) const
TypeTensorColumn< T > & operator=(const TypeVector< T > &rhs)
const T & slice(const unsigned int i) const
TypeTensor< T > * _tensor
T & operator()(const unsigned int i)
void add_scaled(const TypeTensor< T2 > &, const T)
std::tuple< unsigned int, unsigned int > index_type
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< typename CompareTypes< T, Scalar >::supertype > >::type operator/(const Scalar) const
TypeTensor< T > transpose() const
ConstTypeTensorColumn(const TypeTensor< T > &tensor, unsigned int j)
void solve(const TypeVector< T > &b, TypeVector< T > &x) const
TypeVector< T > row(const unsigned int r) const
static const Real TOLERANCE
T _coords[LIBMESH_DIM *LIBMESH_DIM]
const TypeTensor< T > * _tensor
TypeTensor< typename CompareTypes< T, T2 >::supertype > operator+(const TypeTensor< T2 > &) const
TypeTensorColumn(TypeTensor< T > &tensor, unsigned int j)
auto operator*(const Scalar scalar) const -> typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< decltype(T() *scalar)>>::type
TypeTensor< T > operator-() const
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeNTensor< N, typename CompareTypes< Scalar, T >::supertype > >::type operator*(const Scalar &, const TypeNTensor< N, T > &)
const TypeTensor< T > & operator/=(const T)
void subtract(const TypeTensor< T2 > &)
void assign(const TypeTensor< T2 > &)
void add(const TypeTensor< T2 > &)
CompareTypes< T, T2 >::supertype contract(const TypeTensor< T2 > &) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const TypeTensor< T > & operator*=(const Scalar factor)
TypeTensor< T > inverse() const
static PetscErrorCode Mat * A
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor & >::type operator=(const Scalar &libmesh_dbg_var(p))
const TypeTensor< T > & operator-=(const TypeTensor< T2 > &)
ConstTypeTensorColumn< T > slice(const unsigned int i) const
const TypeTensor< T > & operator+=(const TypeTensor< T2 > &)
T & slice(const unsigned int i)
bool operator==(const TypeTensor< T > &rhs) const
bool operator>(const TypeTensor< T > &rhs) const
OStreamProxy out(std::cout)
void print(std::ostream &os=libMesh::out) const
const T & operator()(const unsigned int i, const unsigned int j) const
void subtract_scaled(const TypeTensor< T2 > &, const T)