|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_tensor_h
17 #define dealii_tensor_h
31 #ifdef DEAL_II_WITH_ADOLC
32 # include <adolc/adouble.h>
44 template <
int dim,
typename Number>
46 template <
int rank_,
int dim,
typename Number =
double>
48 template <
typename Number>
50 template <
typename number>
92 template <
int dim,
typename Number>
96 static_assert(dim >= 0,
97 "Tensors must have a dimension greater than or equal to one.");
112 static constexpr
unsigned int rank = 0;
157 template <
typename OtherNumber>
166 template <
typename OtherNumber>
168 Tensor(
const OtherNumber &initializer);
223 template <
typename OtherNumber>
227 #ifdef __INTEL_COMPILER
246 template <
typename OtherNumber>
253 template <
typename OtherNumber>
260 template <
typename OtherNumber>
269 template <
typename OtherNumber>
278 template <
typename OtherNumber>
287 template <
typename OtherNumber>
296 template <
typename OtherNumber>
344 template <
class Archive>
346 serialize(Archive &ar,
const unsigned int version);
363 template <
typename OtherNumber>
366 unsigned int & start_index)
const;
369 template <
int,
int,
typename>
449 template <
int rank_,
int dim,
typename Number>
453 static_assert(rank_ >= 0,
454 "Tensors must have a rank greater than or equal to one.");
455 static_assert(dim >= 0,
456 "Tensors must have a dimension greater than or equal to one.");
470 static constexpr
unsigned int rank = rank_;
522 template <
typename OtherNumber>
529 template <
typename OtherNumber>
536 template <
typename OtherNumber>
538 operator Tensor<1, dim,
Tensor<rank_ - 1, dim, OtherNumber>>()
const;
598 template <
typename OtherNumber>
614 template <
typename OtherNumber>
621 template <
typename OtherNumber>
630 template <
typename OtherNumber>
639 template <
typename OtherNumber>
649 template <
typename OtherNumber>
658 template <
typename OtherNumber>
713 template <
typename OtherNumber>
715 unroll(Vector<OtherNumber> &result)
const;
736 static constexpr std::size_t
743 template <
class Archive>
745 serialize(Archive &ar,
const unsigned int version);
764 template <
typename OtherNumber>
767 unsigned int & start_index)
const;
775 template <
typename ArrayLike, std::size_t... Indices>
777 Tensor(
const ArrayLike &initializer, std_cxx14::index_sequence<Indices...>);
780 template <
int,
int,
typename>
797 template <
int rank,
int dim,
typename T,
typename U>
798 struct ProductTypeImpl<
Tensor<rank, dim,
T>, std::complex<U>>
804 template <
int rank,
int dim,
typename T,
typename U>
805 struct ProductTypeImpl<
Tensor<rank, dim, std::complex<T>>, std::complex<U>>
811 template <
typename T,
int rank,
int dim,
typename U>
812 struct ProductTypeImpl<std::complex<T>,
Tensor<rank, dim, U>>
818 template <
int rank,
int dim,
typename T,
typename U>
819 struct ProductTypeImpl<std::complex<T>,
Tensor<rank, dim, std::complex<U>>>
830 template <
int rank,
int dim,
typename T>
831 struct NumberType<
Tensor<rank, dim,
T>>
853 template <
int dim,
typename Number>
863 template <
int dim,
typename Number>
864 template <
typename OtherNumber>
872 template <
int dim,
typename Number>
873 template <
typename OtherNumber>
881 template <
int dim,
typename Number>
885 return std::addressof(
value);
890 template <
int dim,
typename Number>
891 inline const Number *
894 return std::addressof(
value);
899 template <
int dim,
typename Number>
908 template <
int dim,
typename Number>
917 template <
int dim,
typename Number>
922 # ifndef __CUDA_ARCH__
924 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
930 template <
int dim,
typename Number>
935 # ifndef __CUDA_ARCH__
937 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
943 template <
int dim,
typename Number>
944 template <
typename OtherNumber>
954 # ifdef __INTEL_COMPILER
955 template <
int dim,
typename Number>
966 template <
int dim,
typename Number>
967 template <
typename OtherNumber>
977 template <
int dim,
typename Number>
978 template <
typename OtherNumber>
982 # if defined(DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING)
986 "The Tensor equality operator for ADOL-C taped numbers has not yet "
987 "been extended to support advanced branching."));
994 template <
int dim,
typename Number>
995 template <
typename OtherNumber>
999 return !((*this) == p);
1003 template <
int dim,
typename Number>
1004 template <
typename OtherNumber>
1014 template <
int dim,
typename Number>
1015 template <
typename OtherNumber>
1028 namespace ComplexWorkaround
1030 template <
typename Number,
typename OtherNumber>
1032 multiply_assign_scalar(Number &val,
const OtherNumber &s)
1037 # ifdef __CUDA_ARCH__
1038 template <
typename Number,
typename OtherNumber>
1040 multiply_assign_scalar(std::complex<Number> &,
const OtherNumber &)
1042 printf(
"This function is not implemented for std::complex<Number>!\n");
1050 template <
int dim,
typename Number>
1051 template <
typename OtherNumber>
1056 internal::ComplexWorkaround::multiply_assign_scalar(
value, s);
1062 template <
int dim,
typename Number>
1063 template <
typename OtherNumber>
1072 template <
int dim,
typename Number>
1080 template <
int dim,
typename Number>
1085 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
1090 template <
int dim,
typename Number>
1096 # ifndef __CUDA_ARCH__
1098 ExcMessage(
"Cannot access an object of type Tensor<0,0,Number>"));
1104 template <
int dim,
typename Number>
1105 template <
typename OtherNumber>
1108 unsigned int & index)
const
1111 ExcMessage(
"Cannot unroll an object of type Tensor<0,0,Number>"));
1112 result[index] =
value;
1117 template <
int dim,
typename Number>
1127 template <
int dim,
typename Number>
1128 template <
class Archive>
1136 template <
int dim,
typename Number>
1142 template <
int rank_,
int dim,
typename Number>
1143 template <
typename ArrayLike, std::size_t... indices>
1146 std_cxx14::index_sequence<indices...>)
1147 : values{
Tensor<rank_ - 1, dim, Number>(initializer[indices])...}
1149 static_assert(
sizeof...(indices) == dim,
1150 "dim should match the number of indices");
1154 template <
int rank_,
int dim,
typename Number>
1161 template <
int rank_,
int dim,
typename Number>
1162 template <
typename OtherNumber>
1170 template <
int rank_,
int dim,
typename Number>
1171 template <
typename OtherNumber>
1179 template <
int rank_,
int dim,
typename Number>
1180 template <
typename OtherNumber>
1182 operator
Tensor<1, dim,
Tensor<rank_ - 1, dim, OtherNumber>>()
const
1184 return Tensor<1, dim,
Tensor<rank_ - 1, dim, Number>>(values);
1191 namespace TensorSubscriptor
1193 template <
typename ArrayElementType,
int dim>
1196 subscript(ArrayElementType * values,
1197 const unsigned int i,
1198 std::integral_constant<int, dim>)
1201 # ifndef __CUDA_ARCH__
1212 template <
typename ArrayElementType>
1213 struct Uninitialized
1215 static ArrayElementType
value;
1218 template <
typename Type>
1221 template <
typename ArrayElementType>
1224 subscript(ArrayElementType *,
1226 std::integral_constant<int, 0>)
1229 # ifndef __CUDA_ARCH__
1233 "Cannot access elements of an object of type Tensor<rank,0,Number>."));
1241 template <
int rank_,
int dim,
typename Number>
1246 return ::internal::TensorSubscriptor::subscript(
1247 values, i, std::integral_constant<int, dim>());
1251 template <
int rank_,
int dim,
typename Number>
1260 template <
int rank_,
int dim,
typename Number>
1266 ExcMessage(
"Cannot access an object of type Tensor<rank_,0,Number>"));
1268 return TensorAccessors::extract<rank_>(*
this, indices);
1273 template <
int rank_,
int dim,
typename Number>
1278 ExcMessage(
"Cannot access an object of type Tensor<rank_,0,Number>"));
1280 return TensorAccessors::extract<rank_>(*
this, indices);
1285 template <
int rank_,
int dim,
typename Number>
1289 return std::addressof(
1290 this->
operator[](this->unrolled_to_component_indices(0)));
1295 template <
int rank_,
int dim,
typename Number>
1296 inline const Number *
1299 return std::addressof(
1300 this->
operator[](this->unrolled_to_component_indices(0)));
1305 template <
int rank_,
int dim,
typename Number>
1309 return begin_raw() + n_independent_components;
1314 template <
int rank_,
int dim,
typename Number>
1315 inline const Number *
1318 return begin_raw() + n_independent_components;
1323 template <
int rank_,
int dim,
typename Number>
1324 template <
typename OtherNumber>
1330 for (
unsigned int i = 0; i < dim; ++i)
1336 template <
int rank_,
int dim,
typename Number>
1341 ExcMessage(
"Only assignment with zero is allowed"));
1344 for (
unsigned int i = 0; i < dim; ++i)
1350 template <
int rank_,
int dim,
typename Number>
1351 template <
typename OtherNumber>
1356 for (
unsigned int i = 0; i < dim; ++i)
1357 if (values[i] != p.
values[i])
1377 template <
int rank_,
int dim,
typename Number>
1378 template <
typename OtherNumber>
1383 return !((*this) == p);
1387 template <
int rank_,
int dim,
typename Number>
1388 template <
typename OtherNumber>
1394 for (
unsigned int i = 0; i < dim; ++i)
1395 values[i] += p.
values[i];
1400 template <
int rank_,
int dim,
typename Number>
1401 template <
typename OtherNumber>
1407 for (
unsigned int i = 0; i < dim; ++i)
1408 values[i] -= p.
values[i];
1413 template <
int rank_,
int dim,
typename Number>
1414 template <
typename OtherNumber>
1419 for (
unsigned int i = 0; i < dim; ++i)
1427 namespace TensorImplementation
1432 typename OtherNumber,
1433 typename std::enable_if<
1440 const OtherNumber &factor)
1442 const Number inverse_factor = Number(1.) / factor;
1444 for (
unsigned int d = 0;
d < dim; ++
d)
1445 t[
d] *= inverse_factor;
1452 typename OtherNumber,
1453 typename std::enable_if<
1460 const OtherNumber &factor)
1463 for (
unsigned int d = 0;
d < dim; ++
d)
1470 template <
int rank_,
int dim,
typename Number>
1471 template <
typename OtherNumber>
1481 template <
int rank_,
int dim,
typename Number>
1488 for (
unsigned int i = 0; i < dim; ++i)
1489 tmp.
values[i] = -values[i];
1495 template <
int rank_,
int dim,
typename Number>
1503 template <
int rank_,
int dim,
typename Number>
1510 for (
unsigned int i = 0; i < dim; ++i)
1511 s += values[i].norm_square();
1517 template <
int rank_,
int dim,
typename Number>
1518 template <
typename OtherNumber>
1523 (Utilities::fixed_power<rank_, unsigned int>(dim)));
1525 unsigned int index = 0;
1526 unroll_recursion(result, index);
1530 template <
int rank_,
int dim,
typename Number>
1531 template <
typename OtherNumber>
1534 unsigned int & index)
const
1536 for (
unsigned int i = 0; i < dim; ++i)
1537 values[i].unroll_recursion(result, index);
1541 template <
int rank_,
int dim,
typename Number>
1546 unsigned int index = 0;
1547 for (
int r = 0; r < rank_; ++r)
1548 index = index * dim + indices[r];
1562 mod(
const unsigned int x)
1569 mod<0>(
const unsigned int x)
1577 div(
const unsigned int x)
1584 div<0>(
const unsigned int x)
1594 template <
int rank_,
int dim,
typename Number>
1602 unsigned int remainder = i;
1603 for (
int r = rank_ - 1; r >= 0; --r)
1605 indices[r] = internal::mod<dim>(remainder);
1606 remainder = internal::div<dim>(remainder);
1614 template <
int rank_,
int dim,
typename Number>
1618 for (
unsigned int i = 0; i < dim; ++i)
1623 template <
int rank_,
int dim,
typename Number>
1624 constexpr std::size_t
1631 template <
int rank_,
int dim,
typename Number>
1632 template <
class Archive>
1640 template <
int rank_,
int dim,
typename Number>
1659 template <
int rank_,
int dim,
typename Number>
1660 inline std::ostream &
1663 for (
unsigned int i = 0; i < dim; ++i)
1680 template <
int dim,
typename Number>
1681 inline std::ostream &
1684 out << static_cast<const Number &>(p);
1706 template <
int dim,
typename Number,
typename Other>
1711 return object *
static_cast<const Number &
>(t);
1726 template <
int dim,
typename Number,
typename Other>
1731 return static_cast<const Number &
>(t) *
object;
1746 template <
int dim,
typename Number,
typename OtherNumber>
1752 return static_cast<const Number &
>(src1) *
1753 static_cast<const OtherNumber &
>(src2);
1764 template <
int dim,
typename Number,
typename OtherNumber>
1772 return static_cast<const Number &
>(t) / factor;
1783 template <
int dim,
typename Number,
typename OtherNumber>
1789 return static_cast<const Number &
>(p) +
static_cast<const OtherNumber &
>(q);
1800 template <
int dim,
typename Number,
typename OtherNumber>
1806 return static_cast<const Number &
>(p) -
static_cast<const OtherNumber &
>(q);
1822 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1832 for (
unsigned int d = 0;
d < dim; ++
d)
1833 tt[
d] = t[
d] * factor;
1850 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1865 namespace TensorImplementation
1870 typename OtherNumber,
1871 typename std::enable_if<
1878 const OtherNumber & factor)
1881 const Number inverse_factor = Number(1.) / factor;
1883 for (
unsigned int d = 0;
d < dim; ++
d)
1884 tt[
d] = t[
d] * inverse_factor;
1892 typename OtherNumber,
1893 typename std::enable_if<
1900 const OtherNumber & factor)
1904 for (
unsigned int d = 0;
d < dim; ++
d)
1905 tt[
d] = t[
d] / factor;
1921 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1942 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1950 for (
unsigned int i = 0; i < dim; ++i)
1966 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
1974 for (
unsigned int i = 0; i < dim; ++i)
1986 template <
int dim,
typename Number,
typename OtherNumber>
2015 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2023 for (
unsigned int i = 0; i < dim; ++i)
2060 template <
int rank_1,
2064 typename OtherNumber,
2065 typename =
typename std::enable_if<rank_1 >= 1 && rank_2 >= 1>::type>
2067 typename Tensor<rank_1 + rank_2 - 2,
2073 typename Tensor<rank_1 + rank_2 - 2,
2078 TensorAccessors::internal::
2079 ReorderedIndexView<0, rank_2, const Tensor<rank_2, dim, OtherNumber>>
2080 reordered = TensorAccessors::reordered_index_view<0, rank_2>(src2);
2081 TensorAccessors::contract<1, rank_1, rank_2, dim>(result, src1, reordered);
2116 template <
int index_1,
2122 typename OtherNumber>
2124 typename Tensor<rank_1 + rank_2 - 2,
2130 Assert(0 <= index_1 && index_1 < rank_1,
2132 "The specified index_1 must lie within the range [0,rank_1)"));
2133 Assert(0 <= index_2 && index_2 < rank_2,
2135 "The specified index_2 must lie within the range [0,rank_2)"));
2142 reord_01 = reordered_index_view<index_1, rank_1>(src1);
2146 reord_02 = reordered_index_view<index_2, rank_2>(src2);
2148 typename Tensor<rank_1 + rank_2 - 2,
2152 TensorAccessors::contract<1, rank_1, rank_2, dim>(result, reord_01, reord_02);
2188 template <
int index_1,
2196 typename OtherNumber>
2198 typename Tensor<rank_1 + rank_2 - 4,
2204 Assert(0 <= index_1 && index_1 < rank_1,
2206 "The specified index_1 must lie within the range [0,rank_1)"));
2207 Assert(0 <= index_3 && index_3 < rank_1,
2209 "The specified index_3 must lie within the range [0,rank_1)"));
2210 Assert(index_1 != index_3,
2211 ExcMessage(
"index_1 and index_3 must not be the same"));
2212 Assert(0 <= index_2 && index_2 < rank_2,
2214 "The specified index_2 must lie within the range [0,rank_2)"));
2215 Assert(0 <= index_4 && index_4 < rank_2,
2217 "The specified index_4 must lie within the range [0,rank_2)"));
2218 Assert(index_2 != index_4,
2219 ExcMessage(
"index_2 and index_4 must not be the same"));
2226 reord_1 = TensorAccessors::reordered_index_view<index_1, rank_1>(src1);
2230 reord_2 = TensorAccessors::reordered_index_view<index_2, rank_2>(src2);
2236 (index_3 < index_1 ? index_3 : index_3 - 1),
2248 (index_4 < index_2 ? index_4 : index_4 - 1),
2256 typename Tensor<rank_1 + rank_2 - 4,
2260 TensorAccessors::contract<2, rank_1, rank_2, dim>(result, reord_3, reord_4);
2278 template <
int rank,
int dim,
typename Number,
typename OtherNumber>
2285 TensorAccessors::contract<rank, rank, rank, dim>(result, left, right);
2308 template <
template <
int,
int,
typename>
class TensorT1,
2309 template <
int,
int,
typename>
class TensorT2,
2310 template <
int,
int,
typename>
class TensorT3,
2320 const TensorT2<rank_1 + rank_2, dim, T2> &middle,
2321 const TensorT3<rank_2, dim, T3> & right)
2325 return TensorAccessors::contract3<rank_1, rank_2, dim, return_type>(left,
2342 template <
int rank_1,
2346 typename OtherNumber>
2352 typename Tensor<rank_1 + rank_2,
2356 TensorAccessors::contract<0, rank_1, rank_2, dim>(result, src1, src2);
2379 template <
int dim,
typename Number>
2388 result[1] = -src[0];
2404 template <
int dim,
typename Number1,
typename Number2>
2415 constexpr
int s0 = 0 % dim;
2416 constexpr
int s1 = 1 % dim;
2417 constexpr
int s2 = 2 % dim;
2419 result[s0] = src1[s1] * src2[s2] - src1[s2] * src2[s1];
2420 result[s1] = src1[s2] * src2[s0] - src1[s0] * src2[s2];
2421 result[s2] = src1[s0] * src2[s1] - src1[s1] * src2[s0];
2440 template <
int dim,
typename Number>
2448 for (
unsigned int k = 0; k < dim; ++k)
2450 Tensor<2, dim - 1, Number> minor;
2451 for (
unsigned int i = 0; i < dim - 1; ++i)
2452 for (
unsigned int j = 0; j < dim - 1; ++j)
2453 minor[i][j] = t[i][j < k ? j : j + 1];
2460 return ((dim % 2 == 0) ? 1. : -1.) * det;
2468 template <
typename Number>
2483 template <
int dim,
typename Number>
2488 for (
unsigned int i = 1; i < dim; ++i)
2503 template <
int dim,
typename Number>
2507 Number return_tensor[dim][dim];
2520 template <
typename Number>
2528 return return_tensor;
2532 template <
typename Number>
2541 1.0 / (t[0][0] * t[1][1] - t[1][0] * t[0][1]));
2542 return_tensor[0][0] = t[1][1];
2543 return_tensor[0][1] = -t[0][1];
2544 return_tensor[1][0] = -t[1][0];
2545 return_tensor[1][1] = t[0][0];
2546 return_tensor *= inv_det_t;
2548 return return_tensor;
2552 template <
typename Number>
2565 1.0 / (t4 * t[2][2] - t6 * t[2][1] - t8 * t[2][2] +
2566 t00 * t[2][1] + t01 * t[1][2] - t04 * t[1][1]));
2575 return_tensor[1][1] =
2577 return_tensor[1][2] = t00 - t6;
2580 return_tensor[2][1] =
2583 return_tensor *= inv_det_t;
2585 return return_tensor;
2597 template <
int dim,
typename Number>
2602 for (
unsigned int i = 0; i < dim; ++i)
2605 for (
unsigned int j = i + 1; j < dim; ++j)
2629 template <
int dim,
typename Number>
2651 template <
int dim,
typename Number>
2669 template <
int dim,
typename Number>
2679 matrix.copy_from(tensor);
2690 matrix.copy_to(output_tensor);
2691 return output_tensor;
2703 template <
int dim,
typename Number>
2708 for (
unsigned int j = 0; j < dim; ++j)
2711 for (
unsigned int i = 0; i < dim; ++i)
2730 template <
int dim,
typename Number>
2735 for (
unsigned int i = 0; i < dim; ++i)
2738 for (
unsigned int j = 0; j < dim; ++j)
2754 # ifdef DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING
2763 for (
unsigned int j = 0; j < dim; ++j)
2766 for (
unsigned int i = 0; i < dim; ++i)
2781 for (
unsigned int i = 0; i < dim; ++i)
2784 for (
unsigned int j = 0; j < dim; ++j)
2793 # endif // DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING
constexpr bool operator==(const Tensor< rank_, dim, OtherNumber > &) const
void copy_from(const MatrixType &)
typename numbers::NumberTraits< Number >::real_type real_type
static constexpr unsigned int n_independent_components
constexpr Tensor< 2, dim, Number > adjugate(const Tensor< 2, dim, Number > &t)
constexpr Number trace(const Tensor< 2, dim, Number > &d)
constexpr Tensor & operator*=(const OtherNumber &factor)
constexpr Tensor operator-() const
static ::ExceptionBase & ExcNotImplemented()
constexpr 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)
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator+(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
static constexpr const T & value(const T &t)
constexpr value_type & operator[](const unsigned int i)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
constexpr Tensor< rank, dim, typename ProductType< Number, OtherNumber >::type > division_operator(const Tensor< rank, dim, Number > &t, const OtherNumber &factor)
const LAPACKFullMatrix< number > & get_svd_vt() const
#define AssertIndexRange(index, range)
static constexpr std::enable_if< std::is_same< Dummy, number >::value &&is_cuda_compatible< Dummy >::value, real_type >::type abs_square(const number &x)
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type
constexpr 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)
constexpr Tensor< 2, dim, Number > invert(const Tensor< 2, dim, Number > &)
constexpr ProductType< Number, OtherNumber >::type scalar_product(const Tensor< rank, dim, Number > &left, const Tensor< rank, dim, OtherNumber > &right)
#define DEAL_II_ALWAYS_INLINE
constexpr bool values_are_equal(const Number1 &value_1, const Number2 &value_2)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static constexpr unsigned int dimension
#define DEAL_II_CONSTEXPR
Expression fabs(const Expression &x)
constexpr Tensor< 0, dim, typename ProductType< Number, OtherNumber >::type > schur_product(const Tensor< 0, dim, Number > &src1, const Tensor< 0, dim, OtherNumber > &src2)
constexpr Number determinant(const Tensor< 2, dim, Number > &t)
static constexpr TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
constexpr Tensor & operator+=(const Tensor< rank_, dim, OtherNumber > &)
Point< dim, typename ProductType< Number, typename EnableIfScalar< OtherNumber >::type >::type > operator*(const OtherNumber) const
numbers::NumberTraits< Number >::real_type norm() const
static ::ExceptionBase & ExcMessage(std::string arg1)
constexpr Tensor & operator-=(const Tensor< rank_, dim, OtherNumber > &)
static constexpr std::size_t memory_consumption()
const LAPACKFullMatrix< number > & get_svd_u() const
#define DEAL_II_NAMESPACE_OPEN
void unroll(Vector< OtherNumber > &result) const
static constexpr unsigned int rank
@ matrix
Contents is actually a matrix.
Number l1_norm(const Tensor< 2, dim, Number > &t)
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
T sum(const T &t, const MPI_Comm &mpi_communicator)
#define AssertDimension(dim1, dim2)
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator-(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
static ::ExceptionBase & ExcInternalError()
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 OtherNumber &factor)
inline ::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &x)
#define Assert(cond, exc)
Tensor< rank_, dim, Number > tensor_type
constexpr bool operator!=(const Tensor< rank_, dim, OtherNumber > &) const
void unroll_recursion(Vector< OtherNumber > &result, unsigned int &start_index) const
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
std::enable_if< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typename ProductType< std::complex< T >, std::complex< U > >::type >::type operator/(const std::complex< T > &left, const std::complex< U > &right)
static real_type abs(const number &x)
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 internal::ReorderedIndexView< index, rank, T > reordered_index_view(T &t)
typename Tensor< rank_ - 1, dim, Tensor< 1, spacedim, VectorizedArrayType > >::tensor_type value_type
void mmult(LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const
constexpr Tensor & operator=(const Tensor< rank_, dim, OtherNumber > &rhs)
constexpr Tensor< 2, dim, Number > cofactor(const Tensor< 2, dim, Number > &t)
#define DEAL_II_CUDA_HOST_DEV
#define DEAL_II_NAMESPACE_CLOSE
static constexpr unsigned int component_to_unrolled_index(const TableIndices< rank_ > &indices)
constexpr Tensor< 2, dim, Number > transpose(const Tensor< 2, dim, Number > &t)
constexpr T1 contract3(const T2 &left, const T3 &middle, const T4 &right)
Tensor< 2, dim, Number > project_onto_orthogonal_tensors(const Tensor< 2, dim, Number > &tensor)
#define AssertThrow(cond, exc)
Tensor< rank_ - 1, dim, Number > values[(dim !=0) ? dim :1]
constexpr bool value_is_zero(const Number &value)
decltype(std::declval< T >() *std::declval< U >()) type
typename Tensor< rank_ - 1, dim, Tensor< 1, spacedim, VectorizedArrayType > >::array_type[(dim !=0) ? dim :1] array_type
T max(const T &t, const MPI_Comm &mpi_communicator)
Number linfty_norm(const Tensor< 2, dim, Number > &t)
void serialize(Archive &ar, const unsigned int version)