16#ifndef dealii_tensor_h
17#define dealii_tensor_h
28#ifdef DEAL_II_WITH_ADOLC
29# include <adolc/adouble.h>
42template <
typename ElementType,
typename MemorySpace>
44template <
int dim,
typename Number>
46template <
int rank_,
int dim,
typename Number =
double>
48template <
typename Number>
50template <
typename number>
91template <
int dim,
typename Number>
95 static_assert(dim >= 0,
96 "Tensors must have a dimension greater than or equal to one.");
111 static constexpr unsigned int rank = 0;
156 template <
typename OtherNumber>
165 template <
typename OtherNumber>
169#ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
236 template <
typename OtherNumber>
240#if defined(__INTEL_COMPILER) || defined(DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG)
253#ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
267 template <
typename OtherNumber>
274 template <
typename OtherNumber>
281 template <
typename OtherNumber>
290 template <
typename OtherNumber>
299 template <
typename OtherNumber>
308 template <
typename OtherNumber>
317 template <
typename OtherNumber>
366 template <
class Archive>
385 template <
typename OtherNumber>
388 unsigned int & start_index)
const;
391 template <
int,
int,
typename>
470template <
int rank_,
int dim,
typename Number>
474 static_assert(rank_ >= 1,
475 "Tensors must have a rank greater than or equal to one.");
476 static_assert(dim >= 0,
477 "Tensors must have a dimension greater than or equal to one.");
491 static constexpr unsigned int rank = rank_;
543 template <
typename ElementType,
typename MemorySpace>
554 template <
typename OtherNumber>
561 template <
typename OtherNumber>
568 template <
typename OtherNumber>
572#ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
640 template <
typename OtherNumber>
653#ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
670 template <
typename OtherNumber>
677 template <
typename OtherNumber>
686 template <
typename OtherNumber>
695 template <
typename OtherNumber>
705 template <
typename OtherNumber>
714 template <
typename OtherNumber>
769 template <
typename OtherNumber>
777 static constexpr unsigned int
792 static constexpr std::size_t
800 template <
class Archive>
821 template <
typename OtherNumber>
824 unsigned int & start_index)
const;
832 template <
typename ArrayLike, std::size_t... Indices>
834 Tensor(
const ArrayLike &initializer, std::index_sequence<Indices...>);
837 template <
int,
int,
typename>
842 friend class Point<dim, Number>;
853 template <
int rank,
int dim,
typename T,
typename U>
854 struct ProductTypeImpl<
Tensor<rank, dim,
T>,
std::complex<U>>
860 template <
int rank,
int dim,
typename T,
typename U>
861 struct ProductTypeImpl<
Tensor<rank, dim,
std::complex<T>>, std::complex<U>>
867 template <
typename T,
int rank,
int dim,
typename U>
868 struct ProductTypeImpl<
std::complex<T>,
Tensor<rank, dim, U>>
874 template <
int rank,
int dim,
typename T,
typename U>
875 struct ProductTypeImpl<
std::complex<T>,
Tensor<rank, dim, std::complex<U>>>
886 template <
int rank,
int dim,
typename T>
887 struct NumberType<
Tensor<rank, dim,
T>>
909template <
int dim,
typename Number>
919template <
int dim,
typename Number>
920template <
typename OtherNumber>
928template <
int dim,
typename Number>
929template <
typename OtherNumber>
936# ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
937template <
int dim,
typename Number>
945template <
int dim,
typename Number>
948 :
value{std::move(other.value)}
953template <
int dim,
typename Number>
957 return std::addressof(value);
962template <
int dim,
typename Number>
966 return std::addressof(value);
971template <
int dim,
typename Number>
980template <
int dim,
typename Number>
989template <
int dim,
typename Number>
994# ifndef __CUDA_ARCH__
996 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
1002template <
int dim,
typename Number>
1007# ifndef __CUDA_ARCH__
1009 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
1015template <
int dim,
typename Number>
1016template <
typename OtherNumber>
1026# if defined(__INTEL_COMPILER) || defined(DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG)
1027template <
int dim,
typename Number>
1037# ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
1038template <
int dim,
typename Number>
1042 value = std::move(other.value);
1049template <
int dim,
typename Number>
1050template <
typename OtherNumber>
1060template <
int dim,
typename Number>
1061template <
typename OtherNumber>
1062constexpr inline bool
1065# if defined(DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING)
1066 Assert(!(std::is_same<Number, adouble>::value ||
1067 std::is_same<OtherNumber, adouble>::value),
1069 "The Tensor equality operator for ADOL-C taped numbers has not yet "
1070 "been extended to support advanced branching."));
1077template <
int dim,
typename Number>
1078template <
typename OtherNumber>
1082 return !((*this) == p);
1086template <
int dim,
typename Number>
1087template <
typename OtherNumber>
1097template <
int dim,
typename Number>
1098template <
typename OtherNumber>
1111 namespace ComplexWorkaround
1113 template <
typename Number,
typename OtherNumber>
1115 multiply_assign_scalar(Number &val,
const OtherNumber &s)
1120# ifdef __CUDA_ARCH__
1121 template <
typename Number,
typename OtherNumber>
1123 multiply_assign_scalar(std::complex<Number> &,
const OtherNumber &)
1125 printf(
"This function is not implemented for std::complex<Number>!\n");
1133template <
int dim,
typename Number>
1134template <
typename OtherNumber>
1139 internal::ComplexWorkaround::multiply_assign_scalar(value, s);
1145template <
int dim,
typename Number>
1146template <
typename OtherNumber>
1155template <
int dim,
typename Number>
1163template <
int dim,
typename Number>
1168 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
1173template <
int dim,
typename Number>
1179# ifndef __CUDA_ARCH__
1181 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
1187template <
int dim,
typename Number>
1188template <
typename OtherNumber>
1191 unsigned int & index)
const
1194 ExcMessage(
"Cannot unroll an object of type Tensor<0,0,Number>"));
1195 result[index] =
value;
1200template <
int dim,
typename Number>
1201constexpr inline void
1210template <
int dim,
typename Number>
1211template <
class Archive>
1219template <
int dim,
typename Number>
1225template <
int rank_,
int dim,
typename Number>
1226template <
typename ArrayLike, std::size_t... indices>
1229 std::index_sequence<indices...>)
1230 :
values{
Tensor<rank_ - 1, dim, Number>(initializer[indices])...}
1232 static_assert(
sizeof...(indices) == dim,
1233 "dim should match the number of indices");
1238template <
int rank_,
int dim,
typename Number>
1248template <
int rank_,
int dim,
typename Number>
1251 :
Tensor(initializer,
std::make_index_sequence<dim>{})
1256template <
int rank_,
int dim,
typename Number>
1257template <
typename ElementType,
typename MemorySpace>
1264 for (
unsigned int i = 0; i < n_independent_components; ++i)
1265 (*
this)[unrolled_to_component_indices(i)] = initializer[i];
1270template <
int rank_,
int dim,
typename Number>
1271template <
typename OtherNumber>
1275 :
Tensor(initializer,
std::make_index_sequence<dim>{})
1280template <
int rank_,
int dim,
typename Number>
1281template <
typename OtherNumber>
1285 :
Tensor(initializer,
std::make_index_sequence<dim>{})
1290template <
int rank_,
int dim,
typename Number>
1291template <
typename OtherNumber>
1293 operator
Tensor<1, dim,
Tensor<rank_ - 1, dim, OtherNumber>>()
const
1299# ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
1300template <
int rank_,
int dim,
typename Number>
1304 for (
unsigned int i = 0; i < dim; ++i)
1310template <
int rank_,
int dim,
typename Number>
1314 for (
unsigned int i = 0; i < dim; ++i)
1322 namespace TensorSubscriptor
1324 template <
typename ArrayElementType,
int dim>
1327 subscript(ArrayElementType *
values,
1328 const unsigned int i,
1329 std::integral_constant<int, dim>)
1332# ifndef __CUDA_ARCH__
1343 template <
typename ArrayElementType>
1344 struct Uninitialized
1346 static ArrayElementType
value;
1349 template <
typename Type>
1350 Type Uninitialized<Type>::value;
1352 template <
typename ArrayElementType>
1355 subscript(ArrayElementType *,
1357 std::integral_constant<int, 0>)
1360# ifndef __CUDA_ARCH__
1364 "Cannot access elements of an object of type Tensor<rank,0,Number>."));
1366 return Uninitialized<ArrayElementType>::value;
1372template <
int rank_,
int dim,
typename Number>
1377 return ::internal::TensorSubscriptor::subscript(
1378 values, i, std::integral_constant<int, dim>());
1382template <
int rank_,
int dim,
typename Number>
1387# ifndef DEAL_II_COMPILER_CUDA_AWARE
1395template <
int rank_,
int dim,
typename Number>
1400# ifndef DEAL_II_COMPILER_CUDA_AWARE
1402 ExcMessage(
"Cannot access an object of type Tensor<rank_,0,Number>"));
1405 return TensorAccessors::extract<rank_>(*
this, indices);
1410template <
int rank_,
int dim,
typename Number>
1414# ifndef DEAL_II_COMPILER_CUDA_AWARE
1416 ExcMessage(
"Cannot access an object of type Tensor<rank_,0,Number>"));
1419 return TensorAccessors::extract<rank_>(*
this, indices);
1424template <
int rank_,
int dim,
typename Number>
1428 return std::addressof(
1429 this->
operator[](this->unrolled_to_component_indices(0)));
1434template <
int rank_,
int dim,
typename Number>
1435inline const Number *
1438 return std::addressof(
1439 this->
operator[](this->unrolled_to_component_indices(0)));
1444template <
int rank_,
int dim,
typename Number>
1448 return begin_raw() + n_independent_components;
1453template <
int rank_,
int dim,
typename Number>
1454inline const Number *
1457 return begin_raw() + n_independent_components;
1462template <
int rank_,
int dim,
typename Number>
1463template <
typename OtherNumber>
1469 for (
unsigned int i = 0; i < dim; ++i)
1476template <
int rank_,
int dim,
typename Number>
1483 for (
unsigned int i = 0; i < dim; ++i)
1489# ifdef DEAL_II_DELETED_MOVE_CONSTRUCTOR_BUG
1490template <
int rank_,
int dim,
typename Number>
1494 for (
unsigned int i = 0; i < dim; ++i)
1501template <
int rank_,
int dim,
typename Number>
1506 for (
unsigned int i = 0; i < dim; ++i)
1513template <
int rank_,
int dim,
typename Number>
1514template <
typename OtherNumber>
1515constexpr inline bool
1519 for (
unsigned int i = 0; i < dim; ++i)
1533constexpr inline bool
1540template <
int rank_,
int dim,
typename Number>
1541template <
typename OtherNumber>
1546 return !((*this) == p);
1550template <
int rank_,
int dim,
typename Number>
1551template <
typename OtherNumber>
1557 for (
unsigned int i = 0; i < dim; ++i)
1563template <
int rank_,
int dim,
typename Number>
1564template <
typename OtherNumber>
1570 for (
unsigned int i = 0; i < dim; ++i)
1576template <
int rank_,
int dim,
typename Number>
1577template <
typename OtherNumber>
1582 for (
unsigned int i = 0; i < dim; ++i)
1590 namespace TensorImplementation
1595 typename OtherNumber,
1596 typename std::enable_if<
1599 !std::is_same<Number, Differentiation::SD::Expression>::value,
1603 const OtherNumber &factor)
1605 const Number inverse_factor = Number(1.) / factor;
1607 for (
unsigned int d = 0;
d < dim; ++
d)
1608 t[
d] *= inverse_factor;
1615 typename OtherNumber,
1616 typename std::enable_if<
1619 std::is_same<Number, Differentiation::SD::Expression>::value,
1623 const OtherNumber &factor)
1626 for (
unsigned int d = 0;
d < dim; ++
d)
1633template <
int rank_,
int dim,
typename Number>
1634template <
typename OtherNumber>
1644template <
int rank_,
int dim,
typename Number>
1651 for (
unsigned int i = 0; i < dim; ++i)
1658template <
int rank_,
int dim,
typename Number>
1666template <
int rank_,
int dim,
typename Number>
1673 for (
unsigned int i = 0; i < dim; ++i)
1674 s +=
values[i].norm_square();
1680template <
int rank_,
int dim,
typename Number>
1681template <
typename OtherNumber>
1686 (Utilities::fixed_power<rank_, unsigned int>(dim)));
1688 unsigned int index = 0;
1689 unroll_recursion(result, index);
1693template <
int rank_,
int dim,
typename Number>
1694template <
typename OtherNumber>
1697 unsigned int & index)
const
1699 for (
unsigned int i = 0; i < dim; ++i)
1700 values[i].unroll_recursion(result, index);
1704template <
int rank_,
int dim,
typename Number>
1705constexpr inline unsigned int
1709 unsigned int index = 0;
1710 for (
int r = 0; r < rank_; ++r)
1711 index = index * dim + indices[r];
1724 inline constexpr unsigned int
1725 mod(
const unsigned int x)
1732 mod<0>(
const unsigned int x)
1739 inline constexpr unsigned int
1740 div(
const unsigned int x)
1747 div<0>(
const unsigned int x)
1757template <
int rank_,
int dim,
typename Number>
1765 unsigned int remainder = i;
1766 for (
int r = rank_ - 1; r >= 0; --r)
1768 indices[r] = internal::mod<dim>(remainder);
1769 remainder = internal::div<dim>(remainder);
1777template <
int rank_,
int dim,
typename Number>
1778constexpr inline void
1781 for (
unsigned int i = 0; i < dim; ++i)
1786template <
int rank_,
int dim,
typename Number>
1787constexpr std::size_t
1794template <
int rank_,
int dim,
typename Number>
1795template <
class Archive>
1803template <
int rank_,
int dim,
typename Number>
1822template <
int rank_,
int dim,
typename Number>
1823inline std::ostream &
1826 for (
unsigned int i = 0; i < dim; ++i)
1843template <
int dim,
typename Number>
1844inline std::ostream &
1847 out << static_cast<const Number &>(p);
1869template <
int dim,
typename Number,
typename Other>
1874 return object *
static_cast<const Number &
>(t);
1889template <
int dim,
typename Number,
typename Other>
1894 return static_cast<const Number &
>(t) *
object;
1909template <
int dim,
typename Number,
typename OtherNumber>
1915 return static_cast<const Number &
>(src1) *
1916 static_cast<const OtherNumber &
>(src2);
1927template <
int dim,
typename Number,
typename OtherNumber>
1935 return static_cast<const Number &
>(t) / factor;
1946template <
int dim,
typename Number,
typename OtherNumber>
1952 return static_cast<const Number &
>(p) +
static_cast<const OtherNumber &
>(q);
1963template <
int dim,
typename Number,
typename OtherNumber>
1969 return static_cast<const Number &
>(p) -
static_cast<const OtherNumber &
>(q);
1985template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1995 for (
unsigned int d = 0;
d < dim; ++
d)
1996 tt[
d] = t[
d] * factor;
2013template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2028 namespace TensorImplementation
2033 typename OtherNumber,
2034 typename std::enable_if<
2041 const OtherNumber & factor)
2044 const Number inverse_factor = Number(1.) / factor;
2046 for (
unsigned int d = 0;
d < dim; ++
d)
2047 tt[
d] = t[
d] * inverse_factor;
2055 typename OtherNumber,
2056 typename std::enable_if<
2063 const OtherNumber & factor)
2067 for (
unsigned int d = 0;
d < dim; ++
d)
2068 tt[
d] = t[
d] / factor;
2084template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2105template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2113 for (
unsigned int i = 0; i < dim; ++i)
2129template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2137 for (
unsigned int i = 0; i < dim; ++i)
2149template <
int dim,
typename Number,
typename OtherNumber>
2178template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2186 for (
unsigned int i = 0; i < dim; ++i)
2222template <
int rank_1,
2226 typename OtherNumber,
2227 typename =
typename std::enable_if<rank_1 >= 1 && rank_2 >= 1>::type>
2229 typename Tensor<rank_1 + rank_2 - 2,
2235 typename Tensor<rank_1 + rank_2 - 2,
2240 TensorAccessors::internal::
2241 ReorderedIndexView<0, rank_2, const Tensor<rank_2, dim, OtherNumber>>
2242 reordered = TensorAccessors::reordered_index_view<0, rank_2>(src2);
2243 TensorAccessors::contract<1, rank_1, rank_2, dim>(result, src1, reordered);
2277template <
int index_1,
2283 typename OtherNumber>
2285 typename Tensor<rank_1 + rank_2 - 2,
2291 Assert(0 <= index_1 && index_1 < rank_1,
2293 "The specified index_1 must lie within the range [0,rank_1)"));
2294 Assert(0 <= index_2 && index_2 < rank_2,
2296 "The specified index_2 must lie within the range [0,rank_2)"));
2303 reord_01 = reordered_index_view<index_1, rank_1>(src1);
2309 reord_02 = reordered_index_view<index_2, rank_2>(src2);
2311 typename Tensor<rank_1 + rank_2 - 2,
2315 TensorAccessors::contract<1, rank_1, rank_2, dim>(result, reord_01, reord_02);
2350template <
int index_1,
2358 typename OtherNumber>
2360 typename Tensor<rank_1 + rank_2 - 4,
2366 Assert(0 <= index_1 && index_1 < rank_1,
2368 "The specified index_1 must lie within the range [0,rank_1)"));
2369 Assert(0 <= index_3 && index_3 < rank_1,
2371 "The specified index_3 must lie within the range [0,rank_1)"));
2372 Assert(index_1 != index_3,
2373 ExcMessage(
"index_1 and index_3 must not be the same"));
2374 Assert(0 <= index_2 && index_2 < rank_2,
2376 "The specified index_2 must lie within the range [0,rank_2)"));
2377 Assert(0 <= index_4 && index_4 < rank_2,
2379 "The specified index_4 must lie within the range [0,rank_2)"));
2380 Assert(index_2 != index_4,
2381 ExcMessage(
"index_2 and index_4 must not be the same"));
2388 reord_1 = TensorAccessors::reordered_index_view<index_1, rank_1>(src1);
2392 reord_2 = TensorAccessors::reordered_index_view<index_2, rank_2>(src2);
2398 (index_3 < index_1 ? index_3 : index_3 - 1),
2410 (index_4 < index_2 ? index_4 : index_4 - 1),
2418 typename Tensor<rank_1 + rank_2 - 4,
2422 TensorAccessors::contract<2, rank_1, rank_2, dim>(result, reord_3, reord_4);
2439template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2446 TensorAccessors::contract<rank, rank, rank, dim>(result, left, right);
2468template <
template <
int,
int,
typename>
class TensorT1,
2469 template <
int,
int,
typename>
class TensorT2,
2470 template <
int,
int,
typename>
class TensorT3,
2480 const TensorT2<rank_1 + rank_2, dim, T2> &middle,
2481 const TensorT3<rank_2, dim, T3> & right)
2485 return TensorAccessors::contract3<rank_1, rank_2, dim, return_type>(left,
2501template <
int rank_1,
2505 typename OtherNumber>
2511 typename Tensor<rank_1 + rank_2,
2515 TensorAccessors::contract<0, rank_1, rank_2, dim>(result, src1, src2);
2537template <
int dim,
typename Number>
2546 result[1] = -src[0];
2561template <
int dim,
typename Number1,
typename Number2>
2572 constexpr int s0 = 0 % dim;
2573 constexpr int s1 = 1 % dim;
2574 constexpr int s2 = 2 % dim;
2576 result[s0] = src1[s1] * src2[s2] - src1[s2] * src2[s1];
2577 result[s1] = src1[s2] * src2[s0] - src1[s0] * src2[s2];
2578 result[s2] = src1[s0] * src2[s1] - src1[s1] * src2[s0];
2596template <
int dim,
typename Number>
2604 for (
unsigned int k = 0; k < dim; ++k)
2606 Tensor<2, dim - 1, Number> minor;
2607 for (
unsigned int i = 0; i < dim - 1; ++i)
2608 for (
unsigned int j = 0; j < dim - 1; ++j)
2609 minor[i][j] = t[i][j < k ? j : j + 1];
2616 return ((dim % 2 == 0) ? 1. : -1.) * det;
2624template <
typename Number>
2636template <
typename Number>
2641 return t[0][0] * t[1][1] - t[1][0] * t[0][1];
2649template <
typename Number>
2660 return t[0][0] * C0 + t[0][1] * C1 + t[0][2] * C2;
2670template <
int dim,
typename Number>
2675 for (
unsigned int i = 1; i < dim; ++i)
2689template <
int dim,
typename Number>
2693 Number return_tensor[dim][dim];
2706template <
typename Number>
2714 return return_tensor;
2718template <
typename Number>
2725 1.0 / (t[0][0] * t[1][1] - t[1][0] * t[0][1]));
2726 return_tensor[0][0] = t[1][1];
2727 return_tensor[0][1] = -t[0][1];
2728 return_tensor[1][0] = -t[1][0];
2729 return_tensor[1][1] = t[0][0];
2730 return_tensor *= inv_det_t;
2732 return return_tensor;
2736template <
typename Number>
2761 1.0 / (t[0][0] * return_tensor[0][0] + t[0][1] * return_tensor[1][0] +
2762 t[0][2] * return_tensor[2][0]));
2763 return_tensor *= inv_det_t;
2765 return return_tensor;
2776template <
int dim,
typename Number>
2781 for (
unsigned int i = 0; i < dim; ++i)
2784 for (
unsigned int j = i + 1; j < dim; ++j)
2807template <
int dim,
typename Number>
2828template <
int dim,
typename Number>
2899template <
int dim,
typename Number>
2911template <
int dim,
typename Number>
2916 for (
unsigned int j = 0; j < dim; ++j)
2919 for (
unsigned int i = 0; i < dim; ++i)
2937template <
int dim,
typename Number>
2942 for (
unsigned int i = 0; i < dim; ++i)
2945 for (
unsigned int j = 0; j < dim; ++j)
2961# ifdef DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING
2970 for (
unsigned int j = 0; j < dim; ++j)
2973 for (
unsigned int i = 0; i < dim; ++i)
2988 for (
unsigned int i = 0; i < dim; ++i)
2991 for (
unsigned int j = 0; j < dim; ++j)
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &t)
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &t)
constexpr Tensor & operator*=(const OtherNumber &factor)
void serialize(Archive &ar, const unsigned int version)
constexpr Tensor & operator/=(const OtherNumber &factor)
const Number * end_raw() const
constexpr Tensor & operator-=(const Tensor< 0, dim, OtherNumber > &rhs)
constexpr Tensor(const Tensor< 0, dim, OtherNumber > &initializer)
constexpr Tensor(const OtherNumber &initializer)
constexpr Tensor & operator=(const OtherNumber &d)
constexpr real_type norm_square() const
constexpr bool operator!=(const Tensor< 0, dim, OtherNumber > &rhs) const
constexpr Tensor & operator=(const Tensor< 0, dim, OtherNumber > &rhs)
const Number * begin_raw() const
constexpr Tensor & operator+=(const Tensor< 0, dim, OtherNumber > &rhs)
constexpr bool operator==(const Tensor< 0, dim, OtherNumber > &rhs) const
typename numbers::NumberTraits< Number >::real_type real_type
constexpr Tensor operator-() const
void unroll_recursion(Vector< OtherNumber > &result, unsigned int &start_index) const
constexpr Tensor & operator/=(const OtherNumber &factor)
constexpr Tensor(const ArrayView< ElementType, MemorySpace > &initializer)
constexpr Tensor< 2, dim, Number > cofactor(const Tensor< 2, dim, Number > &t)
constexpr bool operator==(const Tensor< rank_, dim, OtherNumber > &) const
constexpr Tensor< rank, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/(const Tensor< rank, dim, Number > &t, const OtherNumber &factor)
constexpr Tensor(const Tensor< 1, dim, Tensor< rank_ - 1, dim, OtherNumber > > &initializer)
constexpr 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)
constexpr Tensor< 2, dim, Number > adjugate(const Tensor< 2, dim, Number > &t)
typename Tensor< rank_ - 1, dim, Number >::array_type[(dim !=0) ? dim :1] array_type
constexpr const Number & operator[](const TableIndices< rank_ > &indices) const
constexpr ProductType< T1, typenameProductType< 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)
constexpr Tensor< 2, dim, Number > transpose(const Tensor< 2, dim, Number > &t)
constexpr Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > operator-(const Tensor< 0, dim, Number > &p, const Tensor< 0, dim, OtherNumber > &q)
constexpr Tensor< 0, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator/(const Tensor< 0, dim, Number > &t, const OtherNumber &factor)
static constexpr unsigned int rank
constexpr Number determinant(const Tensor< 2, dim, Number > &t)
constexpr Tensor(const Tensor< rank_, dim, OtherNumber > &initializer)
numbers::NumberTraits< Number >::real_type norm() const
constexpr Tensor & operator-=(const Tensor< rank_, dim, OtherNumber > &)
void unroll(Vector< OtherNumber > &result) const
static constexpr unsigned int component_to_unrolled_index(const TableIndices< rank_ > &indices)
constexpr ProductType< Number, OtherNumber >::type operator*(const Tensor< 0, dim, Number > &src1, const Tensor< 0, dim, OtherNumber > &src2)
const Number * begin_raw() const
constexpr bool operator!=(const Tensor< rank_, dim, OtherNumber > &) const
constexpr ProductType< Number, OtherNumber >::type scalar_product(const Tensor< rank, dim, Number > &left, const Tensor< rank, dim, OtherNumber > &right)
constexpr Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > schur_product(const Tensor< 0, dim, Number > &src1, const Tensor< 0, dim, OtherNumber > &src2)
constexpr value_type & operator[](const unsigned int i)
typename Tensor< rank_ - 1, dim, Number >::tensor_type value_type
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > operator+(const Tensor< rank, dim, Number > &p, const Tensor< rank, dim, OtherNumber > &q)
Number linfty_norm(const Tensor< 2, dim, Number > &t)
constexpr ProductType< Other, Number >::type operator*(const Other &object, const Tensor< 0, dim, Number > &t)
constexpr Tensor< rank_1+rank_2-4, dim, typenameProductType< 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 > project_onto_orthogonal_tensors(const Tensor< 2, dim, Number > &A)
Number l1_norm(const Tensor< 2, dim, Number > &t)
constexpr Number trace(const Tensor< 2, dim, Number > &d)
static constexpr unsigned int dimension
static constexpr TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
constexpr Tensor< rank, dim, typename ProductType< typename EnableIfScalar< Number >::type, OtherNumber >::type > operator*(const Number &factor, const Tensor< rank, dim, OtherNumber > &t)
static constexpr std::size_t memory_consumption()
constexpr Tensor & operator=(const Number &d)
OtherNumber::type::tensor_type operator*(const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
constexpr Number determinant(const Tensor< 2, 1, Number > &t)
constexpr Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > operator+(const Tensor< 0, dim, Number > &p, const Tensor< 0, dim, OtherNumber > &q)
constexpr Tensor & operator+=(const Tensor< rank_, dim, OtherNumber > &)
const Number * end_raw() const
void unroll_recursion(Vector< OtherNumber > &result, unsigned int &start_index) const
constexpr Number & operator[](const TableIndices< rank_ > &indices)
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > schur_product(const Tensor< rank, dim, Number > &src1, const Tensor< rank, dim, OtherNumber > &src2)
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > operator-(const Tensor< rank, dim, Number > &p, const Tensor< rank, dim, OtherNumber > &q)
Tensor< rank_, dim, Number > tensor_type
constexpr ProductType< Number, Other >::type operator*(const Tensor< 0, dim, Number > &t, const Other &object)
constexpr Tensor(const ArrayLike &initializer, std::index_sequence< Indices... >)
constexpr Number determinant(const Tensor< 2, 2, Number > &t)
constexpr Tensor< rank_1+rank_2-2, dim, typenameProductType< Number, OtherNumber >::type >::tensor_type contract(const Tensor< rank_1, dim, Number > &src1, const Tensor< rank_2, dim, OtherNumber > &src2)
constexpr Tensor< rank, dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const Tensor< rank, dim, Number > &t, const OtherNumber &factor)
constexpr Tensor< 2, dim, Number > invert(const Tensor< 2, dim, Number > &)
void serialize(Archive &ar, const unsigned int version)
constexpr Tensor & operator*=(const OtherNumber &factor)
constexpr const value_type & operator[](const unsigned int i) const
constexpr Number determinant(const Tensor< 2, 3, Number > &t)
constexpr Tensor(const array_type &initializer)
static constexpr unsigned int n_independent_components
constexpr Tensor operator-() const
Tensor< rank_ - 1, dim, Number > values[(dim !=0) ? dim :1]
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
constexpr Tensor & operator=(const Tensor< rank_, dim, OtherNumber > &rhs)
VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &y)
#define DEAL_II_ALWAYS_INLINE
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Expression fabs(const Expression &x)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr internal::ReorderedIndexView< index, rank, T > reordered_index_view(T &t)
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > division_operator(const Tensor< rank, dim, Number > &t, const OtherNumber &factor)
constexpr bool value_is_zero(const Number &value)
constexpr bool values_are_equal(const Number1 &value_1, const Number2 &value_2)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
#define DEAL_II_CUDA_HOST_DEV
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type
static constexpr const T & value(const T &t)
decltype(std::declval< T >() *std::declval< U >()) type
static constexpr std::enable_if< std::is_same< Dummy, number >::value &&is_cuda_compatible< Dummy >::value, real_type >::type abs_square(const number &x)
static real_type abs(const number &x)