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\}}\)
full_matrix.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2019 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_full_matrix_h
17 #define dealii_full_matrix_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/numbers.h>
23 #include <deal.II/base/table.h>
24 #include <deal.II/base/tensor.h>
25 
27 
28 #include <deal.II/lac/exceptions.h>
30 
31 #include <cstring>
32 #include <iomanip>
33 #include <vector>
34 
36 
37 
38 // forward declarations
39 #ifndef DOXYGEN
40 template <typename number>
41 class Vector;
42 template <typename number>
43 class LAPACKFullMatrix;
44 #endif
45 
70 template <typename number>
71 class FullMatrix : public Table<2, number>
72 {
73 public:
74  // The assertion in full_matrix.templates.h for whether or not a number is
75  // finite is not compatible for AD number types.
76  static_assert(
78  "The FullMatrix class does not support auto-differentiable numbers.");
79 
83  using size_type = std::size_t;
84 
89  using value_type = number;
90 
95 
100 
105 
109  using Table<2, number>::end;
110 
121 
126 
136  explicit FullMatrix(const size_type n = 0);
137 
141  FullMatrix(const size_type rows, const size_type cols);
142 
147  FullMatrix(const size_type rows, const size_type cols, const number *entries);
148 
157  FullMatrix(const IdentityMatrix &id);
172  template <typename number2>
175 
184  operator=(const number d);
185 
195  operator=(const IdentityMatrix &id);
196 
201  template <typename number2>
204 
205 
211  template <typename MatrixType>
212  void
213  copy_from(const MatrixType &);
214 
220  template <typename MatrixType>
221  void
222  copy_transposed(const MatrixType &);
223 
231  template <int dim>
232  void
233  copy_from(const Tensor<2, dim> &T,
234  const unsigned int src_r_i = 0,
235  const unsigned int src_r_j = dim - 1,
236  const unsigned int src_c_i = 0,
237  const unsigned int src_c_j = dim - 1,
238  const size_type dst_r = 0,
239  const size_type dst_c = 0);
240 
248  template <int dim>
249  void copy_to(Tensor<2, dim> & T,
250  const size_type src_r_i = 0,
251  const size_type src_r_j = dim - 1,
252  const size_type src_c_i = 0,
253  const size_type src_c_j = dim - 1,
254  const unsigned int dst_r = 0,
255  const unsigned int dst_c = 0) const;
256 
269  template <typename MatrixType, typename index_type>
270  void
271  extract_submatrix_from(const MatrixType & matrix,
272  const std::vector<index_type> &row_index_set,
273  const std::vector<index_type> &column_index_set);
274 
287  template <typename MatrixType, typename index_type>
288  void
289  scatter_matrix_to(const std::vector<index_type> &row_index_set,
290  const std::vector<index_type> &column_index_set,
291  MatrixType & matrix) const;
292 
303  template <typename number2>
304  void
305  fill(const FullMatrix<number2> &src,
306  const size_type dst_offset_i = 0,
307  const size_type dst_offset_j = 0,
308  const size_type src_offset_i = 0,
309  const size_type src_offset_j = 0);
310 
311 
315  template <typename number2>
316  void
317  fill(const number2 *);
318 
330  template <typename number2>
331  void
333  const std::vector<size_type> &p_rows,
334  const std::vector<size_type> &p_cols);
335 
346  void
347  set(const size_type i, const size_type j, const number value);
363  bool
364  operator==(const FullMatrix<number> &) const;
365 
370  size_type
371  m() const;
372 
377  size_type
378  n() const;
379 
385  bool
386  all_zero() const;
387 
403  template <typename number2>
404  number2
405  matrix_norm_square(const Vector<number2> &v) const;
406 
416  template <typename number2>
417  number2
418  matrix_scalar_product(const Vector<number2> &u,
419  const Vector<number2> &v) const;
420 
425  real_type
426  l1_norm() const;
427 
432  real_type
433  linfty_norm() const;
434 
442  real_type
443  frobenius_norm() const;
444 
453  real_type
454  relative_symmetry_norm2() const;
455 
461  number
462  determinant() const;
463 
469  number
470  trace() const;
471 
478  template <class StreamType>
479  void
480  print(StreamType & s,
481  const unsigned int width = 5,
482  const unsigned int precision = 2) const;
483 
506  void
507  print_formatted(std::ostream & out,
508  const unsigned int precision = 3,
509  const bool scientific = true,
510  const unsigned int width = 0,
511  const char * zero_string = " ",
512  const double denominator = 1.,
513  const double threshold = 0.) const;
514 
519  std::size_t
520  memory_consumption() const;
521 
523 
525 
529  iterator
530  begin(const size_type r);
531 
535  iterator
536  end(const size_type r);
537 
542  begin(const size_type r) const;
543 
548  end(const size_type r) const;
549 
551 
553 
557  FullMatrix &
558  operator*=(const number factor);
559 
563  FullMatrix &
564  operator/=(const number factor);
565 
573  template <typename number2>
574  void
575  add(const number a, const FullMatrix<number2> &A);
576 
584  template <typename number2>
585  void
586  add(const number a,
587  const FullMatrix<number2> &A,
588  const number b,
589  const FullMatrix<number2> &B);
590 
599  template <typename number2>
600  void
601  add(const number a,
602  const FullMatrix<number2> &A,
603  const number b,
604  const FullMatrix<number2> &B,
605  const number c,
606  const FullMatrix<number2> &C);
607 
619  template <typename number2>
620  void
621  add(const FullMatrix<number2> &src,
622  const number factor,
623  const size_type dst_offset_i = 0,
624  const size_type dst_offset_j = 0,
625  const size_type src_offset_i = 0,
626  const size_type src_offset_j = 0);
627 
633  template <typename number2>
634  void
635  Tadd(const number s, const FullMatrix<number2> &B);
636 
648  template <typename number2>
649  void
650  Tadd(const FullMatrix<number2> &src,
651  const number factor,
652  const size_type dst_offset_i = 0,
653  const size_type dst_offset_j = 0,
654  const size_type src_offset_i = 0,
655  const size_type src_offset_j = 0);
656 
660  void
661  add(const size_type row, const size_type column, const number value);
662 
672  template <typename number2, typename index_type>
673  void
674  add(const size_type row,
675  const size_type n_cols,
676  const index_type *col_indices,
677  const number2 * values,
678  const bool elide_zero_values = true,
679  const bool col_indices_are_sorted = false);
680 
684  void
685  add_row(const size_type i, const number s, const size_type j);
686 
691  void
692  add_row(const size_type i,
693  const number s,
694  const size_type j,
695  const number t,
696  const size_type k);
697 
701  void
702  add_col(const size_type i, const number s, const size_type j);
703 
708  void
709  add_col(const size_type i,
710  const number s,
711  const size_type j,
712  const number t,
713  const size_type k);
714 
718  void
719  swap_row(const size_type i, const size_type j);
720 
724  void
725  swap_col(const size_type i, const size_type j);
726 
731  void
732  diagadd(const number s);
733 
737  template <typename number2>
738  void
739  equ(const number a, const FullMatrix<number2> &A);
740 
744  template <typename number2>
745  void
746  equ(const number a,
747  const FullMatrix<number2> &A,
748  const number b,
749  const FullMatrix<number2> &B);
750 
754  template <typename number2>
755  void
756  equ(const number a,
757  const FullMatrix<number2> &A,
758  const number b,
759  const FullMatrix<number2> &B,
760  const number c,
761  const FullMatrix<number2> &C);
762 
769  void
770  symmetrize();
771 
786  void
787  gauss_jordan();
788 
795  template <typename number2>
796  void
797  invert(const FullMatrix<number2> &M);
798 
807  template <typename number2>
808  void
810 
815  template <typename number2>
816  void
817  outer_product(const Vector<number2> &V, const Vector<number2> &W);
818 
824  template <typename number2>
825  void
827 
833  template <typename number2>
834  void
836 
838 
840 
859  template <typename number2>
860  void
862  const FullMatrix<number2> &B,
863  const bool adding = false) const;
864 
883  template <typename number2>
884  void
886  const FullMatrix<number2> &B,
887  const bool adding = false) const;
888 
907  template <typename number2>
908  void
910  const FullMatrix<number2> &B,
911  const bool adding = false) const;
912 
932  template <typename number2>
933  void
935  const FullMatrix<number2> &B,
936  const bool adding = false) const;
937 
948  void
950  const FullMatrix<number> &B,
951  const FullMatrix<number> &D,
952  const bool transpose_B = false,
953  const bool transpose_D = false,
954  const number scaling = number(1.));
955 
968  template <typename number2>
969  void
970  vmult(Vector<number2> & w,
971  const Vector<number2> &v,
972  const bool adding = false) const;
973 
979  template <typename number2>
980  void
981  vmult_add(Vector<number2> &w, const Vector<number2> &v) const;
982 
996  template <typename number2>
997  void
998  Tvmult(Vector<number2> & w,
999  const Vector<number2> &v,
1000  const bool adding = false) const;
1001 
1008  template <typename number2>
1009  void
1010  Tvmult_add(Vector<number2> &w, const Vector<number2> &v) const;
1011 
1017  template <typename somenumber>
1018  void
1019  precondition_Jacobi(Vector<somenumber> & dst,
1020  const Vector<somenumber> &src,
1021  const number omega = 1.) const;
1022 
1029  template <typename number2, typename number3>
1030  number
1031  residual(Vector<number2> & dst,
1032  const Vector<number2> &x,
1033  const Vector<number3> &b) const;
1034 
1045  template <typename number2>
1046  void
1047  forward(Vector<number2> &dst, const Vector<number2> &src) const;
1048 
1056  template <typename number2>
1057  void
1058  backward(Vector<number2> &dst, const Vector<number2> &src) const;
1059 
1061 
1071 
1076  ExcNotRegular,
1077  number,
1078  << "The maximal pivot is " << arg1
1079  << ", which is below the threshold. The matrix may be singular.");
1084  size_type,
1085  size_type,
1086  size_type,
1087  << "Target region not in matrix: size in this direction="
1088  << arg1 << ", size of new matrix=" << arg2
1089  << ", offset=" << arg3);
1094  "You are attempting an operation on two matrices that "
1095  "are the same object, but the operation requires that the "
1096  "two objects are in fact different.");
1102 };
1103 
1106 #ifndef DOXYGEN
1107 /*-------------------------Inline functions -------------------------------*/
1108 
1109 
1110 
1111 template <typename number>
1112 inline typename FullMatrix<number>::size_type
1113 FullMatrix<number>::m() const
1114 {
1115  return this->n_rows();
1116 }
1117 
1118 
1119 
1120 template <typename number>
1121 inline typename FullMatrix<number>::size_type
1122 FullMatrix<number>::n() const
1123 {
1124  return this->n_cols();
1125 }
1126 
1127 
1128 
1129 template <typename number>
1131 FullMatrix<number>::operator=(const number d)
1132 {
1133  Assert(d == number(0), ExcScalarAssignmentOnlyForZeroValue());
1134  (void)d; // removes -Wunused-parameter warning in optimized mode
1135 
1136  if (this->n_elements() != 0)
1137  this->reset_values();
1138 
1139  return *this;
1140 }
1141 
1142 
1143 
1144 template <typename number>
1145 template <typename number2>
1146 inline void
1147 FullMatrix<number>::fill(const number2 *src)
1148 {
1150 }
1151 
1152 
1153 
1154 template <typename number>
1155 template <typename MatrixType>
1156 void
1157 FullMatrix<number>::copy_from(const MatrixType &M)
1158 {
1159  this->reinit(M.m(), M.n());
1160 
1161  // loop over the elements of the argument matrix row by row, as suggested
1162  // in the documentation of the sparse matrix iterator class, and
1163  // copy them into the current object
1164  for (size_type row = 0; row < M.m(); ++row)
1165  {
1166  const typename MatrixType::const_iterator end_row = M.end(row);
1167  for (typename MatrixType::const_iterator entry = M.begin(row);
1168  entry != end_row;
1169  ++entry)
1170  this->el(row, entry->column()) = entry->value();
1171  }
1172 }
1173 
1174 
1175 
1176 template <typename number>
1177 template <typename MatrixType>
1178 void
1179 FullMatrix<number>::copy_transposed(const MatrixType &M)
1180 {
1181  this->reinit(M.n(), M.m());
1182 
1183  // loop over the elements of the argument matrix row by row, as suggested
1184  // in the documentation of the sparse matrix iterator class, and
1185  // copy them into the current object
1186  for (size_type row = 0; row < M.m(); ++row)
1187  {
1188  const typename MatrixType::const_iterator end_row = M.end(row);
1189  for (typename MatrixType::const_iterator entry = M.begin(row);
1190  entry != end_row;
1191  ++entry)
1192  this->el(entry->column(), row) = entry->value();
1193  }
1194 }
1195 
1196 
1197 
1198 template <typename number>
1199 template <typename MatrixType, typename index_type>
1200 inline void
1202  const MatrixType & matrix,
1203  const std::vector<index_type> &row_index_set,
1204  const std::vector<index_type> &column_index_set)
1205 {
1206  AssertDimension(row_index_set.size(), this->n_rows());
1207  AssertDimension(column_index_set.size(), this->n_cols());
1208 
1209  const size_type n_rows_submatrix = row_index_set.size();
1210  const size_type n_cols_submatrix = column_index_set.size();
1211 
1212  for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1213  for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1214  (*this)(sub_row, sub_col) =
1215  matrix.el(row_index_set[sub_row], column_index_set[sub_col]);
1216 }
1217 
1218 
1219 
1220 template <typename number>
1221 template <typename MatrixType, typename index_type>
1222 inline void
1224  const std::vector<index_type> &row_index_set,
1225  const std::vector<index_type> &column_index_set,
1226  MatrixType & matrix) const
1227 {
1228  AssertDimension(row_index_set.size(), this->n_rows());
1229  AssertDimension(column_index_set.size(), this->n_cols());
1230 
1231  const size_type n_rows_submatrix = row_index_set.size();
1232  const size_type n_cols_submatrix = column_index_set.size();
1233 
1234  for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1235  for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1236  matrix.set(row_index_set[sub_row],
1237  column_index_set[sub_col],
1238  (*this)(sub_row, sub_col));
1239 }
1240 
1241 
1242 template <typename number>
1243 inline void
1245  const size_type j,
1246  const number value)
1247 {
1248  (*this)(i, j) = value;
1249 }
1250 
1251 
1252 
1253 template <typename number>
1254 template <typename number2>
1255 void
1256 FullMatrix<number>::vmult_add(Vector<number2> & w,
1257  const Vector<number2> &v) const
1258 {
1259  vmult(w, v, true);
1260 }
1261 
1262 
1263 template <typename number>
1264 template <typename number2>
1265 void
1266 FullMatrix<number>::Tvmult_add(Vector<number2> & w,
1267  const Vector<number2> &v) const
1268 {
1269  Tvmult(w, v, true);
1270 }
1271 
1272 
1273 //---------------------------------------------------------------------------
1274 template <typename number>
1275 inline typename FullMatrix<number>::iterator
1277 {
1278  AssertIndexRange(r, m());
1279  return iterator(this, r, 0);
1280 }
1281 
1282 
1283 
1284 template <typename number>
1285 inline typename FullMatrix<number>::iterator
1287 {
1288  AssertIndexRange(r, m());
1289  return iterator(this, r + 1, 0);
1290 }
1291 
1292 
1293 
1294 template <typename number>
1295 inline typename FullMatrix<number>::const_iterator
1296 FullMatrix<number>::begin(const size_type r) const
1297 {
1298  AssertIndexRange(r, m());
1299  return const_iterator(this, r, 0);
1300 }
1301 
1302 
1303 
1304 template <typename number>
1305 inline typename FullMatrix<number>::const_iterator
1306 FullMatrix<number>::end(const size_type r) const
1307 {
1308  AssertIndexRange(r, m());
1309  return const_iterator(this, r + 1, 0);
1310 }
1311 
1312 
1313 
1314 template <typename number>
1315 inline void
1316 FullMatrix<number>::add(const size_type r, const size_type c, const number v)
1317 {
1318  AssertIndexRange(r, this->m());
1319  AssertIndexRange(c, this->n());
1320 
1321  this->operator()(r, c) += v;
1322 }
1323 
1324 
1325 
1326 template <typename number>
1327 template <typename number2, typename index_type>
1328 inline void
1330  const size_type n_cols,
1331  const index_type *col_indices,
1332  const number2 * values,
1333  const bool,
1334  const bool)
1335 {
1336  AssertIndexRange(row, this->m());
1337  for (size_type col = 0; col < n_cols; ++col)
1338  {
1339  AssertIndexRange(col_indices[col], this->n());
1340  this->operator()(row, col_indices[col]) += values[col];
1341  }
1342 }
1343 
1344 
1345 template <typename number>
1346 template <class StreamType>
1347 inline void
1348 FullMatrix<number>::print(StreamType & s,
1349  const unsigned int w,
1350  const unsigned int p) const
1351 {
1352  Assert(!this->empty(), ExcEmptyMatrix());
1353 
1354  // save the state of out stream
1355  const std::streamsize old_precision = s.precision(p);
1356  const std::streamsize old_width = s.width(w);
1357 
1358  for (size_type i = 0; i < this->m(); ++i)
1359  {
1360  for (size_type j = 0; j < this->n(); ++j)
1361  {
1362  s.width(w);
1363  s.precision(p);
1364  s << this->el(i, j);
1365  }
1366  s << std::endl;
1367  }
1368 
1369  // reset output format
1370  s.precision(old_precision);
1371  s.width(old_width);
1372 }
1373 
1374 
1375 #endif // DOXYGEN
1376 
1378 
1379 #endif
FullMatrix::end
iterator end(const size_type r)
FullMatrix::Tvmult
void Tvmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
FullMatrix::trace
number trace() const
DeclExceptionMsg
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:496
LinearAlgebraDealII::Vector
Vector< double > Vector
Definition: generic_linear_algebra.h:43
FullMatrix::TmTmult
void TmTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
FullMatrix::equ
void equ(const number a, const FullMatrix< number2 > &A)
FullMatrix::operator=
FullMatrix< number > & operator=(const FullMatrix< number2 > &)
FullMatrix::vmult
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
FullMatrix::copy_from
void copy_from(const MatrixType &)
FullMatrix::ExcNotRegular
static ::ExceptionBase & ExcNotRegular(number arg1)
FullMatrix::ExcEmptyMatrix
static ::ExceptionBase & ExcEmptyMatrix()
FullMatrix::copy_transposed
void copy_transposed(const MatrixType &)
FullMatrix::copy_to
void copy_to(Tensor< 2, dim > &T, const size_type src_r_i=0, const size_type src_r_j=dim - 1, const size_type src_c_i=0, const size_type src_c_j=dim - 1, const unsigned int dst_r=0, const unsigned int dst_c=0) const
FullMatrix::forward
void forward(Vector< number2 > &dst, const Vector< number2 > &src) const
FullMatrix::add
void add(const number a, const FullMatrix< number2 > &A)
FullMatrix::operator*=
FullMatrix & operator*=(const number factor)
FullMatrix::Tmmult
void Tmmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
FullMatrix::set
void set(const size_type i, const size_type j, const number value)
FullMatrix::determinant
number determinant() const
FullMatrix::m
size_type m() const
FullMatrix::add_col
void add_col(const size_type i, const number s, const size_type j)
LAPACKSupport::V
static const char V
Definition: lapack_support.h:175
FullMatrix::fill
void fill(const FullMatrix< number2 > &src, const size_type dst_offset_i=0, const size_type dst_offset_j=0, const size_type src_offset_i=0, const size_type src_offset_j=0)
FullMatrix::ExcSourceEqualsDestination
static ::ExceptionBase & ExcSourceEqualsDestination()
FullMatrix::relative_symmetry_norm2
real_type relative_symmetry_norm2() const
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Physics::Elasticity::Kinematics::C
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
FullMatrix::all_zero
bool all_zero() const
FullMatrix::n
size_type n() const
exceptions.h
FullMatrix::linfty_norm
real_type linfty_norm() const
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
FullMatrix::print
void print(StreamType &s, const unsigned int width=5, const unsigned int precision=2) const
LAPACKFullMatrix
Definition: lapack_full_matrix.h:60
Table
Definition: table.h:699
FullMatrix::print_formatted
void print_formatted(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1., const double threshold=0.) const
FullMatrix::size_type
std::size_t size_type
Definition: full_matrix.h:83
FullMatrix::symmetrize
void symmetrize()
LAPACKSupport::T
static const char T
Definition: lapack_support.h:163
Physics::Elasticity::Kinematics::w
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
tensor.h
FullMatrix::begin
iterator begin(const size_type r)
FullMatrix::residual
number residual(Vector< number2 > &dst, const Vector< number2 > &x, const Vector< number3 > &b) const
FullMatrix::precondition_Jacobi
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
FullMatrix::mTmult
void mTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
TableBase< N, CoefficientType >::size_type
typename AlignedVector< CoefficientType >::size_type size_type
Definition: table.h:420
internal::reinit
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:621
Tensor< 2, dim >
IdentityMatrix
Definition: identity_matrix.h:71
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
FullMatrix< CoefficientType >::iterator
typename Table< 2, CoefficientType >::iterator iterator
Definition: full_matrix.h:94
FullMatrix::swap_col
void swap_col(const size_type i, const size_type j)
DeclException3
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:564
FullMatrix::matrix_scalar_product
number2 matrix_scalar_product(const Vector< number2 > &u, const Vector< number2 > &v) const
TableBase::fill
void fill(InputIterator entries, const bool C_style_indexing=true)
FullMatrix::outer_product
void outer_product(const Vector< number2 > &V, const Vector< number2 > &W)
FullMatrix::FullMatrix
FullMatrix(const size_type n=0)
DeclException1
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:518
numbers::NumberTraits::real_type
number real_type
Definition: numbers.h:437
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
FullMatrix::diagadd
void diagadd(const number s)
FullMatrix::cholesky
void cholesky(const FullMatrix< number2 > &A)
FullMatrix::swap_row
void swap_row(const size_type i, const size_type j)
FullMatrix::triple_product
void triple_product(const FullMatrix< number > &A, const FullMatrix< number > &B, const FullMatrix< number > &D, const bool transpose_B=false, const bool transpose_D=false, const number scaling=number(1.))
FullMatrix< CoefficientType >::real_type
typename numbers::NumberTraits< CoefficientType >::real_type real_type
Definition: full_matrix.h:120
identity_matrix.h
FullMatrix::ExcMatrixNotPositiveDefinite
static ::ExceptionBase & ExcMatrixNotPositiveDefinite()
FullMatrix::frobenius_norm
real_type frobenius_norm() const
LAPACKSupport::A
static const char A
Definition: lapack_support.h:155
FullMatrix::ExcInvalidDestination
static ::ExceptionBase & ExcInvalidDestination(size_type arg1, size_type arg2, size_type arg3)
value
static const bool value
Definition: dof_tools_constraints.cc:433
FullMatrix::gauss_jordan
void gauss_jordan()
LinearAlgebra::CUDAWrappers::kernel::size_type
types::global_dof_index size_type
Definition: cuda_kernels.h:45
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
ad_number_traits.h
FullMatrix::Tadd
void Tadd(const number s, const FullMatrix< number2 > &B)
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
FullMatrix::backward
void backward(Vector< number2 > &dst, const Vector< number2 > &src) const
FullMatrix::scatter_matrix_to
void scatter_matrix_to(const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set, MatrixType &matrix) const
FullMatrix::mmult
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
FullMatrix::invert
void invert(const FullMatrix< number2 > &M)
FullMatrix::left_invert
void left_invert(const FullMatrix< number2 > &M)
FullMatrix< CoefficientType >::const_iterator
typename Table< 2, CoefficientType >::const_iterator const_iterator
Definition: full_matrix.h:99
config.h
StandardExceptions::ExcScalarAssignmentOnlyForZeroValue
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
FullMatrix
Definition: full_matrix.h:71
FullMatrix::memory_consumption
std::size_t memory_consumption() const
FullMatrix::right_invert
void right_invert(const FullMatrix< number2 > &M)
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
FullMatrix::Tvmult_add
void Tvmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
FullMatrix::operator==
bool operator==(const FullMatrix< number > &) const
TableBase< N, number >::values
AlignedVector< number > values
Definition: table.h:672
TableBase< N, CoefficientType >::value_type
CoefficientType value_type
Definition: table.h:415
FullMatrix::add_row
void add_row(const size_type i, const number s, const size_type j)
Differentiation::AD::is_ad_number
Definition: numbers.h:658
FullMatrix::vmult_add
void vmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
FullMatrix::l1_norm
real_type l1_norm() const
numbers.h
table.h
FullMatrix::operator/=
FullMatrix & operator/=(const number factor)
FullMatrix::fill_permutation
void fill_permutation(const FullMatrix< number2 > &src, const std::vector< size_type > &p_rows, const std::vector< size_type > &p_cols)
FullMatrix::extract_submatrix_from
void extract_submatrix_from(const MatrixType &matrix, const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set)
FullMatrix::matrix_norm_square
number2 matrix_norm_square(const Vector< number2 > &v) const