Reference documentation for deal.II version 9.0.0
full_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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_full_matrix_h
17 #define dealii_full_matrix_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/numbers.h>
22 #include <deal.II/base/table.h>
23 #include <deal.II/lac/exceptions.h>
24 #include <deal.II/lac/identity_matrix.h>
25 #include <deal.II/base/tensor.h>
26 
27 #include <vector>
28 #include <iomanip>
29 #include <cstring>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
34 // forward declarations
35 template <typename number> class Vector;
36 template <typename number> class LAPACKFullMatrix;
37 
38 
63 template <typename number>
64 class FullMatrix : public Table<2,number>
65 {
66 public:
70  typedef std::size_t size_type;
71 
76  typedef number value_type;
77 
78 
89 
90 
91  class const_iterator;
92 
96  class Accessor
97  {
98  public:
104  const size_type row,
105  const size_type col);
106 
110  size_type row() const;
111 
115  size_type column() const;
116 
120  number value() const;
121 
122  protected:
127 
132 
137 
138  /*
139  * Make enclosing class a friend.
140  */
141  friend class const_iterator;
142  };
143 
148  {
149  public:
153  const_iterator(const FullMatrix<number> *matrix,
154  const size_type row,
155  const size_type col);
156 
161 
166 
170  const Accessor &operator* () const;
171 
175  const Accessor *operator-> () const;
176 
180  bool operator == (const const_iterator &) const;
184  bool operator != (const const_iterator &) const;
185 
190  bool operator < (const const_iterator &) const;
191 
196  bool operator > (const const_iterator &) const;
197 
198  private:
203  };
208 
218  explicit FullMatrix (const size_type n = 0);
219 
223  FullMatrix (const size_type rows,
224  const size_type cols);
225 
230  FullMatrix (const size_type rows,
231  const size_type cols,
232  const number *entries);
233 
242  FullMatrix (const IdentityMatrix &id);
257  template <typename number2>
260 
269  operator = (const number d);
270 
280  operator = (const IdentityMatrix &id);
281 
286  template <typename number2>
289 
290 
296  template <typename MatrixType>
297  void copy_from (const MatrixType &);
298 
304  template <typename MatrixType>
305  void copy_transposed (const MatrixType &);
306 
314  template <int dim>
315  void
316  copy_from (const Tensor<2,dim> &T,
317  const unsigned int src_r_i=0,
318  const unsigned int src_r_j=dim-1,
319  const unsigned int src_c_i=0,
320  const unsigned int src_c_j=dim-1,
321  const size_type dst_r=0,
322  const size_type dst_c=0);
323 
331  template <int dim>
332  void
334  const size_type src_r_i=0,
335  const size_type src_r_j=dim-1,
336  const size_type src_c_i=0,
337  const size_type src_c_j=dim-1,
338  const unsigned int dst_r=0,
339  const unsigned int dst_c=0) const;
340 
353  template <typename MatrixType, typename index_type>
354  void extract_submatrix_from (const MatrixType &matrix,
355  const std::vector<index_type> &row_index_set,
356  const std::vector<index_type> &column_index_set);
357 
370  template <typename MatrixType, typename index_type>
371  void
372  scatter_matrix_to (const std::vector<index_type> &row_index_set,
373  const std::vector<index_type> &column_index_set,
374  MatrixType &matrix) const;
375 
386  template <typename number2>
387  void fill (const FullMatrix<number2> &src,
388  const size_type dst_offset_i = 0,
389  const size_type dst_offset_j = 0,
390  const size_type src_offset_i = 0,
391  const size_type src_offset_j = 0);
392 
393 
397  template <typename number2>
398  void fill (const number2 *);
399 
411  template <typename number2>
412  void fill_permutation (const FullMatrix<number2> &src,
413  const std::vector<size_type> &p_rows,
414  const std::vector<size_type> &p_cols);
415 
426  void set (const size_type i,
427  const size_type j,
428  const number value);
444  bool operator == (const FullMatrix<number> &) const;
445 
450  size_type m () const;
451 
456  size_type n () const;
457 
463  bool all_zero () const;
464 
480  template <typename number2>
481  number2 matrix_norm_square (const Vector<number2> &v) const;
482 
492  template <typename number2>
493  number2 matrix_scalar_product (const Vector<number2> &u,
494  const Vector<number2> &v) const;
495 
500  real_type l1_norm () const;
501 
506  real_type linfty_norm () const;
507 
515  real_type frobenius_norm () const;
516 
526 
532  number determinant () const;
533 
539  number trace () const;
540 
547  template <class StreamType>
548  void print (StreamType &s,
549  const unsigned int width=5,
550  const unsigned int precision=2) const;
551 
574  void print_formatted (std::ostream &out,
575  const unsigned int precision=3,
576  const bool scientific = true,
577  const unsigned int width = 0,
578  const char *zero_string = " ",
579  const double denominator = 1.,
580  const double threshold = 0.) const;
581 
586  std::size_t memory_consumption () const;
587 
589 
591 
595  const_iterator begin () const;
596 
600  const_iterator end () const;
601 
605  const_iterator begin (const size_type r) const;
606 
610  const_iterator end (const size_type r) const;
611 
613 
615 
619  FullMatrix &operator *= (const number factor);
620 
624  FullMatrix &operator /= (const number factor);
625 
633  template <typename number2>
634  void add (const number a,
635  const FullMatrix<number2> &A);
636 
644  template <typename number2>
645  void add (const number a,
646  const FullMatrix<number2> &A,
647  const number b,
648  const FullMatrix<number2> &B);
649 
658  template <typename number2>
659  void add (const number a,
660  const FullMatrix<number2> &A,
661  const number b,
662  const FullMatrix<number2> &B,
663  const number c,
664  const FullMatrix<number2> &C);
665 
677  template <typename number2>
678  void add (const FullMatrix<number2> &src,
679  const number factor,
680  const size_type dst_offset_i = 0,
681  const size_type dst_offset_j = 0,
682  const size_type src_offset_i = 0,
683  const size_type src_offset_j = 0);
684 
690  template <typename number2>
691  void Tadd (const number s,
692  const FullMatrix<number2> &B);
693 
705  template <typename number2>
706  void Tadd (const FullMatrix<number2> &src,
707  const number factor,
708  const size_type dst_offset_i = 0,
709  const size_type dst_offset_j = 0,
710  const size_type src_offset_i = 0,
711  const size_type src_offset_j = 0);
712 
716  void add (const size_type row,
717  const size_type column,
718  const number value);
719 
729  template <typename number2, typename index_type>
730  void add (const size_type row,
731  const size_type n_cols,
732  const index_type *col_indices,
733  const number2 *values,
734  const bool elide_zero_values = true,
735  const bool col_indices_are_sorted = false);
736 
740  void add_row (const size_type i,
741  const number s,
742  const size_type j);
743 
748  void add_row (const size_type i,
749  const number s, const size_type j,
750  const number t, const size_type k);
751 
755  void add_col (const size_type i,
756  const number s,
757  const size_type j);
758 
763  void add_col (const size_type i,
764  const number s, const size_type j,
765  const number t, const size_type k);
766 
770  void swap_row (const size_type i,
771  const size_type j);
772 
776  void swap_col (const size_type i,
777  const size_type j);
778 
783  void diagadd (const number s);
784 
788  template <typename number2>
789  void equ (const number a,
790  const FullMatrix<number2> &A);
791 
795  template <typename number2>
796  void equ (const number a,
797  const FullMatrix<number2> &A,
798  const number b,
799  const FullMatrix<number2> &B);
800 
804  template <typename number2>
805  void equ (const number a,
806  const FullMatrix<number2> &A,
807  const number b,
808  const FullMatrix<number2> &B,
809  const number c,
810  const FullMatrix<number2> &C);
811 
818  void symmetrize ();
819 
834  void gauss_jordan ();
835 
842  template <typename number2>
843  void invert (const FullMatrix<number2> &M);
844 
853  template <typename number2>
854  void cholesky (const FullMatrix<number2> &A);
855 
860  template <typename number2>
861  void outer_product (const Vector<number2> &V,
862  const Vector<number2> &W);
863 
869  template <typename number2>
870  void left_invert (const FullMatrix<number2> &M);
871 
877  template <typename number2>
878  void right_invert (const FullMatrix<number2> &M);
879 
881 
883 
902  template <typename number2>
903  void mmult (FullMatrix<number2> &C,
904  const FullMatrix<number2> &B,
905  const bool adding=false) const;
906 
925  template <typename number2>
926  void Tmmult (FullMatrix<number2> &C,
927  const FullMatrix<number2> &B,
928  const bool adding=false) const;
929 
948  template <typename number2>
949  void mTmult (FullMatrix<number2> &C,
950  const FullMatrix<number2> &B,
951  const bool adding=false) const;
952 
972  template <typename number2>
973  void TmTmult (FullMatrix<number2> &C,
974  const FullMatrix<number2> &B,
975  const bool adding=false) const;
976 
987  void triple_product(const FullMatrix<number> &A,
988  const FullMatrix<number> &B,
989  const FullMatrix<number> &D,
990  const bool transpose_B = false,
991  const bool transpose_D = false,
992  const number scaling = number(1.));
993 
1006  template <typename number2>
1007  void vmult (Vector<number2> &w,
1008  const Vector<number2> &v,
1009  const bool adding=false) const;
1010 
1016  template <typename number2>
1017  void vmult_add (Vector<number2> &w,
1018  const Vector<number2> &v) const;
1019 
1033  template <typename number2>
1034  void Tvmult (Vector<number2> &w,
1035  const Vector<number2> &v,
1036  const bool adding=false) const;
1037 
1044  template <typename number2>
1045  void Tvmult_add (Vector<number2> &w,
1046  const Vector<number2> &v) const;
1047 
1053  template <typename somenumber>
1054  void precondition_Jacobi (Vector<somenumber> &dst,
1055  const Vector<somenumber> &src,
1056  const number omega = 1.) const;
1057 
1064  template <typename number2, typename number3>
1065  number residual (Vector<number2> &dst,
1066  const Vector<number2> &x,
1067  const Vector<number3> &b) const;
1068 
1079  template <typename number2>
1080  void forward (Vector<number2> &dst,
1081  const Vector<number2> &src) const;
1082 
1090  template <typename number2>
1091  void backward (Vector<number2> &dst,
1092  const Vector<number2> &src) const;
1093 
1095 
1105 
1110  number,
1111  << "The maximal pivot is " << arg1
1112  << ", which is below the threshold. The matrix may be singular.");
1118  << "Target region not in matrix: size in this direction="
1119  << arg1 << ", size of new matrix=" << arg2
1120  << ", offset=" << arg3);
1125  "You are attempting an operation on two matrices that "
1126  "are the same object, but the operation requires that the "
1127  "two objects are in fact different.");
1133 
1134  friend class Accessor;
1135 };
1136 
1139 #ifndef DOXYGEN
1140 /*-------------------------Inline functions -------------------------------*/
1141 
1142 
1143 
1144 
1145 template <typename number>
1146 inline
1148 FullMatrix<number>::m() const
1149 {
1150  return this->n_rows();
1151 }
1152 
1153 
1154 
1155 template <typename number>
1156 inline
1158 FullMatrix<number>::n() const
1159 {
1160  return this->n_cols();
1161 }
1162 
1163 
1164 
1165 template <typename number>
1167 FullMatrix<number>::operator = (const number d)
1168 {
1169  Assert (d==number(0), ExcScalarAssignmentOnlyForZeroValue());
1170  (void)d; // removes -Wunused-parameter warning in optimized mode
1171 
1172  if (this->n_elements() != 0)
1173  this->reset_values();
1174 
1175  return *this;
1176 }
1177 
1178 
1179 
1180 template <typename number>
1181 template <typename number2>
1182 inline
1183 void FullMatrix<number>::fill (const number2 *src)
1184 {
1185  Table<2,number>::fill(src);
1186 }
1187 
1188 
1189 
1190 template <typename number>
1191 template <typename MatrixType>
1192 void
1193 FullMatrix<number>::copy_from (const MatrixType &M)
1194 {
1195  this->reinit (M.m(), M.n());
1196 
1197  // loop over the elements of the argument matrix row by row, as suggested
1198  // in the documentation of the sparse matrix iterator class, and
1199  // copy them into the current object
1200  for (size_type row = 0; row < M.m(); ++row)
1201  {
1202  const typename MatrixType::const_iterator end_row = M.end(row);
1203  for (typename MatrixType::const_iterator entry = M.begin(row);
1204  entry != end_row; ++entry)
1205  this->el(row, entry->column()) = entry->value();
1206  }
1207 }
1208 
1209 
1210 
1211 template <typename number>
1212 template <typename MatrixType>
1213 void
1214 FullMatrix<number>::copy_transposed (const MatrixType &M)
1215 {
1216  this->reinit (M.n(), M.m());
1217 
1218  // loop over the elements of the argument matrix row by row, as suggested
1219  // in the documentation of the sparse matrix iterator class, and
1220  // copy them into the current object
1221  for (size_type row = 0; row < M.m(); ++row)
1222  {
1223  const typename MatrixType::const_iterator end_row = M.end(row);
1224  for (typename MatrixType::const_iterator entry = M.begin(row);
1225  entry != end_row; ++entry)
1226  this->el(entry->column(), row) = entry->value();
1227  }
1228 }
1229 
1230 
1231 
1232 template <typename number>
1233 template <typename MatrixType, typename index_type>
1234 inline
1235 void
1236 FullMatrix<number>::extract_submatrix_from (const MatrixType &matrix,
1237  const std::vector<index_type> &row_index_set,
1238  const std::vector<index_type> &column_index_set)
1239 {
1240  AssertDimension(row_index_set.size(), this->n_rows());
1241  AssertDimension(column_index_set.size(), this->n_cols());
1242 
1243  const size_type n_rows_submatrix = row_index_set.size();
1244  const size_type n_cols_submatrix = column_index_set.size();
1245 
1246  for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1247  for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1248  (*this)(sub_row, sub_col) = matrix.el(row_index_set[sub_row], column_index_set[sub_col]);
1249 }
1250 
1251 
1252 
1253 template <typename number>
1254 template <typename MatrixType, typename index_type>
1255 inline
1256 void
1257 FullMatrix<number>::scatter_matrix_to (const std::vector<index_type> &row_index_set,
1258  const std::vector<index_type> &column_index_set,
1259  MatrixType &matrix) const
1260 {
1261  AssertDimension(row_index_set.size(), this->n_rows());
1262  AssertDimension(column_index_set.size(), this->n_cols());
1263 
1264  const size_type n_rows_submatrix = row_index_set.size();
1265  const size_type n_cols_submatrix = column_index_set.size();
1266 
1267  for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1268  for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1269  matrix.set(row_index_set[sub_row],
1270  column_index_set[sub_col],
1271  (*this)(sub_row, sub_col));
1272 }
1273 
1274 
1275 template <typename number>
1276 inline
1277 void
1278 FullMatrix<number>::set (const size_type i,
1279  const size_type j,
1280  const number value)
1281 {
1282  (*this)(i,j) = value;
1283 }
1284 
1285 
1286 
1287 template <typename number>
1288 template <typename number2>
1289 void
1290 FullMatrix<number>::vmult_add (Vector<number2> &w,
1291  const Vector<number2> &v) const
1292 {
1293  vmult(w, v, true);
1294 }
1295 
1296 
1297 template <typename number>
1298 template <typename number2>
1299 void
1300 FullMatrix<number>::Tvmult_add (Vector<number2> &w,
1301  const Vector<number2> &v) const
1302 {
1303  Tvmult(w, v, true);
1304 }
1305 
1306 
1307 //---------------------------------------------------------------------------
1308 
1309 
1310 template <typename number>
1311 inline
1313 Accessor (const FullMatrix<number> *matrix,
1314  const size_type r,
1315  const size_type c)
1316  :
1317  matrix(matrix),
1318  a_row(r),
1319  a_col(c)
1320 {}
1321 
1322 
1323 template <typename number>
1324 inline
1327 {
1328  return a_row;
1329 }
1330 
1331 
1332 template <typename number>
1333 inline
1336 {
1337  return a_col;
1338 }
1339 
1340 
1341 template <typename number>
1342 inline
1343 number
1345 {
1346  AssertIsFinite(matrix->el(a_row, a_col));
1347  return matrix->el(a_row, a_col);
1348 }
1349 
1350 
1351 template <typename number>
1352 inline
1354 const_iterator(const FullMatrix<number> *matrix,
1355  const size_type r,
1356  const size_type c)
1357  :
1358  accessor(matrix, r, c)
1359 {}
1360 
1361 
1362 template <typename number>
1363 inline
1366 {
1367  Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
1368 
1369  ++accessor.a_col;
1370  if (accessor.a_col >= accessor.matrix->n())
1371  {
1372  accessor.a_col = 0;
1373  accessor.a_row++;
1374  }
1375  return *this;
1376 }
1377 
1378 
1379 template <typename number>
1380 inline
1383 {
1384  const typename FullMatrix<number>::const_iterator current = *this;
1385  ++(*this);
1386 
1387  return current;
1388 }
1389 
1390 
1391 template <typename number>
1392 inline
1393 const typename FullMatrix<number>::Accessor &
1395 {
1396  return accessor;
1397 }
1398 
1399 
1400 template <typename number>
1401 inline
1402 const typename FullMatrix<number>::Accessor *
1404 {
1405  return &accessor;
1406 }
1407 
1408 
1409 template <typename number>
1410 inline
1411 bool
1413 operator == (const const_iterator &other) const
1414 {
1415  return (accessor.row() == other.accessor.row() &&
1416  accessor.column() == other.accessor.column());
1417 }
1418 
1419 
1420 template <typename number>
1421 inline
1422 bool
1424 operator != (const const_iterator &other) const
1425 {
1426  return ! (*this == other);
1427 }
1428 
1429 
1430 template <typename number>
1431 inline
1432 bool
1434 operator < (const const_iterator &other) const
1435 {
1436  return (accessor.row() < other.accessor.row() ||
1437  (accessor.row() == other.accessor.row() &&
1438  accessor.column() < other.accessor.column()));
1439 }
1440 
1441 
1442 template <typename number>
1443 inline
1444 bool
1446 operator > (const const_iterator &other) const
1447 {
1448  return (other < *this);
1449 }
1450 
1451 
1452 template <typename number>
1453 inline
1456 {
1457  return const_iterator(this, 0, 0);
1458 }
1459 
1460 
1461 template <typename number>
1462 inline
1464 FullMatrix<number>::end () const
1465 {
1466  return const_iterator(this, m(), 0);
1467 }
1468 
1469 
1470 template <typename number>
1471 inline
1473 FullMatrix<number>::begin (const size_type r) const
1474 {
1475  AssertIndexRange(r,m());
1476  return const_iterator(this, r, 0);
1477 }
1478 
1479 
1480 
1481 template <typename number>
1482 inline
1484 FullMatrix<number>::end (const size_type r) const
1485 {
1486  AssertIndexRange(r,m());
1487  return const_iterator(this, r+1, 0);
1488 }
1489 
1490 
1491 
1492 template <typename number>
1493 inline
1494 void
1495 FullMatrix<number>::add (const size_type r, const size_type c, const number v)
1496 {
1497  AssertIndexRange(r, this->m());
1498  AssertIndexRange(c, this->n());
1499 
1500  this->operator()(r,c) += v;
1501 }
1502 
1503 
1504 
1505 template <typename number>
1506 template <typename number2, typename index_type>
1507 inline
1508 void
1510  const size_type n_cols,
1511  const index_type *col_indices,
1512  const number2 *values,
1513  const bool,
1514  const bool)
1515 {
1516  AssertIndexRange(row, this->m());
1517  for (size_type col=0; col<n_cols; ++col)
1518  {
1519  AssertIndexRange(col_indices[col], this->n());
1520  this->operator()(row,col_indices[col]) += values[col];
1521  }
1522 }
1523 
1524 
1525 template <typename number>
1526 template <class StreamType>
1527 inline
1528 void
1529 FullMatrix<number>::print (StreamType &s,
1530  const unsigned int w,
1531  const unsigned int p) const
1532 {
1533  Assert (!this->empty(), ExcEmptyMatrix());
1534 
1535  // save the state of out stream
1536  const std::streamsize old_precision = s.precision (p);
1537  const std::streamsize old_width = s.width (w);
1538 
1539  for (size_type i=0; i<this->m(); ++i)
1540  {
1541  for (size_type j=0; j<this->n(); ++j)
1542  {
1543  s.width(w);
1544  s.precision(p);
1545  s << this->el(i,j);
1546  }
1547  s << std::endl;
1548  }
1549 
1550  // reset output format
1551  s.precision(old_precision);
1552  s.width(old_width);
1553 }
1554 
1555 
1556 #endif // DOXYGEN
1557 
1558 DEAL_II_NAMESPACE_CLOSE
1559 
1560 #endif
size_type m() const
number determinant() const
const Accessor & operator*() const
void diagadd(const number s)
bool operator==(const FullMatrix< number > &) const
static ::ExceptionBase & ExcEmptyMatrix()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
void vmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
const_iterator(const FullMatrix< number > *matrix, const size_type row, const size_type col)
number2 matrix_scalar_product(const Vector< number2 > &u, const Vector< number2 > &v) const
Contents is actually a matrix.
FullMatrix(const size_type n=0)
bool operator==(const const_iterator &) const
FullMatrix & operator/=(const number factor)
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
void gauss_jordan()
static ::ExceptionBase & ExcNotRegular(number arg1)
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
void cholesky(const FullMatrix< number2 > &A)
void add_row(const size_type i, const number s, const size_type j)
bool operator>(const const_iterator &) const
const_iterator begin() const
void outer_product(const Vector< number2 > &V, const Vector< number2 > &W)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
void right_invert(const FullMatrix< number2 > &M)
std::size_t size_type
Definition: full_matrix.h:70
void scatter_matrix_to(const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set, MatrixType &matrix) const
void fill_permutation(const FullMatrix< number2 > &src, const std::vector< size_type > &p_rows, const std::vector< size_type > &p_cols)
void left_invert(const FullMatrix< number2 > &M)
void Tadd(const number s, const FullMatrix< number2 > &B)
void copy_transposed(const MatrixType &)
real_type l1_norm() const
FullMatrix & operator*=(const number factor)
const_iterator & operator++()
void set(const size_type i, const size_type j, const number value)
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.))
const Accessor * operator->() const
real_type frobenius_norm() const
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
const FullMatrix< number > * matrix
Definition: full_matrix.h:126
void Tmmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void backward(Vector< number2 > &dst, const Vector< number2 > &src) const
size_type row() const
numbers::NumberTraits< number >::real_type real_type
Definition: full_matrix.h:88
Accessor(const FullMatrix< number > *matrix, const size_type row, const size_type col)
size_type n() const
void swap_row(const size_type i, const size_type j)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:346
real_type linfty_norm() const
void invert(const FullMatrix< number2 > &M)
#define Assert(cond, exc)
Definition: exceptions.h:1142
bool operator!=(const const_iterator &) const
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:335
bool operator<(const const_iterator &) const
#define DeclException0(Exception0)
Definition: exceptions.h:323
AlignedVector< T >::reference el(const TableIndices< N > &indices)
void extract_submatrix_from(const MatrixType &matrix, const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set)
void Tvmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
number value() const
const_iterator end() const
real_type relative_symmetry_norm2() const
number residual(Vector< number2 > &dst, const Vector< number2 > &x, const Vector< number3 > &b) const
FullMatrix< number > & operator=(const FullMatrix< number2 > &)
void TmTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
number2 matrix_norm_square(const Vector< number2 > &v) const
void add_col(const size_type i, const number s, const size_type j)
static ::ExceptionBase & ExcIteratorPastEnd()
number trace() const
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
std::size_t memory_consumption() const
void Tvmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void forward(Vector< number2 > &dst, const Vector< number2 > &src) const
number value_type
Definition: full_matrix.h:76
static ::ExceptionBase & ExcSourceEqualsDestination()
void swap_col(const size_type i, const size_type j)
void copy_from(const MatrixType &)
void add(const number a, const FullMatrix< number2 > &A)
static ::ExceptionBase & ExcInvalidDestination(size_type arg1, size_type arg2, size_type arg3)
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
Definition: table.h:34
void print(StreamType &s, const unsigned int width=5, const unsigned int precision=2) const
bool empty() const
void symmetrize()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:370
static ::ExceptionBase & ExcMatrixNotPositiveDefinite()
AlignedVector< T >::size_type size_type
Definition: table.h:413
AlignedVector< T >::reference operator()(const TableIndices< N > &indices)
AlignedVector< T > values
Definition: table.h:648
void equ(const number a, const FullMatrix< number2 > &A)
#define AssertIsFinite(number)
Definition: exceptions.h:1299
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)
void fill(InputIterator entries, const bool C_style_indexing=true)
void mTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
size_type column() const
bool all_zero() const