Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
full_matrix.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1999 - 2023 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
28
29#include <cstring>
30#include <iomanip>
31#include <vector>
32
34
35
36// forward declarations
37#ifndef DOXYGEN
38template <typename number>
39class Vector;
40template <typename number>
42#endif
43
77template <typename number>
78class FullMatrix : public Table<2, number>
79{
80public:
88 static_assert(
89 std::is_arithmetic<
91 "The FullMatrix class only supports basic numeric types. In particular, it "
92 "does not support automatically differentiated numbers.");
93
94
98 using size_type = std::size_t;
99
104 using value_type = number;
105
110
115
119 using Table<2, number>::begin;
120
124 using Table<2, number>::end;
125
136
151 explicit FullMatrix(const size_type n = 0);
152
156 FullMatrix(const size_type rows, const size_type cols);
157
162 FullMatrix(const size_type rows, const size_type cols, const number *entries);
163
185 template <typename number2>
188
197 operator=(const number d);
198
209
214 template <typename number2>
217
218
224 template <typename MatrixType>
225 void
226 copy_from(const MatrixType &);
227
233 template <typename MatrixType>
234 void
235 copy_transposed(const MatrixType &);
236
244 template <int dim>
245 void
247 const unsigned int src_r_i = 0,
248 const unsigned int src_r_j = dim - 1,
249 const unsigned int src_c_i = 0,
250 const unsigned int src_c_j = dim - 1,
251 const size_type dst_r = 0,
252 const size_type dst_c = 0);
253
261 template <int dim>
262 void
264 const size_type src_r_i = 0,
265 const size_type src_r_j = dim - 1,
266 const size_type src_c_i = 0,
267 const size_type src_c_j = dim - 1,
268 const unsigned int dst_r = 0,
269 const unsigned int dst_c = 0) const;
270
283 template <typename MatrixType, typename index_type>
284 void
285 extract_submatrix_from(const MatrixType & matrix,
286 const std::vector<index_type> &row_index_set,
287 const std::vector<index_type> &column_index_set);
288
301 template <typename MatrixType, typename index_type>
302 void
303 scatter_matrix_to(const std::vector<index_type> &row_index_set,
304 const std::vector<index_type> &column_index_set,
305 MatrixType & matrix) const;
306
317 template <typename number2>
318 void
320 const size_type dst_offset_i = 0,
321 const size_type dst_offset_j = 0,
322 const size_type src_offset_i = 0,
323 const size_type src_offset_j = 0);
324
325
329 template <typename number2>
330 void
331 fill(const number2 *);
332
344 template <typename number2>
345 void
347 const std::vector<size_type> &p_rows,
348 const std::vector<size_type> &p_cols);
349
360 void
361 set(const size_type i, const size_type j, const number value);
375 bool
377
383 m() const;
384
390 n() const;
391
397 bool
398 all_zero() const;
399
415 template <typename number2>
416 number2
418
428 template <typename number2>
429 number2
431 const Vector<number2> &v) const;
432
438 l1_norm() const;
439
445 linfty_norm() const;
446
456
467
473 number
474 determinant() const;
475
481 number
482 trace() const;
483
490 template <class StreamType>
491 void
492 print(StreamType & s,
493 const unsigned int width = 5,
494 const unsigned int precision = 2) const;
495
518 void
519 print_formatted(std::ostream & out,
520 const unsigned int precision = 3,
521 const bool scientific = true,
522 const unsigned int width = 0,
523 const char * zero_string = " ",
524 const double denominator = 1.,
525 const double threshold = 0.) const;
526
531 std::size_t
533
544 begin(const size_type r);
545
550 end(const size_type r);
551
556 begin(const size_type r) const;
557
562 end(const size_type r) const;
563
573 FullMatrix &
574 operator*=(const number factor);
575
579 FullMatrix &
580 operator/=(const number factor);
581
589 template <typename number2>
590 void
591 add(const number a, const FullMatrix<number2> &A);
592
600 template <typename number2>
601 void
602 add(const number a,
603 const FullMatrix<number2> &A,
604 const number b,
605 const FullMatrix<number2> &B);
606
615 template <typename number2>
616 void
617 add(const number a,
618 const FullMatrix<number2> &A,
619 const number b,
620 const FullMatrix<number2> &B,
621 const number c,
622 const FullMatrix<number2> &C);
623
635 template <typename number2>
636 void
638 const number factor,
639 const size_type dst_offset_i = 0,
640 const size_type dst_offset_j = 0,
641 const size_type src_offset_i = 0,
642 const size_type src_offset_j = 0);
643
649 template <typename number2>
650 void
651 Tadd(const number s, const FullMatrix<number2> &B);
652
664 template <typename number2>
665 void
667 const number factor,
668 const size_type dst_offset_i = 0,
669 const size_type dst_offset_j = 0,
670 const size_type src_offset_i = 0,
671 const size_type src_offset_j = 0);
672
676 void
677 add(const size_type row, const size_type column, const number value);
678
688 template <typename number2, typename index_type>
689 void
690 add(const size_type row,
691 const size_type n_cols,
692 const index_type *col_indices,
693 const number2 * values,
694 const bool elide_zero_values = true,
695 const bool col_indices_are_sorted = false);
696
700 void
701 add_row(const size_type i, const number s, const size_type j);
702
707 void
709 const number s,
710 const size_type j,
711 const number t,
712 const size_type k);
713
717 void
718 add_col(const size_type i, const number s, const size_type j);
719
724 void
726 const number s,
727 const size_type j,
728 const number t,
729 const size_type k);
730
734 void
735 swap_row(const size_type i, const size_type j);
736
740 void
741 swap_col(const size_type i, const size_type j);
742
747 void
748 diagadd(const number s);
749
753 template <typename number2>
754 void
755 equ(const number a, const FullMatrix<number2> &A);
756
760 template <typename number2>
761 void
762 equ(const number a,
763 const FullMatrix<number2> &A,
764 const number b,
765 const FullMatrix<number2> &B);
766
770 template <typename number2>
771 void
772 equ(const number a,
773 const FullMatrix<number2> &A,
774 const number b,
775 const FullMatrix<number2> &B,
776 const number c,
777 const FullMatrix<number2> &C);
778
785 void
787
802 void
804
811 template <typename number2>
812 void
814
823 template <typename number2>
824 void
826
831 template <typename number2>
832 void
834
840 template <typename number2>
841 void
843
849 template <typename number2>
850 void
852
877 template <typename number2>
878 void
880 const FullMatrix<number2> &B,
881 const bool adding = false) const;
882
901 template <typename number2>
902 void
904 const FullMatrix<number2> &B,
905 const bool adding = false) const;
906
925 template <typename number2>
926 void
928 const FullMatrix<number2> &B,
929 const bool adding = false) const;
930
950 template <typename number2>
951 void
953 const FullMatrix<number2> &B,
954 const bool adding = false) const;
955
966 void
968 const FullMatrix<number> &B,
969 const FullMatrix<number> &D,
970 const bool transpose_B = false,
971 const bool transpose_D = false,
972 const number scaling = number(1.));
973
986 template <typename number2>
987 void
989 const Vector<number2> &v,
990 const bool adding = false) const;
991
997 template <typename number2>
998 void
1000
1014 template <typename number2>
1015 void
1017 const Vector<number2> &v,
1018 const bool adding = false) const;
1019
1026 template <typename number2>
1027 void
1029
1035 template <typename somenumber>
1036 void
1038 const Vector<somenumber> &src,
1039 const number omega = 1.) const;
1040
1047 template <typename number2, typename number3>
1048 number
1050 const Vector<number2> &x,
1051 const Vector<number3> &b) const;
1052
1063 template <typename number2>
1064 void
1065 forward(Vector<number2> &dst, const Vector<number2> &src) const;
1066
1074 template <typename number2>
1075 void
1077
1089
1095 number,
1096 << "The maximal pivot is " << arg1
1097 << ", which is below the threshold. The matrix may be singular.");
1102 size_type,
1103 size_type,
1104 size_type,
1105 << "Target region not in matrix: size in this direction="
1106 << arg1 << ", size of new matrix=" << arg2
1107 << ", offset=" << arg3);
1112 "You are attempting an operation on two vectors that "
1113 "are the same object, but the operation requires that the "
1114 "two objects are in fact different.");
1120};
1121
1124#ifndef DOXYGEN
1125/*-------------------------Inline functions -------------------------------*/
1126
1127
1128
1129template <typename number>
1130inline typename FullMatrix<number>::size_type
1132{
1133 return this->n_rows();
1134}
1135
1136
1137
1138template <typename number>
1139inline typename FullMatrix<number>::size_type
1141{
1142 return this->n_cols();
1143}
1144
1145
1146
1147template <typename number>
1149FullMatrix<number>::operator=(const number d)
1150{
1152 (void)d; // removes -Wunused-parameter warning in optimized mode
1153
1154 if (this->n_elements() != 0)
1155 this->reset_values();
1156
1157 return *this;
1158}
1159
1160
1161
1162template <typename number>
1163template <typename number2>
1164inline void
1165FullMatrix<number>::fill(const number2 *src)
1166{
1168}
1169
1170
1171
1172template <typename number>
1173template <typename MatrixType>
1174void
1175FullMatrix<number>::copy_from(const MatrixType &M)
1176{
1177 this->reinit(M.m(), M.n());
1178
1179 // loop over the elements of the argument matrix row by row, as suggested
1180 // in the documentation of the sparse matrix iterator class, and
1181 // copy them into the current object
1182 for (size_type row = 0; row < M.m(); ++row)
1183 {
1184 const typename MatrixType::const_iterator end_row = M.end(row);
1185 for (typename MatrixType::const_iterator entry = M.begin(row);
1186 entry != end_row;
1187 ++entry)
1188 this->el(row, entry->column()) = entry->value();
1189 }
1190}
1191
1192
1193
1194template <typename number>
1195template <int dim>
1196void
1198 const unsigned int src_r_i,
1199 const unsigned int src_r_j,
1200 const unsigned int src_c_i,
1201 const unsigned int src_c_j,
1202 const size_type dst_r,
1203 const size_type dst_c)
1204{
1205 Assert(!this->empty(), ExcEmptyMatrix());
1206 AssertIndexRange(src_r_j - src_r_i, this->m() - dst_r);
1207 AssertIndexRange(src_c_j - src_c_i, this->n() - dst_c);
1208 AssertIndexRange(src_r_j, dim);
1209 AssertIndexRange(src_c_j, dim);
1210 AssertIndexRange(src_r_i, src_r_j + 1);
1211 AssertIndexRange(src_c_i, src_c_j + 1);
1212
1213 for (size_type i = 0; i < src_r_j - src_r_i + 1; ++i)
1214 for (size_type j = 0; j < src_c_j - src_c_i + 1; ++j)
1215 {
1216 const unsigned int src_r_index = static_cast<unsigned int>(i + src_r_i);
1217 const unsigned int src_c_index = static_cast<unsigned int>(j + src_c_i);
1218 (*this)(i + dst_r, j + dst_c) = number(T[src_r_index][src_c_index]);
1219 }
1220}
1221
1222
1223
1224template <typename number>
1225template <int dim>
1226void
1228 const size_type src_r_i,
1229 const size_type src_r_j,
1230 const size_type src_c_i,
1231 const size_type src_c_j,
1232 const unsigned int dst_r,
1233 const unsigned int dst_c) const
1234{
1235 Assert(!this->empty(), ExcEmptyMatrix());
1236 AssertIndexRange(src_r_j - src_r_i, dim - dst_r);
1237 AssertIndexRange(src_c_j - src_c_i, dim - dst_c);
1238 AssertIndexRange(src_r_j, this->m());
1239 AssertIndexRange(src_r_j, this->n());
1240 AssertIndexRange(src_r_i, src_r_j + 1);
1241 AssertIndexRange(src_c_j, src_c_j + 1);
1242
1243 for (size_type i = 0; i < src_r_j - src_r_i + 1; ++i)
1244 for (size_type j = 0; j < src_c_j - src_c_i + 1; ++j)
1245 {
1246 const unsigned int dst_r_index = static_cast<unsigned int>(i + dst_r);
1247 const unsigned int dst_c_index = static_cast<unsigned int>(j + dst_c);
1248 T[dst_r_index][dst_c_index] = double((*this)(i + src_r_i, j + src_c_i));
1249 }
1250}
1251
1252
1253
1254template <typename number>
1255template <typename MatrixType>
1256void
1257FullMatrix<number>::copy_transposed(const MatrixType &M)
1258{
1259 this->reinit(M.n(), M.m());
1260
1261 // loop over the elements of the argument matrix row by row, as suggested
1262 // in the documentation of the sparse matrix iterator class, and
1263 // copy them into the current object
1264 for (size_type row = 0; row < M.m(); ++row)
1265 {
1266 const typename MatrixType::const_iterator end_row = M.end(row);
1267 for (typename MatrixType::const_iterator entry = M.begin(row);
1268 entry != end_row;
1269 ++entry)
1270 this->el(entry->column(), row) = entry->value();
1271 }
1272}
1273
1274
1275
1276template <typename number>
1277template <typename MatrixType, typename index_type>
1278inline void
1280 const MatrixType & matrix,
1281 const std::vector<index_type> &row_index_set,
1282 const std::vector<index_type> &column_index_set)
1283{
1284 AssertDimension(row_index_set.size(), this->n_rows());
1285 AssertDimension(column_index_set.size(), this->n_cols());
1286
1287 const size_type n_rows_submatrix = row_index_set.size();
1288 const size_type n_cols_submatrix = column_index_set.size();
1289
1290 for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1291 for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1292 (*this)(sub_row, sub_col) =
1293 matrix.el(row_index_set[sub_row], column_index_set[sub_col]);
1294}
1295
1296
1297
1298template <typename number>
1299template <typename MatrixType, typename index_type>
1300inline void
1302 const std::vector<index_type> &row_index_set,
1303 const std::vector<index_type> &column_index_set,
1304 MatrixType & matrix) const
1305{
1306 AssertDimension(row_index_set.size(), this->n_rows());
1307 AssertDimension(column_index_set.size(), this->n_cols());
1308
1309 const size_type n_rows_submatrix = row_index_set.size();
1310 const size_type n_cols_submatrix = column_index_set.size();
1311
1312 for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1313 for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1314 matrix.set(row_index_set[sub_row],
1315 column_index_set[sub_col],
1316 (*this)(sub_row, sub_col));
1317}
1318
1319
1320template <typename number>
1321inline void
1322FullMatrix<number>::set(const size_type i,
1323 const size_type j,
1324 const number value)
1325{
1326 (*this)(i, j) = value;
1327}
1328
1329
1330
1331template <typename number>
1332template <typename number2>
1333void
1335 const Vector<number2> &v) const
1336{
1337 vmult(w, v, true);
1338}
1339
1340
1341template <typename number>
1342template <typename number2>
1343void
1345 const Vector<number2> &v) const
1346{
1347 Tvmult(w, v, true);
1348}
1349
1350
1351//---------------------------------------------------------------------------
1352template <typename number>
1353inline typename FullMatrix<number>::iterator
1354FullMatrix<number>::begin(const size_type r)
1355{
1356 AssertIndexRange(r, m());
1357 return iterator(this, r, 0);
1358}
1359
1360
1361
1362template <typename number>
1363inline typename FullMatrix<number>::iterator
1364FullMatrix<number>::end(const size_type r)
1365{
1366 AssertIndexRange(r, m());
1367 return iterator(this, r + 1, 0);
1368}
1369
1370
1371
1372template <typename number>
1374FullMatrix<number>::begin(const size_type r) const
1375{
1376 AssertIndexRange(r, m());
1377 return const_iterator(this, r, 0);
1378}
1379
1380
1381
1382template <typename number>
1384FullMatrix<number>::end(const size_type r) const
1385{
1386 AssertIndexRange(r, m());
1387 return const_iterator(this, r + 1, 0);
1388}
1389
1390
1391
1392template <typename number>
1393inline void
1394FullMatrix<number>::add(const size_type r, const size_type c, const number v)
1395{
1396 AssertIndexRange(r, this->m());
1397 AssertIndexRange(c, this->n());
1398
1399 this->operator()(r, c) += v;
1400}
1401
1402
1403
1404template <typename number>
1405template <typename number2, typename index_type>
1406inline void
1407FullMatrix<number>::add(const size_type row,
1408 const size_type n_cols,
1409 const index_type *col_indices,
1410 const number2 * values,
1411 const bool,
1412 const bool)
1413{
1414 AssertIndexRange(row, this->m());
1415 for (size_type col = 0; col < n_cols; ++col)
1416 {
1417 AssertIndexRange(col_indices[col], this->n());
1418 this->operator()(row, col_indices[col]) += values[col];
1419 }
1420}
1421
1422
1423template <typename number>
1424template <class StreamType>
1425inline void
1426FullMatrix<number>::print(StreamType & s,
1427 const unsigned int w,
1428 const unsigned int p) const
1429{
1430 Assert(!this->empty(), ExcEmptyMatrix());
1431
1432 // save the state of out stream
1433 const std::streamsize old_precision = s.precision(p);
1434 const std::streamsize old_width = s.width(w);
1435
1436 for (size_type i = 0; i < this->m(); ++i)
1437 {
1438 for (size_type j = 0; j < this->n(); ++j)
1439 {
1440 s.width(w);
1441 s.precision(p);
1442 s << this->el(i, j);
1443 }
1444 s << std::endl;
1445 }
1446
1447 // reset output format
1448 s.precision(old_precision);
1449 s.width(old_width);
1450}
1451
1452
1453#endif // DOXYGEN
1454
1456
1457#endif
typename numbers::NumberTraits< number >::real_type real_type
typename Table< 2, number >::const_iterator const_iterator
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:98
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
typename Table< 2, number >::iterator iterator
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
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
#define DeclException0(Exception0)
Definition exceptions.h:465
static ::ExceptionBase & ExcEmptyMatrix()
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcNotRegular(number arg1)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:488
static ::ExceptionBase & ExcMatrixNotPositiveDefinite()
static ::ExceptionBase & ExcSourceEqualsDestination()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition exceptions.h:556
static ::ExceptionBase & ExcInvalidDestination(size_type arg1, size_type arg2, size_type arg3)
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:510
@ matrix
Contents is actually a matrix.
types::global_dof_index size_type
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)