Reference documentation for deal.II version 9.3.3
\(\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 - 2020 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
23#include <deal.II/base/table.h>
24#include <deal.II/base/tensor.h>
25
27
30
31#include <cstring>
32#include <iomanip>
33#include <vector>
34
36
37
38// forward declarations
39#ifndef DOXYGEN
40template <typename number>
41class Vector;
42template <typename number>
44#endif
45
68template <typename number>
69class FullMatrix : public Table<2, number>
70{
71public:
72 // The assertion in full_matrix.templates.h for whether or not a number is
73 // finite is not compatible for AD number types.
74 static_assert(
76 "The FullMatrix class does not support auto-differentiable numbers.");
77
81 using size_type = std::size_t;
82
87 using value_type = number;
88
93
98
102 using Table<2, number>::begin;
103
107 using Table<2, number>::end;
108
119
124
134 explicit FullMatrix(const size_type n = 0);
135
139 FullMatrix(const size_type rows, const size_type cols);
140
145 FullMatrix(const size_type rows, const size_type cols, const number *entries);
146
170 template <typename number2>
173
182 operator=(const number d);
183
194
199 template <typename number2>
202
203
209 template <typename MatrixType>
210 void
211 copy_from(const MatrixType &);
212
218 template <typename MatrixType>
219 void
220 copy_transposed(const MatrixType &);
221
229 template <int dim>
230 void
232 const unsigned int src_r_i = 0,
233 const unsigned int src_r_j = dim - 1,
234 const unsigned int src_c_i = 0,
235 const unsigned int src_c_j = dim - 1,
236 const size_type dst_r = 0,
237 const size_type dst_c = 0);
238
246 template <int dim>
248 const size_type src_r_i = 0,
249 const size_type src_r_j = dim - 1,
250 const size_type src_c_i = 0,
251 const size_type src_c_j = dim - 1,
252 const unsigned int dst_r = 0,
253 const unsigned int dst_c = 0) const;
254
267 template <typename MatrixType, typename index_type>
268 void
269 extract_submatrix_from(const MatrixType & matrix,
270 const std::vector<index_type> &row_index_set,
271 const std::vector<index_type> &column_index_set);
272
285 template <typename MatrixType, typename index_type>
286 void
287 scatter_matrix_to(const std::vector<index_type> &row_index_set,
288 const std::vector<index_type> &column_index_set,
289 MatrixType & matrix) const;
290
301 template <typename number2>
302 void
304 const size_type dst_offset_i = 0,
305 const size_type dst_offset_j = 0,
306 const size_type src_offset_i = 0,
307 const size_type src_offset_j = 0);
308
309
313 template <typename number2>
314 void
315 fill(const number2 *);
316
328 template <typename number2>
329 void
331 const std::vector<size_type> &p_rows,
332 const std::vector<size_type> &p_cols);
333
344 void
345 set(const size_type i, const size_type j, const number value);
361 bool
363
369 m() const;
370
376 n() const;
377
383 bool
384 all_zero() const;
385
401 template <typename number2>
402 number2
404
414 template <typename number2>
415 number2
417 const Vector<number2> &v) const;
418
424 l1_norm() const;
425
431 linfty_norm() const;
432
442
453
459 number
460 determinant() const;
461
467 number
468 trace() const;
469
476 template <class StreamType>
477 void
478 print(StreamType & s,
479 const unsigned int width = 5,
480 const unsigned int precision = 2) const;
481
504 void
505 print_formatted(std::ostream & out,
506 const unsigned int precision = 3,
507 const bool scientific = true,
508 const unsigned int width = 0,
509 const char * zero_string = " ",
510 const double denominator = 1.,
511 const double threshold = 0.) const;
512
517 std::size_t
519
521
523
528 begin(const size_type r);
529
534 end(const size_type r);
535
540 begin(const size_type r) const;
541
546 end(const size_type r) const;
547
549
551
555 FullMatrix &
556 operator*=(const number factor);
557
561 FullMatrix &
562 operator/=(const number factor);
563
571 template <typename number2>
572 void
573 add(const number a, const FullMatrix<number2> &A);
574
582 template <typename number2>
583 void
584 add(const number a,
585 const FullMatrix<number2> &A,
586 const number b,
587 const FullMatrix<number2> &B);
588
597 template <typename number2>
598 void
599 add(const number a,
600 const FullMatrix<number2> &A,
601 const number b,
602 const FullMatrix<number2> &B,
603 const number c,
604 const FullMatrix<number2> &C);
605
617 template <typename number2>
618 void
620 const number factor,
621 const size_type dst_offset_i = 0,
622 const size_type dst_offset_j = 0,
623 const size_type src_offset_i = 0,
624 const size_type src_offset_j = 0);
625
631 template <typename number2>
632 void
633 Tadd(const number s, const FullMatrix<number2> &B);
634
646 template <typename number2>
647 void
649 const number factor,
650 const size_type dst_offset_i = 0,
651 const size_type dst_offset_j = 0,
652 const size_type src_offset_i = 0,
653 const size_type src_offset_j = 0);
654
658 void
659 add(const size_type row, const size_type column, const number value);
660
670 template <typename number2, typename index_type>
671 void
672 add(const size_type row,
673 const size_type n_cols,
674 const index_type *col_indices,
675 const number2 * values,
676 const bool elide_zero_values = true,
677 const bool col_indices_are_sorted = false);
678
682 void
683 add_row(const size_type i, const number s, const size_type j);
684
689 void
691 const number s,
692 const size_type j,
693 const number t,
694 const size_type k);
695
699 void
700 add_col(const size_type i, const number s, const size_type j);
701
706 void
708 const number s,
709 const size_type j,
710 const number t,
711 const size_type k);
712
716 void
717 swap_row(const size_type i, const size_type j);
718
722 void
723 swap_col(const size_type i, const size_type j);
724
729 void
730 diagadd(const number s);
731
735 template <typename number2>
736 void
737 equ(const number a, const FullMatrix<number2> &A);
738
742 template <typename number2>
743 void
744 equ(const number a,
745 const FullMatrix<number2> &A,
746 const number b,
747 const FullMatrix<number2> &B);
748
752 template <typename number2>
753 void
754 equ(const number a,
755 const FullMatrix<number2> &A,
756 const number b,
757 const FullMatrix<number2> &B,
758 const number c,
759 const FullMatrix<number2> &C);
760
767 void
769
784 void
786
793 template <typename number2>
794 void
796
805 template <typename number2>
806 void
808
813 template <typename number2>
814 void
816
822 template <typename number2>
823 void
825
831 template <typename number2>
832 void
834
836
838
857 template <typename number2>
858 void
860 const FullMatrix<number2> &B,
861 const bool adding = false) const;
862
881 template <typename number2>
882 void
884 const FullMatrix<number2> &B,
885 const bool adding = false) const;
886
905 template <typename number2>
906 void
908 const FullMatrix<number2> &B,
909 const bool adding = false) const;
910
930 template <typename number2>
931 void
933 const FullMatrix<number2> &B,
934 const bool adding = false) const;
935
946 void
948 const FullMatrix<number> &B,
949 const FullMatrix<number> &D,
950 const bool transpose_B = false,
951 const bool transpose_D = false,
952 const number scaling = number(1.));
953
966 template <typename number2>
967 void
969 const Vector<number2> &v,
970 const bool adding = false) const;
971
977 template <typename number2>
978 void
980
994 template <typename number2>
995 void
997 const Vector<number2> &v,
998 const bool adding = false) const;
999
1006 template <typename number2>
1007 void
1009
1015 template <typename somenumber>
1016 void
1018 const Vector<somenumber> &src,
1019 const number omega = 1.) const;
1020
1027 template <typename number2, typename number3>
1028 number
1030 const Vector<number2> &x,
1031 const Vector<number3> &b) const;
1032
1043 template <typename number2>
1044 void
1045 forward(Vector<number2> &dst, const Vector<number2> &src) const;
1046
1054 template <typename number2>
1055 void
1057
1059
1069
1075 number,
1076 << "The maximal pivot is " << arg1
1077 << ", which is below the threshold. The matrix may be singular.");
1082 size_type,
1083 size_type,
1084 size_type,
1085 << "Target region not in matrix: size in this direction="
1086 << arg1 << ", size of new matrix=" << arg2
1087 << ", offset=" << arg3);
1092 "You are attempting an operation on two matrices that "
1093 "are the same object, but the operation requires that the "
1094 "two objects are in fact different.");
1100};
1101
1104#ifndef DOXYGEN
1105/*-------------------------Inline functions -------------------------------*/
1106
1107
1108
1109template <typename number>
1110inline typename FullMatrix<number>::size_type
1112{
1113 return this->n_rows();
1114}
1115
1116
1117
1118template <typename number>
1119inline typename FullMatrix<number>::size_type
1121{
1122 return this->n_cols();
1123}
1124
1125
1126
1127template <typename number>
1129FullMatrix<number>::operator=(const number d)
1130{
1132 (void)d; // removes -Wunused-parameter warning in optimized mode
1133
1134 if (this->n_elements() != 0)
1135 this->reset_values();
1136
1137 return *this;
1138}
1139
1140
1141
1142template <typename number>
1143template <typename number2>
1144inline void
1145FullMatrix<number>::fill(const number2 *src)
1146{
1148}
1149
1150
1151
1152template <typename number>
1153template <typename MatrixType>
1154void
1155FullMatrix<number>::copy_from(const MatrixType &M)
1156{
1157 this->reinit(M.m(), M.n());
1158
1159 // loop over the elements of the argument matrix row by row, as suggested
1160 // in the documentation of the sparse matrix iterator class, and
1161 // copy them into the current object
1162 for (size_type row = 0; row < M.m(); ++row)
1163 {
1164 const typename MatrixType::const_iterator end_row = M.end(row);
1165 for (typename MatrixType::const_iterator entry = M.begin(row);
1166 entry != end_row;
1167 ++entry)
1168 this->el(row, entry->column()) = entry->value();
1169 }
1170}
1171
1172
1173
1174template <typename number>
1175template <typename MatrixType>
1176void
1177FullMatrix<number>::copy_transposed(const MatrixType &M)
1178{
1179 this->reinit(M.n(), M.m());
1180
1181 // loop over the elements of the argument matrix row by row, as suggested
1182 // in the documentation of the sparse matrix iterator class, and
1183 // copy them into the current object
1184 for (size_type row = 0; row < M.m(); ++row)
1185 {
1186 const typename MatrixType::const_iterator end_row = M.end(row);
1187 for (typename MatrixType::const_iterator entry = M.begin(row);
1188 entry != end_row;
1189 ++entry)
1190 this->el(entry->column(), row) = entry->value();
1191 }
1192}
1193
1194
1195
1196template <typename number>
1197template <typename MatrixType, typename index_type>
1198inline void
1200 const MatrixType & matrix,
1201 const std::vector<index_type> &row_index_set,
1202 const std::vector<index_type> &column_index_set)
1203{
1204 AssertDimension(row_index_set.size(), this->n_rows());
1205 AssertDimension(column_index_set.size(), this->n_cols());
1206
1207 const size_type n_rows_submatrix = row_index_set.size();
1208 const size_type n_cols_submatrix = column_index_set.size();
1209
1210 for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1211 for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1212 (*this)(sub_row, sub_col) =
1213 matrix.el(row_index_set[sub_row], column_index_set[sub_col]);
1214}
1215
1216
1217
1218template <typename number>
1219template <typename MatrixType, typename index_type>
1220inline void
1222 const std::vector<index_type> &row_index_set,
1223 const std::vector<index_type> &column_index_set,
1224 MatrixType & matrix) const
1225{
1226 AssertDimension(row_index_set.size(), this->n_rows());
1227 AssertDimension(column_index_set.size(), this->n_cols());
1228
1229 const size_type n_rows_submatrix = row_index_set.size();
1230 const size_type n_cols_submatrix = column_index_set.size();
1231
1232 for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1233 for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1234 matrix.set(row_index_set[sub_row],
1235 column_index_set[sub_col],
1236 (*this)(sub_row, sub_col));
1237}
1238
1239
1240template <typename number>
1241inline void
1243 const size_type j,
1244 const number value)
1245{
1246 (*this)(i, j) = value;
1247}
1248
1249
1250
1251template <typename number>
1252template <typename number2>
1253void
1255 const Vector<number2> &v) const
1256{
1257 vmult(w, v, true);
1258}
1259
1260
1261template <typename number>
1262template <typename number2>
1263void
1265 const Vector<number2> &v) const
1266{
1267 Tvmult(w, v, true);
1268}
1269
1270
1271//---------------------------------------------------------------------------
1272template <typename number>
1273inline typename FullMatrix<number>::iterator
1275{
1276 AssertIndexRange(r, m());
1277 return iterator(this, r, 0);
1278}
1279
1280
1281
1282template <typename number>
1283inline typename FullMatrix<number>::iterator
1285{
1286 AssertIndexRange(r, m());
1287 return iterator(this, r + 1, 0);
1288}
1289
1290
1291
1292template <typename number>
1295{
1296 AssertIndexRange(r, m());
1297 return const_iterator(this, r, 0);
1298}
1299
1300
1301
1302template <typename number>
1305{
1306 AssertIndexRange(r, m());
1307 return const_iterator(this, r + 1, 0);
1308}
1309
1310
1311
1312template <typename number>
1313inline void
1314FullMatrix<number>::add(const size_type r, const size_type c, const number v)
1315{
1316 AssertIndexRange(r, this->m());
1317 AssertIndexRange(c, this->n());
1318
1319 this->operator()(r, c) += v;
1320}
1321
1322
1323
1324template <typename number>
1325template <typename number2, typename index_type>
1326inline void
1328 const size_type n_cols,
1329 const index_type *col_indices,
1330 const number2 * values,
1331 const bool,
1332 const bool)
1333{
1334 AssertIndexRange(row, this->m());
1335 for (size_type col = 0; col < n_cols; ++col)
1336 {
1337 AssertIndexRange(col_indices[col], this->n());
1338 this->operator()(row, col_indices[col]) += values[col];
1339 }
1340}
1341
1342
1343template <typename number>
1344template <class StreamType>
1345inline void
1346FullMatrix<number>::print(StreamType & s,
1347 const unsigned int w,
1348 const unsigned int p) const
1349{
1350 Assert(!this->empty(), ExcEmptyMatrix());
1351
1352 // save the state of out stream
1353 const std::streamsize old_precision = s.precision(p);
1354 const std::streamsize old_width = s.width(w);
1355
1356 for (size_type i = 0; i < this->m(); ++i)
1357 {
1358 for (size_type j = 0; j < this->n(); ++j)
1359 {
1360 s.width(w);
1361 s.precision(p);
1362 s << this->el(i, j);
1363 }
1364 s << std::endl;
1365 }
1366
1367 // reset output format
1368 s.precision(old_precision);
1369 s.width(old_width);
1370}
1371
1372
1373#endif // DOXYGEN
1374
1376
1377#endif
typename numbers::NumberTraits< number >::real_type real_type
Definition: full_matrix.h:118
typename Table< 2, number >::const_iterator const_iterator
Definition: full_matrix.h:97
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< number > & operator=(const number d)
void TmTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
FullMatrix(const size_type rows, const size_type cols)
number residual(Vector< number2 > &dst, const Vector< number2 > &x, const Vector< number3 > &b) const
std::size_t memory_consumption() const
void diagadd(const number s)
void fill_permutation(const FullMatrix< number2 > &src, const std::vector< size_type > &p_rows, const std::vector< size_type > &p_cols)
void add_row(const size_type i, const number s, const size_type j)
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
real_type relative_symmetry_norm2() const
void add(const size_type row, const size_type column, const number value)
void copy_from(const Tensor< 2, dim > &T, const unsigned int src_r_i=0, const unsigned int src_r_j=dim - 1, const unsigned int src_c_i=0, const unsigned int src_c_j=dim - 1, const size_type dst_r=0, const size_type dst_c=0)
void symmetrize()
void equ(const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B)
number trace() const
void right_invert(const FullMatrix< number2 > &M)
FullMatrix< number > & operator=(const FullMatrix< number2 > &)
FullMatrix & operator/=(const number factor)
std::size_t size_type
Definition: full_matrix.h:81
void set(const size_type i, const size_type j, const number value)
size_type n() const
void add_row(const size_type i, const number s, const size_type j, const number t, const size_type k)
FullMatrix< number > & operator=(const IdentityMatrix &id)
number value_type
Definition: full_matrix.h:87
typename Table< 2, number >::iterator iterator
Definition: full_matrix.h:92
void add(const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B, const number c, const FullMatrix< number2 > &C)
void Tadd(const number s, const FullMatrix< number2 > &B)
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
bool all_zero() const
void add_col(const size_type i, const number s, const size_type j, const number t, const size_type k)
void scatter_matrix_to(const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set, MatrixType &matrix) const
void swap_row(const size_type i, const size_type j)
void add(const FullMatrix< number2 > &src, const number factor, 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 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
number2 matrix_norm_square(const Vector< number2 > &v) const
void equ(const number a, const FullMatrix< number2 > &A)
const_iterator end(const size_type r) const
FullMatrix(const size_type rows, const size_type cols, const number *entries)
bool operator==(const FullMatrix< number > &) const
void Tvmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void Tmmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
number2 matrix_scalar_product(const Vector< number2 > &u, const Vector< number2 > &v) const
void gauss_jordan()
void add(const size_type row, const size_type n_cols, const index_type *col_indices, const number2 *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
void add_col(const size_type i, const number s, const size_type j)
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
FullMatrix(const IdentityMatrix &id)
void vmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
void invert(const FullMatrix< number2 > &M)
void cholesky(const FullMatrix< number2 > &A)
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 add(const number a, const FullMatrix< number2 > &A)
number determinant() const
void equ(const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B, const number c, const FullMatrix< number2 > &C)
void fill(const number2 *)
void copy_transposed(const MatrixType &)
void print(StreamType &s, const unsigned int width=5, const unsigned int precision=2) const
FullMatrix< number > & operator=(const LAPACKFullMatrix< number2 > &)
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 & operator*=(const number factor)
void Tvmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
void mTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void left_invert(const FullMatrix< number2 > &M)
iterator begin(const size_type r)
iterator end(const size_type r)
void add(const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B)
void outer_product(const Vector< number2 > &V, const Vector< number2 > &W)
real_type frobenius_norm() const
void forward(Vector< number2 > &dst, const Vector< number2 > &src) const
void Tadd(const FullMatrix< number2 > &src, const number factor, 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(const size_type n=0)
size_type m() const
void swap_col(const size_type i, const size_type j)
void extract_submatrix_from(const MatrixType &matrix, const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set)
void backward(Vector< number2 > &dst, const Vector< number2 > &src) const
void copy_from(const MatrixType &)
const_iterator begin(const size_type r) const
real_type l1_norm() const
real_type linfty_norm() const
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcEmptyMatrix()
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
static ::ExceptionBase & ExcNotRegular(number arg1)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcMatrixNotPositiveDefinite()
static ::ExceptionBase & ExcSourceEqualsDestination()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:561
static ::ExceptionBase & ExcInvalidDestination(size_type arg1, size_type arg2, size_type arg3)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
@ matrix
Contents is actually a matrix.
static const char A
static const char T
static const char V
types::global_dof_index size_type
Definition: cuda_kernels.h:45
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618