Reference documentation for deal.II version 9.0.0
tensor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_tensor_h
17 #define dealii_tensor_h
18 
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>
26 
27 #ifdef DEAL_II_WITH_ADOLC
28 #include <adolc/adouble.h> // Taped double
29 #endif
30 
31 #include <cmath>
32 #include <ostream>
33 #include <vector>
34 
35 
36 DEAL_II_NAMESPACE_OPEN
37 
38 // Forward declarations:
39 
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;
43 template <typename Number> class VectorizedArray;
44 
45 #ifndef DOXYGEN
46 // Overload invalid tensor types of negative rank that come up during
47 // overload resolution of operator* and related contraction variants.
48 template <int dim, typename Number>
49 class Tensor<-2, dim, Number>
50 {
51 };
52 
53 template <int dim, typename Number>
54 class Tensor<-1, dim, Number>
55 {
56 };
57 #endif /* DOXYGEN */
58 
59 
90 template <int dim, typename Number>
91 class Tensor<0,dim,Number>
92 {
93 public:
102  static const unsigned int dimension = dim;
103 
107  static const unsigned int rank = 0;
108 
112  static const unsigned int n_independent_components = 1;
113 
123 
128  typedef Number value_type;
129 
135  typedef Number array_type;
136 
142  DEAL_II_CUDA_HOST_DEV Tensor ();
143 
149  template <typename OtherNumber>
150  Tensor (const Tensor<0,dim,OtherNumber> &initializer);
151 
155  template <typename OtherNumber>
156  Tensor (const OtherNumber &initializer);
157 
161  Number *
162  begin_raw();
163 
167  const Number *
168  begin_raw() const;
169 
173  Number *
174  end_raw();
175 
180  const Number *
181  end_raw() const;
182 
192  DEAL_II_CUDA_HOST_DEV operator Number &();
193 
202  DEAL_II_CUDA_HOST_DEV operator const Number &() const;
203 
209  template <typename OtherNumber>
211 
212 #ifdef __INTEL_COMPILER
213 
220 #endif
221 
226  template <typename OtherNumber>
227  Tensor &operator = (const OtherNumber &d);
228 
232  template <typename OtherNumber>
233  bool operator == (const Tensor<0,dim,OtherNumber> &rhs) const;
234 
238  template <typename OtherNumber>
239  bool operator != (const Tensor<0,dim,OtherNumber> &rhs) const;
240 
244  template <typename OtherNumber>
246 
250  template <typename OtherNumber>
252 
258  template <typename OtherNumber>
259  DEAL_II_CUDA_HOST_DEV Tensor &operator *= (const OtherNumber &factor);
260 
264  template <typename OtherNumber>
265  Tensor &operator /= (const OtherNumber &factor);
266 
270  Tensor operator - () const;
271 
284  void clear ();
285 
291  real_type norm () const;
292 
299  DEAL_II_CUDA_HOST_DEV real_type norm_square () const;
300 
305  template <class Archive>
306  void serialize(Archive &ar, const unsigned int version);
307 
312  typedef Number tensor_type;
313 
314 private:
318  Number value;
319 
323  template <typename OtherNumber>
324  void unroll_recursion(Vector<OtherNumber> &result,
325  unsigned int &start_index) const;
326 
330  template <int, int, typename> friend class Tensor;
331 };
332 
333 
334 
376 template <int rank_, int dim, typename Number>
377 class Tensor
378 {
379 public:
388  static const unsigned int dimension = dim;
389 
393  static const unsigned int rank = rank_;
394 
399  static const unsigned int
401 
407  typedef typename Tensor<rank_-1,dim,Number>::tensor_type value_type;
408 
413  typedef typename Tensor<rank_-1,dim,Number>::array_type
414  array_type[(dim != 0) ? dim : 1];
415  // ... avoid a compiler warning in case of dim == 0 and ensure that the
416  // array always has positive size.
417 
423  DEAL_II_CUDA_HOST_DEV Tensor ();
424 
428  explicit Tensor (const array_type &initializer);
429 
435  template <typename OtherNumber>
436  Tensor (const Tensor<rank_,dim,OtherNumber> &initializer);
437 
441  template <typename OtherNumber>
442  Tensor (const Tensor<1,dim,Tensor<rank_-1,dim,OtherNumber> > &initializer);
443 
447  template <typename OtherNumber>
448  operator Tensor<1,dim,Tensor<rank_-1,dim,OtherNumber> > () const;
449 
455  DEAL_II_CUDA_HOST_DEV value_type &operator [] (const unsigned int i);
456 
462  DEAL_II_CUDA_HOST_DEV const value_type &operator[](const unsigned int i) const;
463 
467  const Number &operator [] (const TableIndices<rank_> &indices) const;
468 
472  Number &operator [] (const TableIndices<rank_> &indices);
473 
477  Number *
478  begin_raw();
479 
483  const Number *
484  begin_raw() const;
485 
489  Number *
490  end_raw();
491 
495  const Number *
496  end_raw() const;
497 
503  template <typename OtherNumber>
505 
512  Tensor &operator = (const Number &d);
513 
517  template <typename OtherNumber>
518  bool operator == (const Tensor<rank_,dim,OtherNumber> &) const;
519 
523  template <typename OtherNumber>
524  bool operator != (const Tensor<rank_,dim,OtherNumber> &) const;
525 
529  template <typename OtherNumber>
531 
535  template <typename OtherNumber>
537 
544  template <typename OtherNumber>
545  DEAL_II_CUDA_HOST_DEV Tensor &operator *= (const OtherNumber &factor);
546 
550  template <typename OtherNumber>
551  Tensor &operator /= (const OtherNumber &factor);
552 
556  Tensor operator - () const;
557 
570  void clear ();
571 
579 
586  DEAL_II_CUDA_HOST_DEV typename numbers::NumberTraits<Number>::real_type norm_square() const;
587 
595  template <typename OtherNumber>
596  void unroll (Vector<OtherNumber> &result) const;
597 
602  static
603  unsigned int
605 
610  static
612 
617  static std::size_t memory_consumption ();
618 
623  template <class Archive>
624  void serialize(Archive &ar, const unsigned int version);
625 
631 
632 private:
636  Tensor<rank_-1, dim, Number> values[(dim != 0) ? dim : 1];
637  // ... avoid a compiler warning in case of dim == 0 and ensure that the
638  // array always has positive size.
639 
643  template <typename OtherNumber>
644  void unroll_recursion(Vector<OtherNumber> &result,
645  unsigned int &start_index) const;
646 
650  template <int, int, typename> friend class Tensor;
651 
656  friend class Point<dim,Number>;
657 };
658 
659 
660 namespace internal
661 {
672  template <int rank, int dim, typename T>
673  struct NumberType<Tensor<rank,dim,T> >
674  {
675  static const Tensor<rank,dim,T> &value (const Tensor<rank,dim,T> &t)
676  {
677  return t;
678  }
679 
680  static Tensor<rank,dim,T> value (const T &t)
681  {
682  Tensor<rank,dim,T> tmp;
683  tmp=t;
684  return tmp;
685  }
686  };
687 
688  template <int rank, int dim, typename T>
689  struct NumberType<Tensor<rank,dim,VectorizedArray<T> > >
690  {
691  static const Tensor<rank,dim,VectorizedArray<T> > &value (const Tensor<rank,dim,VectorizedArray<T> > &t)
692  {
693  return t;
694  }
695 
696  static Tensor<rank,dim,VectorizedArray<T> > value (const T &t)
697  {
700  return tmp;
701  }
702 
704  {
706  tmp=t;
707  return tmp;
708  }
709  };
710 }
711 
712 
713 /*---------------------- Inline functions: Tensor<0,dim> ---------------------*/
714 
715 
716 template <int dim,typename Number>
717 inline
718 DEAL_II_CUDA_HOST_DEV Tensor<0,dim,Number>::Tensor ()
719 // Some auto-differentiable numbers need explicit
720 // zero initialization.
721  : value(internal::NumberType<Number>::value(0.0))
722 {
723 }
724 
725 
726 
727 template <int dim, typename Number>
728 template <typename OtherNumber>
729 inline
730 Tensor<0,dim,Number>::Tensor (const OtherNumber &initializer)
731 {
732  value = internal::NumberType<Number>::value(initializer);
733 }
734 
735 
736 
737 template <int dim, typename Number>
738 template <typename OtherNumber>
739 inline
741 {
742  value = p.value;
743 }
744 
745 
746 
747 template <int dim, typename Number>
748 inline
749 Number *
751 {
752  return std::addressof(value);
753 }
754 
755 
756 
757 template <int dim, typename Number>
758 inline
759 const Number *
761 {
762  return std::addressof(value);
763 }
764 
765 
766 
767 template <int dim, typename Number>
768 inline
769 Number *
771 {
773 }
774 
775 
776 
777 template <int dim, typename Number>
778 inline
779 const Number *
781 {
783 }
784 
785 
786 
787 template <int dim, typename Number>
788 inline
789 DEAL_II_CUDA_HOST_DEV Tensor<0,dim,Number>::operator Number &()
790 {
791  // We cannot use Assert inside a CUDA kernel
792 #ifndef __CUDA_ARCH__
793  Assert(dim != 0, ExcMessage("Cannot access an object of type Tensor<0,0,Number>"));
794 #endif
795  return value;
796 }
797 
798 
799 template <int dim, typename Number>
800 inline
801 DEAL_II_CUDA_HOST_DEV Tensor<0,dim,Number>::operator const Number &() const
802 {
803  // We cannot use Assert inside a CUDA kernel
804 #ifndef __CUDA_ARCH__
805  Assert(dim != 0, ExcMessage("Cannot access an object of type Tensor<0,0,Number>"));
806 #endif
807  return value;
808 }
809 
810 
811 template <int dim, typename Number>
812 template <typename OtherNumber>
813 inline
815 {
817  return *this;
818 }
819 
820 
821 #ifdef __INTEL_COMPILER
822 template <int dim, typename Number>
823 inline
825 {
826  value = p.value;
827  return *this;
828 }
829 #endif
830 
831 
832 template <int dim, typename Number>
833 template <typename OtherNumber>
834 inline
836 {
838  return *this;
839 }
840 
841 
842 template <int dim, typename Number>
843 template <typename OtherNumber>
844 inline
846 {
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."));
851 #endif
852 
853  return numbers::values_are_equal(value,p.value);
854 }
855 
856 
857 template <int dim, typename Number>
858 template <typename OtherNumber>
859 inline
861 {
862  return !((*this) == p);
863 }
864 
865 
866 template <int dim, typename Number>
867 template <typename OtherNumber>
868 inline
870 {
871  value += p.value;
872  return *this;
873 }
874 
875 
876 template <int dim, typename Number>
877 template <typename OtherNumber>
878 inline
880 {
881  value -= p.value;
882  return *this;
883 }
884 
885 
886 template <int dim, typename Number>
887 template <typename OtherNumber>
888 inline
889 DEAL_II_CUDA_HOST_DEV Tensor<0,dim,Number> &Tensor<0,dim,Number>::operator *= (const OtherNumber &s)
890 {
891  value *= s;
892  return *this;
893 }
894 
895 
896 template <int dim, typename Number>
897 template <typename OtherNumber>
898 inline
900 {
901  value /= s;
902  return *this;
903 }
904 
905 
906 template <int dim, typename Number>
907 inline
909 {
910  return -value;
911 }
912 
913 
914 template <int dim, typename Number>
915 inline
918 {
919  Assert(dim != 0, ExcMessage("Cannot access an object of type Tensor<0,0,Number>"));
920  return numbers::NumberTraits<Number>::abs (value);
921 }
922 
923 
924 template <int dim, typename Number>
925 inline
927 DEAL_II_CUDA_HOST_DEV Tensor<0,dim,Number>::norm_square () const
928 {
929  // We cannot use Assert inside a CUDA kernel
930 #ifndef __CUDA_ARCH__
931  Assert(dim != 0, ExcMessage("Cannot access an object of type Tensor<0,0,Number>"));
932 #endif
934 }
935 
936 
937 template <int dim, typename Number>
938 template <typename OtherNumber>
939 inline
940 void
941 Tensor<0, dim, Number>::unroll_recursion (Vector<OtherNumber> &result,
942  unsigned int &index) const
943 {
944  Assert(dim != 0, ExcMessage("Cannot unroll an object of type Tensor<0,0,Number>"));
945  result[index] = value;
946  ++index;
947 }
948 
949 
950 template <int dim, typename Number>
951 inline
953 {
954  // Some auto-differentiable numbers need explicit
955  // zero initialization.
957 }
958 
959 
960 template <int dim, typename Number>
961 template <class Archive>
962 inline
963 void Tensor<0,dim,Number>::serialize(Archive &ar, const unsigned int)
964 {
965  ar &value;
966 }
967 
968 
969 /*-------------------- Inline functions: Tensor<rank,dim> --------------------*/
970 
971 
972 template <int rank_, int dim, typename Number>
973 inline DEAL_II_ALWAYS_INLINE
974 DEAL_II_CUDA_HOST_DEV Tensor<rank_,dim,Number>::Tensor ()
975 {
976  // All members of the c-style array values are already default initialized
977  // and thus all values are already set to zero recursively.
978 }
979 
980 
981 template <int rank_, int dim, typename Number>
982 inline DEAL_II_ALWAYS_INLINE
983 Tensor<rank_,dim,Number>::Tensor (const array_type &initializer)
984 {
985  for (unsigned int i=0; i<dim; ++i)
986  values[i] = Tensor<rank_-1, dim, Number>(initializer[i]);
987 }
988 
989 
990 template <int rank_, int dim, typename Number>
991 template <typename OtherNumber>
992 inline DEAL_II_ALWAYS_INLINE
994 {
995  for (unsigned int i=0; i!=dim; ++i)
996  values[i] = Tensor<rank_-1,dim,Number>(initializer[i]);
997 }
998 
999 
1000 template <int rank_, int dim, typename Number>
1001 template <typename OtherNumber>
1002 inline DEAL_II_ALWAYS_INLINE
1004 (const Tensor<1,dim,Tensor<rank_-1,dim,OtherNumber> > &initializer)
1005 {
1006  for (unsigned int i=0; i<dim; ++i)
1007  values[i] = initializer[i];
1008 }
1009 
1010 
1011 template <int rank_, int dim, typename Number>
1012 template <typename OtherNumber>
1013 inline DEAL_II_ALWAYS_INLINE
1015 operator Tensor<1,dim,Tensor<rank_-1,dim,OtherNumber> > () const
1016 {
1017  return Tensor<1,dim,Tensor<rank_-1,dim,Number> > (values);
1018 }
1019 
1020 
1021 
1022 namespace internal
1023 {
1024  namespace TensorSubscriptor
1025  {
1026  template <typename ArrayElementType, int dim>
1027  inline DEAL_II_ALWAYS_INLINE
1028  DEAL_II_CUDA_HOST_DEV
1029  ArrayElementType &
1030  subscript (ArrayElementType *values,
1031  const unsigned int i,
1032  std::integral_constant<int, dim>)
1033  {
1034  // We cannot use Assert in a CUDA kernel
1035 #ifndef __CUDA_ARCH__
1036  Assert (i<dim, ExcIndexRange(i, 0, dim));
1037 #endif
1038  return values[i];
1039  }
1040 
1041 
1042  template <typename ArrayElementType>
1043  DEAL_II_CUDA_HOST_DEV
1044  ArrayElementType &
1045  subscript (ArrayElementType *,
1046  const unsigned int,
1047  std::integral_constant<int, 0>)
1048  {
1049  // We cannot use Assert in a CUDA kernel
1050 #ifndef __CUDA_ARCH__
1051  Assert(false, ExcMessage("Cannot access elements of an object of type Tensor<rank,0,Number>."));
1052 #endif
1053  static ArrayElementType t;
1054  return t;
1055  }
1056  }
1057 }
1058 
1059 
1060 template <int rank_, int dim, typename Number>
1061 inline DEAL_II_ALWAYS_INLINE
1062 DEAL_II_CUDA_HOST_DEV
1065 {
1066  return ::internal::TensorSubscriptor::subscript(values, i, std::integral_constant<int, dim>());
1067 }
1068 
1069 
1070 template <int rank_, int dim, typename Number>
1071 inline DEAL_II_ALWAYS_INLINE
1072 DEAL_II_CUDA_HOST_DEV
1073 const typename Tensor<rank_,dim,Number>::value_type &
1074 Tensor<rank_,dim,Number>::operator[] (const unsigned int i) const
1075 {
1076  return ::internal::TensorSubscriptor::subscript(values, i, std::integral_constant<int, dim>());
1077 }
1078 
1079 
1080 template <int rank_, int dim, typename Number>
1081 inline
1082 const Number &
1084 {
1085  Assert(dim != 0, ExcMessage("Cannot access an object of type Tensor<rank_,0,Number>"));
1086 
1087  return TensorAccessors::extract<rank_>(*this, indices);
1088 }
1089 
1090 
1091 
1092 template <int rank_, int dim, typename Number>
1093 inline
1094 Number &
1096 {
1097  Assert(dim != 0, ExcMessage("Cannot access an object of type Tensor<rank_,0,Number>"));
1098 
1099  return TensorAccessors::extract<rank_>(*this, indices);
1100 }
1101 
1102 
1103 
1104 template <int rank_, int dim, typename Number>
1105 inline
1106 Number *
1108 {
1109  return std::addressof(this->operator[](this->unrolled_to_component_indices(0)));
1110 }
1111 
1112 
1113 
1114 template <int rank_, int dim, typename Number>
1115 inline
1116 const Number *
1118 {
1119  return std::addressof(this->operator[](this->unrolled_to_component_indices(0)));
1120 }
1121 
1122 
1123 
1124 template <int rank_, int dim, typename Number>
1125 inline
1126 Number *
1128 {
1129  return begin_raw()+n_independent_components;
1130 }
1131 
1132 
1133 
1134 template <int rank_, int dim, typename Number>
1135 inline
1136 const Number *
1138 {
1139  return begin_raw()+n_independent_components;
1140 }
1141 
1142 
1143 
1144 template <int rank_, int dim, typename Number>
1145 template <typename OtherNumber>
1146 inline DEAL_II_ALWAYS_INLINE
1149 {
1150  if (dim > 0)
1151  std::copy (&t.values[0], &t.values[0]+dim, &values[0]);
1152  return *this;
1153 }
1154 
1155 
1156 template <int rank_, int dim, typename Number>
1157 inline DEAL_II_ALWAYS_INLINE
1160 {
1161  Assert (numbers::value_is_zero(d), ExcMessage ("Only assignment with zero is allowed"));
1162  (void) d;
1163 
1164  for (unsigned int i=0; i<dim; ++i)
1165  values[i] = internal::NumberType<Number>::value(0.0);
1166  return *this;
1167 }
1168 
1169 
1170 template <int rank_, int dim, typename Number>
1171 template <typename OtherNumber>
1172 inline
1173 bool
1175 {
1176  for (unsigned int i=0; i<dim; ++i)
1177  if (values[i] != p.values[i])
1178  return false;
1179  return true;
1180 }
1181 
1182 
1183 // At some places in the library, we have Point<0> for formal reasons
1184 // (e.g., we sometimes have Quadrature<dim-1> for faces, so we have
1185 // Quadrature<0> for dim=1, and then we have Point<0>). To avoid warnings
1186 // in the above function that the loop end check always fails, we
1187 // implement this function here
1188 template <>
1189 template <>
1190 inline
1192 {
1193  return true;
1194 }
1195 
1196 
1197 template <int rank_, int dim, typename Number>
1198 template <typename OtherNumber>
1199 inline
1200 bool
1202 {
1203  return !((*this) == p);
1204 }
1205 
1206 
1207 template <int rank_, int dim, typename Number>
1208 template <typename OtherNumber>
1209 inline
1212 {
1213  for (unsigned int i=0; i<dim; ++i)
1214  values[i] += p.values[i];
1215  return *this;
1216 }
1217 
1218 
1219 template <int rank_, int dim, typename Number>
1220 template <typename OtherNumber>
1221 inline
1224 {
1225  for (unsigned int i=0; i<dim; ++i)
1226  values[i] -= p.values[i];
1227  return *this;
1228 }
1229 
1230 
1231 template <int rank_, int dim, typename Number>
1232 template <typename OtherNumber>
1233 inline
1234 DEAL_II_CUDA_HOST_DEV
1236 Tensor<rank_,dim,Number>::operator *= (const OtherNumber &s)
1237 {
1238  for (unsigned int i=0; i<dim; ++i)
1239  values[i] *= s;
1240  return *this;
1241 }
1242 
1243 
1244 template <int rank_, int dim, typename Number>
1245 template <typename OtherNumber>
1246 inline
1248 Tensor<rank_,dim,Number>::operator /= (const OtherNumber &s)
1249 {
1250  for (unsigned int i=0; i<dim; ++i)
1251  values[i] /= s;
1252  return *this;
1253 }
1254 
1255 
1256 template <int rank_, int dim, typename Number>
1257 inline
1260 {
1262 
1263  for (unsigned int i=0; i<dim; ++i)
1264  tmp.values[i] = -values[i];
1265 
1266  return tmp;
1267 }
1268 
1269 
1270 template <int rank_, int dim, typename Number>
1271 inline
1274 {
1275  return std::sqrt (norm_square());
1276 }
1277 
1278 
1279 template <int rank_, int dim, typename Number>
1280 inline
1281 DEAL_II_CUDA_HOST_DEV
1284 {
1287  for (unsigned int i=0; i<dim; ++i)
1288  s += values[i].norm_square();
1289 
1290  return s;
1291 }
1292 
1293 
1294 template <int rank_, int dim, typename Number>
1295 template <typename OtherNumber>
1296 inline
1297 void
1298 Tensor<rank_, dim, Number>::unroll (Vector<OtherNumber> &result) const
1299 {
1300  AssertDimension (result.size(),(Utilities::fixed_power<rank_, unsigned int>(dim)));
1301 
1302  unsigned int index = 0;
1303  unroll_recursion (result, index);
1304 }
1305 
1306 
1307 template <int rank_, int dim, typename Number>
1308 template <typename OtherNumber>
1309 inline
1310 void
1312  unsigned int &index) const
1313 {
1314  for (unsigned int i=0; i<dim; ++i)
1315  values[i].unroll_recursion(result, index);
1316 }
1317 
1318 
1319 template <int rank_, int dim, typename Number>
1320 inline
1321 unsigned int
1323 {
1324  unsigned int index = 0;
1325  for (int r = 0; r < rank_; ++r)
1326  index = index * dim + indices[r];
1327 
1328  return index;
1329 }
1330 
1331 
1332 template <int rank_, int dim, typename Number>
1333 inline
1336 {
1337  Assert (i < n_independent_components,
1338  ExcIndexRange (i, 0, n_independent_components));
1339 
1340  TableIndices<rank_> indices;
1341 
1342  unsigned int remainder = i;
1343  for (int r=rank_-1; r>=0; --r)
1344  {
1345  indices[r] = (remainder % dim);
1346  remainder /= dim;
1347  }
1348  Assert (remainder == 0, ExcInternalError());
1349 
1350  return indices;
1351 }
1352 
1353 
1354 template <int rank_, int dim, typename Number>
1355 inline
1357 {
1358  for (unsigned int i=0; i<dim; ++i)
1359  values[i] = internal::NumberType<Number>::value(0.0);
1360 }
1361 
1362 
1363 template <int rank_, int dim, typename Number>
1364 inline
1365 std::size_t
1367 {
1368  return sizeof(Tensor<rank_,dim,Number>);
1369 }
1370 
1371 
1372 template <int rank_, int dim, typename Number>
1373 template <class Archive>
1374 inline
1375 void
1376 Tensor<rank_,dim,Number>::serialize(Archive &ar, const unsigned int)
1377 {
1378  ar &values;
1379 }
1380 
1381 
1382 /* ----------------- Non-member functions operating on tensors. ------------ */
1383 
1388 
1396 template <int rank_, int dim, typename Number>
1397 inline
1398 std::ostream &operator << (std::ostream &out, const Tensor<rank_,dim,Number> &p)
1399 {
1400  for (unsigned int i = 0; i < dim; ++i)
1401  {
1402  out << p[i];
1403  if (i != dim - 1)
1404  out << ' ';
1405  }
1406 
1407  return out;
1408 }
1409 
1410 
1417 template <int dim, typename Number>
1418 inline
1419 std::ostream &operator << (std::ostream &out, const Tensor<0,dim,Number> &p)
1420 {
1421  out << static_cast<const Number &>(p);
1422  return out;
1423 }
1424 
1425 
1427 
1431 
1432 
1441 template <int dim, typename Number, typename Other>
1442 inline DEAL_II_ALWAYS_INLINE
1443 typename ProductType<Other, Number>::type
1444 operator * (const Other &object,
1445  const Tensor<0,dim,Number> &t)
1446 {
1447  return object * static_cast<const Number &>(t);
1448 }
1449 
1450 
1451 
1460 template <int dim, typename Number, typename Other>
1461 inline DEAL_II_ALWAYS_INLINE
1462 typename ProductType<Number, Other>::type
1464  const Other &object)
1465 {
1466  return static_cast<const Number &>(t) * object;
1467 }
1468 
1469 
1479 template <int dim, typename Number, typename OtherNumber>
1480 inline DEAL_II_ALWAYS_INLINE
1481 typename ProductType<Number, OtherNumber>::type
1483  const Tensor<0, dim, OtherNumber> &src2)
1484 {
1485  return static_cast<const Number &>(src1) *
1486  static_cast<const OtherNumber &>(src2);
1487 }
1488 
1489 
1495 template <int dim, typename Number, typename OtherNumber>
1496 inline DEAL_II_ALWAYS_INLINE
1499  const OtherNumber &factor)
1500 {
1501  return static_cast<const Number &>(t) / factor;
1502 }
1503 
1504 
1510 template <int dim, typename Number, typename OtherNumber>
1511 inline DEAL_II_ALWAYS_INLINE
1514  const Tensor<0,dim,OtherNumber> &q)
1515 {
1516  return static_cast<const Number &>(p) + static_cast<const OtherNumber &>(q);
1517 }
1518 
1519 
1525 template <int dim, typename Number, typename OtherNumber>
1526 inline DEAL_II_ALWAYS_INLINE
1529  const Tensor<0,dim,OtherNumber> &q)
1530 {
1531  return static_cast<const Number &>(p) - static_cast<const OtherNumber &>(q);
1532 }
1533 
1534 
1545 template <int rank, int dim,
1546  typename Number,
1547  typename OtherNumber>
1548 inline DEAL_II_ALWAYS_INLINE
1551  const OtherNumber &factor)
1552 {
1553  // recurse over the base objects
1555  for (unsigned int d=0; d<dim; ++d)
1556  tt[d] = t[d] * factor;
1557  return tt;
1558 }
1559 
1560 
1571 template <int rank, int dim,
1572  typename Number,
1573  typename OtherNumber>
1574 inline DEAL_II_ALWAYS_INLINE
1576 operator * (const Number &factor,
1578 {
1579  // simply forward to the operator above
1580  return t * factor;
1581 }
1582 
1583 
1591 template <int rank, int dim,
1592  typename Number,
1593  typename OtherNumber>
1594 inline
1597  const OtherNumber &factor)
1598 {
1599  // recurse over the base objects
1601  for (unsigned int d=0; d<dim; ++d)
1602  tt[d] = t[d] / factor;
1603  return tt;
1604 }
1605 
1606 
1614 template <int rank, int dim, typename Number, typename OtherNumber>
1615 inline DEAL_II_ALWAYS_INLINE
1619 {
1621 
1622  for (unsigned int i=0; i<dim; ++i)
1623  tmp[i] += q[i];
1624 
1625  return tmp;
1626 }
1627 
1628 
1636 template <int rank, int dim, typename Number, typename OtherNumber>
1637 inline DEAL_II_ALWAYS_INLINE
1641 {
1643 
1644  for (unsigned int i=0; i<dim; ++i)
1645  tmp[i] -= q[i];
1646 
1647  return tmp;
1648 }
1649 
1650 
1652 
1656 
1657 
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
1687 {
1688  typename Tensor<rank_1 + rank_2 - 2, dim, typename ProductType<Number, OtherNumber>::type>::tensor_type result;
1689 
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);
1693 
1694  return result;
1695 }
1696 
1697 
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
1734 {
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)"));
1739 
1740  using namespace TensorAccessors;
1741  using namespace TensorAccessors::internal;
1742 
1743  // Reorder index_1 to the end of src1:
1744  ReorderedIndexView<index_1, rank_1, const Tensor<rank_1, dim, Number> >
1745  reord_01 = reordered_index_view<index_1, rank_1>(src1);
1746 
1747  // Reorder index_2 to the end of src2:
1748  ReorderedIndexView<index_2, rank_2, const Tensor<rank_2, dim, OtherNumber> >
1749  reord_02 = reordered_index_view<index_2, rank_2>(src2);
1750 
1751  typename Tensor<rank_1 + rank_2 - 2, dim, typename ProductType<Number, OtherNumber>::type>::tensor_type
1752  result;
1753  TensorAccessors::contract<1, rank_1, rank_2, dim>(result, reord_01, reord_02);
1754  return result;
1755 }
1756 
1757 
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>
1792 inline
1793 typename Tensor<rank_1 + rank_2 - 4, dim, typename ProductType<Number, OtherNumber>::type>::tensor_type
1796 {
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"));
1809 
1810  using namespace TensorAccessors;
1811  using namespace TensorAccessors::internal;
1812 
1813  // Reorder index_1 to the end of src1:
1814  ReorderedIndexView<index_1, rank_1, const Tensor<rank_1, dim, Number> >
1815  reord_1 = TensorAccessors::reordered_index_view<index_1, rank_1>(src1);
1816 
1817  // Reorder index_2 to the end of src2:
1818  ReorderedIndexView<index_2, rank_2, const Tensor<rank_2, dim, OtherNumber> >
1819  reord_2 = TensorAccessors::reordered_index_view<index_2, rank_2>(src2);
1820 
1821  // Now, reorder index_3 to the end of src1. We have to make sure to
1822  // preserve the orginial ordering: index_1 has been removed. If
1823  // index_3 > index_1, we have to use (index_3 - 1) instead:
1824  ReorderedIndexView<(index_3 < index_1 ? index_3 : index_3 - 1), rank_1, ReorderedIndexView<index_1, rank_1, const Tensor<rank_1, dim, Number> > >
1825  reord_3 = TensorAccessors::reordered_index_view<index_3 < index_1 ? index_3 : index_3 - 1, rank_1>(reord_1);
1826 
1827  // Now, reorder index_4 to the end of src2. We have to make sure to
1828  // preserve the orginial ordering: index_2 has been removed. If
1829  // index_4 > index_2, we have to use (index_4 - 1) instead:
1830  ReorderedIndexView<(index_4 < index_2 ? index_4 : index_4 - 1), rank_2, ReorderedIndexView<index_2, rank_2, const Tensor<rank_2, dim, OtherNumber> > >
1831  reord_4 = TensorAccessors::reordered_index_view<index_4 < index_2 ? index_4 : index_4 - 1, rank_2>(reord_2);
1832 
1833  typename Tensor<rank_1 + rank_2 - 4, dim, typename ProductType<Number, OtherNumber>::type>::tensor_type
1834  result;
1835  TensorAccessors::contract<2, rank_1, rank_2, dim>(result, reord_3, reord_4);
1836  return result;
1837 }
1838 
1839 
1853 template <int rank, int dim, typename Number, typename OtherNumber>
1854 inline DEAL_II_ALWAYS_INLINE
1855 typename ProductType<Number, OtherNumber>::type
1857  const Tensor<rank, dim, OtherNumber> &right)
1858 {
1859  typename ProductType<Number, OtherNumber>::type result;
1860  TensorAccessors::contract<rank, rank, rank, dim>(result, left, right);
1861  return result;
1862 }
1863 
1864 
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>
1889 contract3 (const TensorT1<rank_1, dim, T1> &left,
1890  const TensorT2<rank_1 + rank_2, dim, T2> &middle,
1891  const TensorT3<rank_2, dim, T3> &right)
1892 {
1894  return_type;
1895  return TensorAccessors::contract3<rank_1, rank_2, dim, return_type>(
1896  left, middle, right);
1897 }
1898 
1899 
1911 template <int rank_1, int rank_2, int dim,
1912  typename Number, typename OtherNumber>
1913 inline DEAL_II_ALWAYS_INLINE
1917 {
1919  TensorAccessors::contract<0, rank_1, rank_2, dim>(result, src1, src2);
1920  return result;
1921 }
1922 
1923 
1925 
1929 
1930 
1942 template <int dim, typename Number>
1943 inline DEAL_II_ALWAYS_INLINE
1946 {
1947  Assert (dim==2, ExcInternalError());
1948 
1949  Tensor<1, dim, Number> result;
1950 
1951  result[0] = src[1];
1952  result[1] = -src[0];
1953 
1954  return result;
1955 }
1956 
1957 
1968 template <int dim, typename Number>
1969 inline DEAL_II_ALWAYS_INLINE
1972  const Tensor<1,dim,Number> &src2)
1973 {
1974  Assert (dim==3, ExcInternalError());
1975 
1976  Tensor<1, dim, Number> result;
1977 
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];
1981 
1982  return result;
1983 }
1984 
1985 
1987 
1991 
1992 
1999 template <int dim, typename Number>
2000 inline
2002 {
2003  // Compute the determinant using the Laplace expansion of the
2004  // determinant. We expand along the last row.
2005  Number det = internal::NumberType<Number>::value(0.0);
2006 
2007  for (unsigned int k=0; k<dim; ++k)
2008  {
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];
2013 
2014  const Number cofactor = ((k % 2 == 0) ? -1. : 1.) * determinant(minor);
2015 
2016  det += t[dim-1][k] * cofactor;
2017  }
2018 
2019  return ((dim % 2 == 0) ? 1. : -1.) * det;
2020 }
2021 
2027 template <typename Number>
2028 inline
2030 {
2031  return t[0][0];
2032 }
2033 
2034 
2042 template <int dim, typename Number>
2043 inline DEAL_II_ALWAYS_INLINE
2044 Number trace (const Tensor<2,dim,Number> &d)
2045 {
2046  Number t=d[0][0];
2047  for (unsigned int i=1; i<dim; ++i)
2048  t += d[i][i];
2049  return t;
2050 }
2051 
2052 
2062 template <int dim, typename Number>
2063 inline
2066 {
2067  Number return_tensor [dim][dim];
2068 
2069  // if desired, take over the
2070  // inversion of a 4x4 tensor
2071  // from the FullMatrix
2072  AssertThrow (false, ExcNotImplemented());
2073 
2074  return Tensor<2,dim,Number>(return_tensor);
2075 }
2076 
2077 
2078 #ifndef DOXYGEN
2079 
2080 template <typename Number>
2081 inline
2083 invert (const Tensor<2,1,Number> &t)
2084 {
2085  Number return_tensor [1][1];
2086 
2087  return_tensor[0][0] = internal::NumberType<Number>::value(1.0/t[0][0]);
2088 
2089  return Tensor<2,1,Number>(return_tensor);
2090 }
2091 
2092 
2093 template <typename Number>
2094 inline
2096 invert (const Tensor<2,2,Number> &t)
2097 {
2098  Tensor<2,2,Number> return_tensor;
2099 
2100  // this is Maple output,
2101  // thus a bit unstructured
2102  const Number inv_det_t = internal::NumberType<Number>::value(1.0/(t[0][0]*t[1][1]-t[1][0]*t[0][1]));
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;
2108 
2109  return return_tensor;
2110 }
2111 
2112 
2113 template <typename Number>
2114 inline
2116 invert (const Tensor<2,3,Number> &t)
2117 {
2118  Tensor<2,3,Number> return_tensor;
2119 
2120  const Number t4 = internal::NumberType<Number>::value(t[0][0]*t[1][1]),
2121  t6 = internal::NumberType<Number>::value(t[0][0]*t[1][2]),
2122  t8 = internal::NumberType<Number>::value(t[0][1]*t[1][0]),
2123  t00 = internal::NumberType<Number>::value(t[0][2]*t[1][0]),
2124  t01 = internal::NumberType<Number>::value(t[0][1]*t[2][0]),
2125  t04 = internal::NumberType<Number>::value(t[0][2]*t[2][0]),
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]));
2129  return_tensor[0][0] = internal::NumberType<Number>::value(t[1][1]*t[2][2])-internal::NumberType<Number>::value(t[1][2]*t[2][1]);
2130  return_tensor[0][1] = internal::NumberType<Number>::value(t[0][2]*t[2][1])-internal::NumberType<Number>::value(t[0][1]*t[2][2]);
2131  return_tensor[0][2] = internal::NumberType<Number>::value(t[0][1]*t[1][2])-internal::NumberType<Number>::value(t[0][2]*t[1][1]);
2132  return_tensor[1][0] = internal::NumberType<Number>::value(t[1][2]*t[2][0])-internal::NumberType<Number>::value(t[1][0]*t[2][2]);
2133  return_tensor[1][1] = internal::NumberType<Number>::value(t[0][0]*t[2][2])-t04;
2134  return_tensor[1][2] = t00-t6;
2135  return_tensor[2][0] = internal::NumberType<Number>::value(t[1][0]*t[2][1])-internal::NumberType<Number>::value(t[1][1]*t[2][0]);
2136  return_tensor[2][1] = t01-internal::NumberType<Number>::value(t[0][0]*t[2][1]);
2137  return_tensor[2][2] = internal::NumberType<Number>::value(t4-t8);
2138  return_tensor *= inv_det_t;
2139 
2140  return return_tensor;
2141 }
2142 
2143 #endif /* DOXYGEN */
2144 
2145 
2152 template <int dim, typename Number>
2153 inline DEAL_II_ALWAYS_INLINE
2156 {
2158  for (unsigned int i=0; i<dim; ++i)
2159  {
2160  tt[i][i] = t[i][i];
2161  for (unsigned int j=i+1; j<dim; ++j)
2162  {
2163  tt[i][j] = t[j][i];
2164  tt[j][i] = t[i][j];
2165  };
2166  }
2167  return tt;
2168 }
2169 
2170 
2184 template <int dim, typename Number>
2185 inline
2188 {
2189  return determinant(t)*invert(t);
2190 }
2191 
2192 
2207 template <int dim, typename Number>
2208 inline
2211 {
2212  return transpose(adjugate(t));
2213 }
2214 
2215 
2223 template <int dim, typename Number>
2224 inline
2225 Number
2227 {
2228  Number max = internal::NumberType<Number>::value(0.0);
2229  for (unsigned int j=0; j<dim; ++j)
2230  {
2231  Number sum = internal::NumberType<Number>::value(0.0);
2232  for (unsigned int i=0; i<dim; ++i)
2233  sum += std::fabs(t[i][j]);
2234 
2235  if (sum > max)
2236  max = sum;
2237  }
2238 
2239  return max;
2240 }
2241 
2242 
2250 template <int dim, typename Number>
2251 inline
2252 Number
2254 {
2255  Number max = internal::NumberType<Number>::value(0.0);
2256  for (unsigned int i=0; i<dim; ++i)
2257  {
2258  Number sum = internal::NumberType<Number>::value(0.0);
2259  for (unsigned int j=0; j<dim; ++j)
2260  sum += std::fabs(t[i][j]);
2261 
2262  if (sum > max)
2263  max = sum;
2264  }
2265 
2266  return max;
2267 }
2268 
2270 
2271 
2272 #ifndef DOXYGEN
2273 
2274 
2275 #ifdef DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING
2276 
2277 // Specialization of functions for Adol-C number types when
2278 // the advanced branching feature is used
2279 template <int dim>
2280 inline
2281 adouble
2282 l1_norm (const Tensor<2,dim,adouble> &t)
2283 {
2284  adouble max = internal::NumberType<adouble>::value(0.0);
2285  for (unsigned int j=0; j<dim; ++j)
2286  {
2287  adouble sum = internal::NumberType<adouble>::value(0.0);
2288  for (unsigned int i=0; i<dim; ++i)
2289  sum += std::fabs(t[i][j]);
2290 
2291  condassign(max, (sum > max), sum, max);
2292  }
2293 
2294  return max;
2295 }
2296 
2297 
2298 template <int dim>
2299 inline
2300 adouble
2301 linfty_norm (const Tensor<2,dim,adouble> &t)
2302 {
2304  for (unsigned int i=0; i<dim; ++i)
2305  {
2307  for (unsigned int j=0; j<dim; ++j)
2308  sum += std::fabs(t[i][j]);
2309 
2310  condassign(max, (sum > max), sum, max);
2311  }
2312 
2313  return max;
2314 }
2315 
2316 #endif // DEAL_II_ADOLC_WITH_ADVANCED_BRANCHING
2317 
2318 
2319 #endif // DOXYGEN
2320 
2321 DEAL_II_NAMESPACE_CLOSE
2322 
2323 // include deprecated non-member functions operating on Tensor
2324 #include <deal.II/base/tensor_deprecated.h>
2325 
2326 #endif
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
Definition: tensor.h:122
Number trace(const Tensor< 2, dim, Number > &d)
Definition: tensor.h:2044
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
Number determinant(const SymmetricTensor< 2, dim, Number > &)
static unsigned int component_to_unrolled_index(const TableIndices< rank_ > &indices)
Definition: tensor.h:1322
static const unsigned int n_independent_components
Definition: tensor.h:400
bool value_is_zero(const Number &value)
Definition: numbers.h:784
bool operator==(const Tensor< rank_, dim, OtherNumber > &) const
Definition: tensor.h:1174
Number l1_norm(const Tensor< 2, dim, Number > &t)
Definition: tensor.h:2226
static std::size_t memory_consumption()
Definition: tensor.h:1366
Tensor & operator/=(const OtherNumber &factor)
bool operator!=(const Tensor< rank_, dim, OtherNumber > &) const
Definition: tensor.h:1201
Number linfty_norm(const Tensor< 2, dim, Number > &t)
Definition: tensor.h:2253
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)
Definition: tensor.h:1794
Tensor< 2, dim, Number > adjugate(const Tensor< 2, dim, Number > &t)
Definition: tensor.h:2187
static const unsigned int rank
Definition: tensor.h:393
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static real_type abs(const number &x)
Definition: numbers.h:465
internal::ReorderedIndexView< index, rank, T > reordered_index_view(T &t)
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:1273
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
Definition: tensor.h:1311
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Tensor & operator=(const Tensor< rank_, dim, OtherNumber > &rhs)
Definition: point.h:104
bool values_are_equal(const Number1 &value_1, const Number2 &value_2)
Definition: numbers.h:768
Tensor & operator*=(const OtherNumber &factor)
Tensor< 2, dim, Number > cofactor(const Tensor< 2, dim, Number > &t)
Definition: tensor.h:2210
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)
Definition: tensor.h:2155
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)
Definition: tensor.h:1889
T sum(const T &t, const MPI_Comm &mpi_communicator)
Tensor< 2, dim, Number > invert(const Tensor< 2, dim, Number > &)
Definition: tensor.h:2065
static real_type abs_square(const number &x)
Definition: numbers.h:456
Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition: tensor.h:1945
#define Assert(cond, exc)
Definition: exceptions.h:1142
ProductType< Number, OtherNumber >::type scalar_product(const Tensor< rank, dim, Number > &left, const Tensor< rank, dim, OtherNumber > &right)
Definition: tensor.h:1856
Tensor< rank_, dim, Number > tensor_type
Definition: tensor.h:630
void serialize(Archive &ar, const unsigned int version)
Definition: tensor.h:1376
numbers::NumberTraits< Number >::real_type norm_square() const
Definition: tensor.h:1283
Number * end_raw()
Definition: tensor.h:1127
value_type & operator[](const unsigned int i)
Definition: tensor.h:1064
Tensor()
Definition: tensor.h:974
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]
Definition: tensor.h:414
static TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
Definition: tensor.h:1335
Tensor< rank_-1, dim, Number > values[(dim !=0) ? dim :1]
Definition: tensor.h:636
Tensor< 1, dim, Number > cross_product_3d(const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2)
Definition: tensor.h:1971
Number determinant(const Tensor< 2, dim, Number > &t)
Definition: tensor.h:2001
static const unsigned int dimension
Definition: tensor.h:388
Definition: mpi.h:53
Tensor< rank_-1, dim, Number >::tensor_type value_type
Definition: tensor.h:407
Number * begin_raw()
Definition: tensor.h:1107
static ::ExceptionBase & ExcNotImplemented()
void unroll(Vector< OtherNumber > &result) const
Definition: tensor.h:1298
Tensor operator-() const
Definition: tensor.h:1259
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)
Definition: tensor.h:1732
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)
Definition: tensor.h:1915
Tensor & operator+=(const Tensor< rank_, dim, OtherNumber > &)
void clear()
Definition: tensor.h:1356
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)
Definition: tensor.h:2029