16 #ifndef dealii_tensor_h 17 #define dealii_tensor_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/exceptions.h> 21 #include <deal.II/base/numbers.h> 22 #include <deal.II/base/table_indices.h> 23 #include <deal.II/base/tensor_accessors.h> 24 #include <deal.II/base/template_constraints.h> 25 #include <deal.II/base/utilities.h> 27 #ifdef DEAL_II_WITH_ADOLC 28 #include <adolc/adouble.h> 36 DEAL_II_NAMESPACE_OPEN
40 template <
int dim,
typename Number>
class Point;
41 template <
int rank_,
int dim,
typename Number =
double>
class Tensor;
42 template <
typename Number>
class Vector;
48 template <
int dim,
typename Number>
49 class Tensor<-2, dim, Number>
53 template <
int dim,
typename Number>
54 class Tensor<-1, dim, Number>
90 template <
int dim,
typename Number>
107 static const unsigned int rank = 0;
142 DEAL_II_CUDA_HOST_DEV
Tensor ();
149 template <
typename OtherNumber>
155 template <
typename OtherNumber>
156 Tensor (
const OtherNumber &initializer);
192 DEAL_II_CUDA_HOST_DEV
operator Number &();
202 DEAL_II_CUDA_HOST_DEV
operator const Number &()
const;
209 template <
typename OtherNumber>
212 #ifdef __INTEL_COMPILER 226 template <
typename OtherNumber>
232 template <
typename OtherNumber>
238 template <
typename OtherNumber>
244 template <
typename OtherNumber>
250 template <
typename OtherNumber>
258 template <
typename OtherNumber>
264 template <
typename OtherNumber>
305 template <
class Archive>
306 void serialize(Archive &ar,
const unsigned int version);
323 template <
typename OtherNumber>
325 unsigned int &start_index)
const;
330 template <
int,
int,
typename>
friend class Tensor;
376 template <
int rank_,
int dim,
typename Number>
393 static const unsigned int rank = rank_;
399 static const unsigned int 423 DEAL_II_CUDA_HOST_DEV
Tensor ();
435 template <
typename OtherNumber>
441 template <
typename OtherNumber>
447 template <
typename OtherNumber>
448 operator Tensor<1,dim,
Tensor<rank_-1,dim,OtherNumber> > ()
const;
503 template <
typename OtherNumber>
517 template <
typename OtherNumber>
523 template <
typename OtherNumber>
529 template <
typename OtherNumber>
535 template <
typename OtherNumber>
544 template <
typename OtherNumber>
550 template <
typename OtherNumber>
595 template <
typename OtherNumber>
596 void unroll (Vector<OtherNumber> &result)
const;
623 template <
class Archive>
624 void serialize(Archive &ar,
const unsigned int version);
643 template <
typename OtherNumber>
645 unsigned int &start_index)
const;
650 template <
int,
int,
typename>
friend class Tensor;
672 template <
int rank,
int dim,
typename T>
688 template <
int rank,
int dim,
typename T>
716 template <
int dim,
typename Number>
721 : value(
internal::NumberType<Number>::value(0.0))
727 template <
int dim,
typename Number>
728 template <
typename OtherNumber>
737 template <
int dim,
typename Number>
738 template <
typename OtherNumber>
747 template <
int dim,
typename Number>
752 return std::addressof(value);
757 template <
int dim,
typename Number>
762 return std::addressof(value);
767 template <
int dim,
typename Number>
777 template <
int dim,
typename Number>
787 template <
int dim,
typename Number>
792 #ifndef __CUDA_ARCH__ 793 Assert(dim != 0,
ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
799 template <
int dim,
typename Number>
804 #ifndef __CUDA_ARCH__ 805 Assert(dim != 0,
ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
811 template <
int dim,
typename Number>
812 template <
typename OtherNumber>
821 #ifdef __INTEL_COMPILER 822 template <
int dim,
typename Number>
832 template <
int dim,
typename Number>
833 template <
typename OtherNumber>
842 template <
int dim,
typename Number>
843 template <
typename OtherNumber>
847 #ifdef DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING 848 Assert(!(std::is_same<Number,adouble>::value || std::is_same<OtherNumber,adouble>::value),
849 ExcMessage(
"The Tensor equality operator for Adol-C taped numbers has not yet " 850 "been extended to support advanced branching."));
857 template <
int dim,
typename Number>
858 template <
typename OtherNumber>
862 return !((*this) == p);
866 template <
int dim,
typename Number>
867 template <
typename OtherNumber>
876 template <
int dim,
typename Number>
877 template <
typename OtherNumber>
886 template <
int dim,
typename Number>
887 template <
typename OtherNumber>
896 template <
int dim,
typename Number>
897 template <
typename OtherNumber>
906 template <
int dim,
typename Number>
914 template <
int dim,
typename Number>
919 Assert(dim != 0,
ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
924 template <
int dim,
typename Number>
930 #ifndef __CUDA_ARCH__ 931 Assert(dim != 0,
ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
937 template <
int dim,
typename Number>
938 template <
typename OtherNumber>
942 unsigned int &index)
const 944 Assert(dim != 0,
ExcMessage(
"Cannot unroll an object of type Tensor<0,0,Number>"));
945 result[index] = value;
950 template <
int dim,
typename Number>
960 template <
int dim,
typename Number>
961 template <
class Archive>
972 template <
int rank_,
int dim,
typename Number>
973 inline DEAL_II_ALWAYS_INLINE
981 template <
int rank_,
int dim,
typename Number>
982 inline DEAL_II_ALWAYS_INLINE
985 for (
unsigned int i=0; i<dim; ++i)
990 template <
int rank_,
int dim,
typename Number>
991 template <
typename OtherNumber>
992 inline DEAL_II_ALWAYS_INLINE
995 for (
unsigned int i=0; i!=dim; ++i)
1000 template <
int rank_,
int dim,
typename Number>
1001 template <
typename OtherNumber>
1002 inline DEAL_II_ALWAYS_INLINE
1006 for (
unsigned int i=0; i<dim; ++i)
1007 values[i] = initializer[i];
1011 template <
int rank_,
int dim,
typename Number>
1012 template <
typename OtherNumber>
1013 inline DEAL_II_ALWAYS_INLINE
1017 return Tensor<1,dim,
Tensor<rank_-1,dim,Number> > (values);
1024 namespace TensorSubscriptor
1026 template <
typename ArrayElementType,
int dim>
1027 inline DEAL_II_ALWAYS_INLINE
1028 DEAL_II_CUDA_HOST_DEV
1030 subscript (ArrayElementType *values,
1031 const unsigned int i,
1032 std::integral_constant<int, dim>)
1035 #ifndef __CUDA_ARCH__ 1042 template <
typename ArrayElementType>
1043 DEAL_II_CUDA_HOST_DEV
1045 subscript (ArrayElementType *,
1047 std::integral_constant<int, 0>)
1050 #ifndef __CUDA_ARCH__ 1051 Assert(
false,
ExcMessage(
"Cannot access elements of an object of type Tensor<rank,0,Number>."));
1053 static ArrayElementType t;
1060 template <
int rank_,
int dim,
typename Number>
1061 inline DEAL_II_ALWAYS_INLINE
1062 DEAL_II_CUDA_HOST_DEV
1066 return ::internal::TensorSubscriptor::subscript(values, i, std::integral_constant<int, dim>());
1070 template <
int rank_,
int dim,
typename Number>
1071 inline DEAL_II_ALWAYS_INLINE
1072 DEAL_II_CUDA_HOST_DEV
1076 return ::internal::TensorSubscriptor::subscript(values, i, std::integral_constant<int, dim>());
1080 template <
int rank_,
int dim,
typename Number>
1085 Assert(dim != 0,
ExcMessage(
"Cannot access an object of type Tensor<rank_,0,Number>"));
1087 return TensorAccessors::extract<rank_>(*
this, indices);
1092 template <
int rank_,
int dim,
typename Number>
1097 Assert(dim != 0,
ExcMessage(
"Cannot access an object of type Tensor<rank_,0,Number>"));
1099 return TensorAccessors::extract<rank_>(*
this, indices);
1104 template <
int rank_,
int dim,
typename Number>
1109 return std::addressof(this->
operator[](this->unrolled_to_component_indices(0)));
1114 template <
int rank_,
int dim,
typename Number>
1119 return std::addressof(this->
operator[](this->unrolled_to_component_indices(0)));
1124 template <
int rank_,
int dim,
typename Number>
1129 return begin_raw()+n_independent_components;
1134 template <
int rank_,
int dim,
typename Number>
1139 return begin_raw()+n_independent_components;
1144 template <
int rank_,
int dim,
typename Number>
1145 template <
typename OtherNumber>
1146 inline DEAL_II_ALWAYS_INLINE
1151 std::copy (&t.
values[0], &t.
values[0]+dim, &values[0]);
1156 template <
int rank_,
int dim,
typename Number>
1157 inline DEAL_II_ALWAYS_INLINE
1164 for (
unsigned int i=0; i<dim; ++i)
1170 template <
int rank_,
int dim,
typename Number>
1171 template <
typename OtherNumber>
1176 for (
unsigned int i=0; i<dim; ++i)
1177 if (values[i] != p.
values[i])
1197 template <
int rank_,
int dim,
typename Number>
1198 template <
typename OtherNumber>
1203 return !((*this) == p);
1207 template <
int rank_,
int dim,
typename Number>
1208 template <
typename OtherNumber>
1213 for (
unsigned int i=0; i<dim; ++i)
1214 values[i] += p.
values[i];
1219 template <
int rank_,
int dim,
typename Number>
1220 template <
typename OtherNumber>
1225 for (
unsigned int i=0; i<dim; ++i)
1226 values[i] -= p.
values[i];
1231 template <
int rank_,
int dim,
typename Number>
1232 template <
typename OtherNumber>
1234 DEAL_II_CUDA_HOST_DEV
1238 for (
unsigned int i=0; i<dim; ++i)
1244 template <
int rank_,
int dim,
typename Number>
1245 template <
typename OtherNumber>
1250 for (
unsigned int i=0; i<dim; ++i)
1256 template <
int rank_,
int dim,
typename Number>
1263 for (
unsigned int i=0; i<dim; ++i)
1264 tmp.
values[i] = -values[i];
1270 template <
int rank_,
int dim,
typename Number>
1275 return std::sqrt (norm_square());
1279 template <
int rank_,
int dim,
typename Number>
1281 DEAL_II_CUDA_HOST_DEV
1287 for (
unsigned int i=0; i<dim; ++i)
1288 s += values[i].norm_square();
1294 template <
int rank_,
int dim,
typename Number>
1295 template <
typename OtherNumber>
1300 AssertDimension (result.size(),(Utilities::fixed_power<rank_, unsigned int>(dim)));
1302 unsigned int index = 0;
1303 unroll_recursion (result, index);
1307 template <
int rank_,
int dim,
typename Number>
1308 template <
typename OtherNumber>
1312 unsigned int &index)
const 1314 for (
unsigned int i=0; i<dim; ++i)
1315 values[i].unroll_recursion(result, index);
1319 template <
int rank_,
int dim,
typename Number>
1324 unsigned int index = 0;
1325 for (
int r = 0; r < rank_; ++r)
1326 index = index * dim + indices[r];
1332 template <
int rank_,
int dim,
typename Number>
1337 Assert (i < n_independent_components,
1342 unsigned int remainder = i;
1343 for (
int r=rank_-1; r>=0; --r)
1345 indices[r] = (remainder % dim);
1354 template <
int rank_,
int dim,
typename Number>
1358 for (
unsigned int i=0; i<dim; ++i)
1363 template <
int rank_,
int dim,
typename Number>
1372 template <
int rank_,
int dim,
typename Number>
1373 template <
class Archive>
1396 template <
int rank_,
int dim,
typename Number>
1398 std::ostream &operator << (std::ostream &out, const Tensor<rank_,dim,Number> &p)
1400 for (
unsigned int i = 0; i < dim; ++i)
1417 template <
int dim,
typename Number>
1419 std::ostream &operator << (std::ostream &out, const Tensor<0,dim,Number> &p)
1421 out << static_cast<const Number &>(p);
1441 template <
int dim,
typename Number,
typename Other>
1442 inline DEAL_II_ALWAYS_INLINE
1443 typename ProductType<Other, Number>::type
1447 return object *
static_cast<const Number &
>(t);
1460 template <
int dim,
typename Number,
typename Other>
1461 inline DEAL_II_ALWAYS_INLINE
1462 typename ProductType<Number, Other>::type
1464 const Other &
object)
1466 return static_cast<const Number &
>(t) *
object;
1479 template <
int dim,
typename Number,
typename OtherNumber>
1480 inline DEAL_II_ALWAYS_INLINE
1481 typename ProductType<Number, OtherNumber>::type
1485 return static_cast<const Number &
>(src1) *
1486 static_cast<const OtherNumber &>(src2);
1495 template <
int dim,
typename Number,
typename OtherNumber>
1496 inline DEAL_II_ALWAYS_INLINE
1499 const OtherNumber &factor)
1501 return static_cast<const Number &
>(t) / factor;
1510 template <
int dim,
typename Number,
typename OtherNumber>
1511 inline DEAL_II_ALWAYS_INLINE
1516 return static_cast<const Number &
>(p) + static_cast<const OtherNumber &>(q);
1525 template <
int dim,
typename Number,
typename OtherNumber>
1526 inline DEAL_II_ALWAYS_INLINE
1531 return static_cast<const Number &
>(p) - static_cast<const OtherNumber &>(q);
1545 template <
int rank,
int dim,
1547 typename OtherNumber>
1548 inline DEAL_II_ALWAYS_INLINE
1551 const OtherNumber &factor)
1555 for (
unsigned int d=0; d<dim; ++d)
1556 tt[d] = t[d] * factor;
1571 template <
int rank,
int dim,
1573 typename OtherNumber>
1574 inline DEAL_II_ALWAYS_INLINE
1591 template <
int rank,
int dim,
1593 typename OtherNumber>
1597 const OtherNumber &factor)
1601 for (
unsigned int d=0; d<dim; ++d)
1602 tt[d] = t[d] / factor;
1614 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1615 inline DEAL_II_ALWAYS_INLINE
1622 for (
unsigned int i=0; i<dim; ++i)
1636 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1637 inline DEAL_II_ALWAYS_INLINE
1644 for (
unsigned int i=0; i<dim; ++i)
1681 template <
int rank_1,
int rank_2,
int dim,
1682 typename Number,
typename OtherNumber>
1683 inline DEAL_II_ALWAYS_INLINE
1684 typename Tensor<rank_1 + rank_2 - 2, dim,
typename ProductType<Number, OtherNumber>::type>::tensor_type
1688 typename Tensor<rank_1 + rank_2 - 2, dim,
typename ProductType<Number, OtherNumber>::type>
::tensor_type result;
1690 TensorAccessors::internal::ReorderedIndexView<0, rank_2, const Tensor<rank_2, dim, OtherNumber> >
1691 reordered = TensorAccessors::reordered_index_view<0, rank_2>(src2);
1692 TensorAccessors::contract<1, rank_1, rank_2, dim>(result, src1, reordered);
1727 template <
int index_1,
int index_2,
1728 int rank_1,
int rank_2,
int dim,
1729 typename Number,
typename OtherNumber>
1730 inline DEAL_II_ALWAYS_INLINE
1731 typename Tensor<rank_1 + rank_2 - 2, dim,
typename ProductType<Number, OtherNumber>::type>::tensor_type
1735 Assert(0 <= index_1 && index_1 < rank_1,
1736 ExcMessage(
"The specified index_1 must lie within the range [0,rank_1)"));
1737 Assert(0 <= index_2 && index_2 < rank_2,
1738 ExcMessage(
"The specified index_2 must lie within the range [0,rank_2)"));
1744 ReorderedIndexView<index_1, rank_1, const Tensor<rank_1, dim, Number> >
1745 reord_01 = reordered_index_view<index_1, rank_1>(src1);
1748 ReorderedIndexView<index_2, rank_2, const Tensor<rank_2, dim, OtherNumber> >
1749 reord_02 = reordered_index_view<index_2, rank_2>(src2);
1751 typename Tensor<rank_1 + rank_2 - 2, dim,
typename ProductType<Number, OtherNumber>::type>
::tensor_type 1753 TensorAccessors::contract<1, rank_1, rank_2, dim>(result, reord_01, reord_02);
1789 template <
int index_1,
int index_2,
int index_3,
int index_4,
1790 int rank_1,
int rank_2,
int dim,
1791 typename Number,
typename OtherNumber>
1793 typename Tensor<rank_1 + rank_2 - 4, dim,
typename ProductType<Number, OtherNumber>::type>::tensor_type
1797 Assert(0 <= index_1 && index_1 < rank_1,
1798 ExcMessage(
"The specified index_1 must lie within the range [0,rank_1)"));
1799 Assert(0 <= index_3 && index_3 < rank_1,
1800 ExcMessage(
"The specified index_3 must lie within the range [0,rank_1)"));
1801 Assert(index_1 != index_3,
1802 ExcMessage(
"index_1 and index_3 must not be the same"));
1803 Assert(0 <= index_2 && index_2 < rank_2,
1804 ExcMessage(
"The specified index_2 must lie within the range [0,rank_2)"));
1805 Assert(0 <= index_4 && index_4 < rank_2,
1806 ExcMessage(
"The specified index_4 must lie within the range [0,rank_2)"));
1807 Assert(index_2 != index_4,
1808 ExcMessage(
"index_2 and index_4 must not be the same"));
1814 ReorderedIndexView<index_1, rank_1, const Tensor<rank_1, dim, Number> >
1815 reord_1 = TensorAccessors::reordered_index_view<index_1, rank_1>(src1);
1818 ReorderedIndexView<index_2, rank_2, const Tensor<rank_2, dim, OtherNumber> >
1819 reord_2 = TensorAccessors::reordered_index_view<index_2, rank_2>(src2);
1824 ReorderedIndexView<(index_3 < index_1 ? index_3 : index_3 - 1), rank_1, ReorderedIndexView<index_1, rank_1,
const Tensor<rank_1, dim, Number> > >
1833 typename Tensor<rank_1 + rank_2 - 4, dim,
typename ProductType<Number, OtherNumber>::type>
::tensor_type 1835 TensorAccessors::contract<2, rank_1, rank_2, dim>(result, reord_3, reord_4);
1853 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1854 inline DEAL_II_ALWAYS_INLINE
1855 typename ProductType<Number, OtherNumber>::type
1859 typename ProductType<Number, OtherNumber>::type result;
1860 TensorAccessors::contract<rank, rank, rank, dim>(result, left, right);
1883 template <
template <
int,
int,
typename>
class TensorT1,
1884 template <
int,
int,
typename>
class TensorT2,
1885 template <
int,
int,
typename>
class TensorT3,
1886 int rank_1,
int rank_2,
int dim,
1887 typename T1,
typename T2,
typename T3>
1890 const TensorT2<rank_1 + rank_2, dim, T2> &middle,
1891 const TensorT3<rank_2, dim, T3> &right)
1895 return TensorAccessors::contract3<rank_1, rank_2, dim, return_type>(
1896 left, middle, right);
1911 template <
int rank_1,
int rank_2,
int dim,
1912 typename Number,
typename OtherNumber>
1913 inline DEAL_II_ALWAYS_INLINE
1919 TensorAccessors::contract<0, rank_1, rank_2, dim>(result, src1, src2);
1942 template <
int dim,
typename Number>
1943 inline DEAL_II_ALWAYS_INLINE
1952 result[1] = -src[0];
1968 template <
int dim,
typename Number>
1969 inline DEAL_II_ALWAYS_INLINE
1978 result[0] = src1[1]*src2[2] - src1[2]*src2[1];
1979 result[1] = src1[2]*src2[0] - src1[0]*src2[2];
1980 result[2] = src1[0]*src2[1] - src1[1]*src2[0];
1999 template <
int dim,
typename Number>
2007 for (
unsigned int k=0; k<dim; ++k)
2009 Tensor<2,dim-1,Number> minor;
2010 for (
unsigned int i=0; i<dim-1; ++i)
2011 for (
unsigned int j=0; j<dim-1; ++j)
2012 minor[i][j] = t[i][j<k ? j : j+1];
2014 const Number cofactor = ((k % 2 == 0) ? -1. : 1.) *
determinant(minor);
2016 det += t[dim-1][k] * cofactor;
2019 return ((dim % 2 == 0) ? 1. : -1.) * det;
2027 template <
typename Number>
2042 template <
int dim,
typename Number>
2043 inline DEAL_II_ALWAYS_INLINE
2047 for (
unsigned int i=1; i<dim; ++i)
2062 template <
int dim,
typename Number>
2067 Number return_tensor [dim][dim];
2080 template <
typename Number>
2085 Number return_tensor [1][1];
2093 template <
typename Number>
2103 return_tensor[0][0] = t[1][1];
2104 return_tensor[0][1] = -t[0][1];
2105 return_tensor[1][0] = -t[1][0];
2106 return_tensor[1][1] = t[0][0];
2107 return_tensor *= inv_det_t;
2109 return return_tensor;
2113 template <
typename Number>
2127 1.0/(t4*t[2][2]-t6*t[2][1]-t8*t[2][2]+
2128 t00*t[2][1]+t01*t[1][2]-t04*t[1][1]));
2134 return_tensor[1][2] = t00-t6;
2138 return_tensor *= inv_det_t;
2140 return return_tensor;
2152 template <
int dim,
typename Number>
2153 inline DEAL_II_ALWAYS_INLINE
2158 for (
unsigned int i=0; i<dim; ++i)
2161 for (
unsigned int j=i+1; j<dim; ++j)
2184 template <
int dim,
typename Number>
2207 template <
int dim,
typename Number>
2223 template <
int dim,
typename Number>
2229 for (
unsigned int j=0; j<dim; ++j)
2232 for (
unsigned int i=0; i<dim; ++i)
2233 sum += std::fabs(t[i][j]);
2250 template <
int dim,
typename Number>
2256 for (
unsigned int i=0; i<dim; ++i)
2259 for (
unsigned int j=0; j<dim; ++j)
2260 sum += std::fabs(t[i][j]);
2275 #ifdef DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING 2285 for (
unsigned int j=0; j<dim; ++j)
2288 for (
unsigned int i=0; i<dim; ++i)
2289 sum += std::fabs(t[i][j]);
2291 condassign(max, (sum > max), sum, max);
2304 for (
unsigned int i=0; i<dim; ++i)
2307 for (
unsigned int j=0; j<dim; ++j)
2308 sum += std::fabs(t[i][j]);
2310 condassign(max, (sum > max), sum, max);
2316 #endif // DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING 2321 DEAL_II_NAMESPACE_CLOSE
2324 #include <deal.II/base/tensor_deprecated.h> SymmetricTensor< rank_, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/(const SymmetricTensor< rank_, dim, Number > &t, const OtherNumber &factor)
numbers::NumberTraits< Number >::real_type real_type
Number trace(const Tensor< 2, dim, Number > &d)
#define AssertDimension(dim1, dim2)
Number determinant(const SymmetricTensor< 2, dim, Number > &)
static unsigned int component_to_unrolled_index(const TableIndices< rank_ > &indices)
static const unsigned int n_independent_components
bool value_is_zero(const Number &value)
bool operator==(const Tensor< rank_, dim, OtherNumber > &) const
Number l1_norm(const Tensor< 2, dim, Number > &t)
static std::size_t memory_consumption()
Tensor & operator/=(const OtherNumber &factor)
bool operator!=(const Tensor< rank_, dim, OtherNumber > &) const
Number linfty_norm(const Tensor< 2, dim, Number > &t)
Tensor< rank_1+rank_2 - 4, dim, typename ProductType< Number, OtherNumber >::type >::tensor_type double_contract(const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
Tensor< 2, dim, Number > adjugate(const Tensor< 2, dim, Number > &t)
static const unsigned int rank
#define AssertThrow(cond, exc)
static real_type abs(const number &x)
internal::ReorderedIndexView< index, rank, T > reordered_index_view(T &t)
numbers::NumberTraits< Number >::real_type norm() const
SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator+(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
void unroll_recursion(Vector< OtherNumber > &result, unsigned int &start_index) const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Tensor & operator=(const Tensor< rank_, dim, OtherNumber > &rhs)
bool values_are_equal(const Number1 &value_1, const Number2 &value_2)
Tensor & operator*=(const OtherNumber &factor)
Tensor< 2, dim, Number > cofactor(const Tensor< 2, dim, Number > &t)
SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator-(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
Tensor< 2, dim, Number > transpose(const Tensor< 2, dim, Number > &t)
static ::ExceptionBase & ExcMessage(std::string arg1)
ProductType< T1, typename ProductType< T2, T3 >::type >::type contract3(const TensorT1< rank_1, dim, T1 > &left, const TensorT2< rank_1+rank_2, dim, T2 > &middle, const TensorT3< rank_2, dim, T3 > &right)
T sum(const T &t, const MPI_Comm &mpi_communicator)
Tensor< 2, dim, Number > invert(const Tensor< 2, dim, Number > &)
static real_type abs_square(const number &x)
Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
#define Assert(cond, exc)
ProductType< Number, OtherNumber >::type scalar_product(const Tensor< rank, dim, Number > &left, const Tensor< rank, dim, OtherNumber > &right)
Tensor< rank_, dim, Number > tensor_type
void serialize(Archive &ar, const unsigned int version)
numbers::NumberTraits< Number >::real_type norm_square() const
value_type & operator[](const unsigned int i)
Tensor & operator-=(const Tensor< rank_, dim, OtherNumber > &)
SymmetricTensor< rank_, dim, Number > transpose(const SymmetricTensor< rank_, dim, Number > &t)
Tensor< rank_-1, dim, Number >::array_type array_type[(dim !=0) ? dim :1]
static TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
Tensor< rank_-1, dim, Number > values[(dim !=0) ? dim :1]
Tensor< 1, dim, Number > cross_product_3d(const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2)
Number determinant(const Tensor< 2, dim, Number > &t)
static const unsigned int dimension
Tensor< rank_-1, dim, Number >::tensor_type value_type
static ::ExceptionBase & ExcNotImplemented()
void unroll(Vector< OtherNumber > &result) const
Tensor< rank_1+rank_2 - 2, dim, typename ProductType< Number, OtherNumber >::type >::tensor_type contract(const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
Tensor< rank_1+rank_2, dim, typename ProductType< Number, OtherNumber >::type > outer_product(const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
Tensor & operator+=(const Tensor< rank_, dim, OtherNumber > &)
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber) const
T max(const T &t, const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcInternalError()
Number determinant(const Tensor< 2, 1, Number > &t)