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\}}\)
petsc_matrix_base.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2004 - 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_petsc_matrix_base_h
17# define dealii_petsc_matrix_base_h
18
19
20# include <deal.II/base/config.h>
21
22# ifdef DEAL_II_WITH_PETSC
23
25
31
32# include <petscmat.h>
33
34# include <cmath>
35# include <memory>
36# include <vector>
37
39
40// Forward declarations
41# ifndef DOXYGEN
42template <typename Matrix>
43class BlockMatrixBase;
44# endif
45
46
47namespace PETScWrappers
48{
49 // forward declarations
50 class MatrixBase;
51
52 namespace MatrixIterators
53 {
68 {
69 private:
74 {
75 public:
80
86 const size_type row,
87 const size_type index);
88
93 row() const;
94
99 index() const;
100
105 column() const;
106
110 PetscScalar
111 value() const;
112
121 int,
122 int,
123 int,
124 << "You tried to access row " << arg1
125 << " of a distributed matrix, but only rows " << arg2
126 << " through " << arg3
127 << " are stored locally and can be accessed.");
128
129 private:
134
139
144
157 std::shared_ptr<const std::vector<size_type>> colnum_cache;
158
162 std::shared_ptr<const std::vector<PetscScalar>> value_cache;
163
169 void
171
172 // Make enclosing class a friend.
173 friend class const_iterator;
174 };
175
176 public:
181
187 const size_type row,
188 const size_type index);
189
195
201
205 const Accessor &operator*() const;
206
210 const Accessor *operator->() const;
211
216 bool
221 bool
223
229 bool
230 operator<(const const_iterator &) const;
231
236 int,
237 int,
238 << "Attempt to access element " << arg2 << " of row "
239 << arg1 << " which doesn't have that many elements.");
240
241 private:
246 };
247
248 } // namespace MatrixIterators
249
250
282 class MatrixBase : public Subscriptor
283 {
284 public:
289
294
298 using value_type = PetscScalar;
299
303 MatrixBase();
304
310 MatrixBase(const MatrixBase &) = delete;
311
317 MatrixBase &
318 operator=(const MatrixBase &) = delete;
319
323 virtual ~MatrixBase() override;
324
334 MatrixBase &
335 operator=(const value_type d);
340 void
341 clear();
342
352 void
353 set(const size_type i, const size_type j, const PetscScalar value);
354
374 void
375 set(const std::vector<size_type> & indices,
376 const FullMatrix<PetscScalar> &full_matrix,
377 const bool elide_zero_values = false);
378
384 void
385 set(const std::vector<size_type> & row_indices,
386 const std::vector<size_type> & col_indices,
387 const FullMatrix<PetscScalar> &full_matrix,
388 const bool elide_zero_values = false);
389
404 void
405 set(const size_type row,
406 const std::vector<size_type> & col_indices,
407 const std::vector<PetscScalar> &values,
408 const bool elide_zero_values = false);
409
424 void
425 set(const size_type row,
426 const size_type n_cols,
427 const size_type * col_indices,
428 const PetscScalar *values,
429 const bool elide_zero_values = false);
430
440 void
441 add(const size_type i, const size_type j, const PetscScalar value);
442
462 void
463 add(const std::vector<size_type> & indices,
464 const FullMatrix<PetscScalar> &full_matrix,
465 const bool elide_zero_values = true);
466
472 void
473 add(const std::vector<size_type> & row_indices,
474 const std::vector<size_type> & col_indices,
475 const FullMatrix<PetscScalar> &full_matrix,
476 const bool elide_zero_values = true);
477
492 void
493 add(const size_type row,
494 const std::vector<size_type> & col_indices,
495 const std::vector<PetscScalar> &values,
496 const bool elide_zero_values = true);
497
512 void
513 add(const size_type row,
514 const size_type n_cols,
515 const size_type * col_indices,
516 const PetscScalar *values,
517 const bool elide_zero_values = true,
518 const bool col_indices_are_sorted = false);
519
536 void
537 clear_row(const size_type row, const PetscScalar new_diag_value = 0);
538
547 void
548 clear_rows(const std::vector<size_type> &rows,
549 const PetscScalar new_diag_value = 0);
550
562 void
563 compress(const VectorOperation::values operation);
564
576 PetscScalar
577 operator()(const size_type i, const size_type j) const;
578
586 PetscScalar
587 el(const size_type i, const size_type j) const;
588
598 PetscScalar
599 diag_element(const size_type i) const;
600
605 m() const;
606
611 n() const;
612
622 local_size() const;
623
632 std::pair<size_type, size_type>
633 local_range() const;
634
639 bool
640 in_local_range(const size_type index) const;
641
646 virtual const MPI_Comm &
648
655 n_nonzero_elements() const;
656
661 row_length(const size_type row) const;
662
670 PetscReal
671 l1_norm() const;
672
680 PetscReal
681 linfty_norm() const;
682
687 PetscReal
688 frobenius_norm() const;
689
690
710 PetscScalar
711 matrix_norm_square(const VectorBase &v) const;
712
713
727 PetscScalar
728 matrix_scalar_product(const VectorBase &u, const VectorBase &v) const;
729
734 PetscScalar
735 trace() const;
736
740 MatrixBase &
741 operator*=(const PetscScalar factor);
742
746 MatrixBase &
747 operator/=(const PetscScalar factor);
748
749
754 MatrixBase &
755 add(const PetscScalar factor, const MatrixBase &other);
756
757
764 MatrixBase &
765 add(const MatrixBase &other, const PetscScalar factor);
766
778 void
779 vmult(VectorBase &dst, const VectorBase &src) const;
780
793 void
794 Tvmult(VectorBase &dst, const VectorBase &src) const;
795
807 void
808 vmult_add(VectorBase &dst, const VectorBase &src) const;
809
822 void
823 Tvmult_add(VectorBase &dst, const VectorBase &src) const;
824
837 PetscScalar
838 residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const;
839
846 begin() const;
847
854 end() const;
855
865 begin(const size_type r) const;
866
876 end(const size_type r) const;
877
885 operator Mat() const;
886
892 Mat &
893 petsc_matrix();
894
898 void
899 transpose();
900
905 PetscBool
906 is_symmetric(const double tolerance = 1.e-12);
907
913 PetscBool
914 is_hermitian(const double tolerance = 1.e-12);
915
922 void
923 write_ascii(const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
924
932 void
933 print(std::ostream &out, const bool alternative_output = false) const;
934
938 std::size_t
939 memory_consumption() const;
940
945 "You are attempting an operation on two matrices that "
946 "are the same object, but the operation requires that the "
947 "two objects are in fact different.");
948
953 int,
954 int,
955 << "You tried to do a "
956 << (arg1 == 1 ? "'set'" : (arg1 == 2 ? "'add'" : "???"))
957 << " operation but the matrix is currently in "
958 << (arg2 == 1 ? "'set'" : (arg2 == 2 ? "'add'" : "???"))
959 << " mode. You first have to call 'compress()'.");
960
961 protected:
967
972
978 void
980
986 void
988
999 void
1006 void
1008
1024 void
1025 mmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const;
1026
1043 void
1044 Tmmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const;
1045
1046 private:
1061 mutable std::vector<PetscInt> column_indices;
1062
1071 mutable std::vector<PetscScalar> column_values;
1072
1073
1074 // To allow calling protected prepare_add() and prepare_set().
1075 template <class>
1076 friend class ::BlockMatrixBase;
1077 };
1078
1079
1080
1081# ifndef DOXYGEN
1082 // ---------------------- inline and template functions ---------------------
1083
1084
1085 namespace MatrixIterators
1086 {
1088 const size_type row,
1089 const size_type index)
1090 : matrix(const_cast<MatrixBase *>(matrix))
1091 , a_row(row)
1092 , a_index(index)
1093 {
1094 visit_present_row();
1095 }
1096
1097
1098
1101 {
1102 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1103 return a_row;
1104 }
1105
1106
1109 {
1110 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1111 return (*colnum_cache)[a_index];
1112 }
1113
1114
1117 {
1118 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1119 return a_index;
1120 }
1121
1122
1123 inline PetscScalar
1125 {
1126 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1127 return (*value_cache)[a_index];
1128 }
1129
1130
1131 inline const_iterator::const_iterator(const MatrixBase *matrix,
1132 const size_type row,
1133 const size_type index)
1134 : accessor(matrix, row, index)
1135 {}
1136
1137
1138
1139 inline const_iterator &
1141 {
1142 Assert(accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
1143
1144 ++accessor.a_index;
1145
1146 // if at end of line: do one step, then cycle until we find a
1147 // row with a nonzero number of entries
1148 if (accessor.a_index >= accessor.colnum_cache->size())
1149 {
1150 accessor.a_index = 0;
1151 ++accessor.a_row;
1152
1153 while ((accessor.a_row < accessor.matrix->m()) &&
1154 (accessor.a_row < accessor.matrix->local_range().second) &&
1155 (accessor.matrix->row_length(accessor.a_row) == 0))
1156 ++accessor.a_row;
1157
1158 accessor.visit_present_row();
1159 }
1160 return *this;
1161 }
1162
1163
1164 inline const_iterator
1166 {
1167 const const_iterator old_state = *this;
1168 ++(*this);
1169 return old_state;
1170 }
1171
1172
1173 inline const const_iterator::Accessor &const_iterator::operator*() const
1174 {
1175 return accessor;
1176 }
1177
1178
1179 inline const const_iterator::Accessor *const_iterator::operator->() const
1180 {
1181 return &accessor;
1182 }
1183
1184
1185 inline bool
1186 const_iterator::operator==(const const_iterator &other) const
1187 {
1188 return (accessor.a_row == other.accessor.a_row &&
1189 accessor.a_index == other.accessor.a_index);
1190 }
1191
1192
1193 inline bool
1194 const_iterator::operator!=(const const_iterator &other) const
1195 {
1196 return !(*this == other);
1197 }
1198
1199
1200 inline bool
1201 const_iterator::operator<(const const_iterator &other) const
1202 {
1203 return (accessor.row() < other.accessor.row() ||
1204 (accessor.row() == other.accessor.row() &&
1205 accessor.index() < other.accessor.index()));
1206 }
1207
1208 } // namespace MatrixIterators
1209
1210
1211
1212 // Inline the set() and add()
1213 // functions, since they will be
1214 // called frequently, and the
1215 // compiler can optimize away
1216 // some unnecessary loops when
1217 // the sizes are given at
1218 // compile time.
1219 inline void
1220 MatrixBase::set(const size_type i, const size_type j, const PetscScalar value)
1221 {
1222 AssertIsFinite(value);
1223
1224 set(i, 1, &j, &value, false);
1225 }
1226
1227
1228
1229 inline void
1230 MatrixBase::set(const std::vector<size_type> & indices,
1232 const bool elide_zero_values)
1233 {
1234 Assert(indices.size() == values.m(),
1235 ExcDimensionMismatch(indices.size(), values.m()));
1236 Assert(values.m() == values.n(), ExcNotQuadratic());
1237
1238 for (size_type i = 0; i < indices.size(); ++i)
1239 set(indices[i],
1240 indices.size(),
1241 indices.data(),
1242 &values(i, 0),
1243 elide_zero_values);
1244 }
1245
1246
1247
1248 inline void
1249 MatrixBase::set(const std::vector<size_type> & row_indices,
1250 const std::vector<size_type> & col_indices,
1252 const bool elide_zero_values)
1253 {
1254 Assert(row_indices.size() == values.m(),
1255 ExcDimensionMismatch(row_indices.size(), values.m()));
1256 Assert(col_indices.size() == values.n(),
1257 ExcDimensionMismatch(col_indices.size(), values.n()));
1258
1259 for (size_type i = 0; i < row_indices.size(); ++i)
1260 set(row_indices[i],
1261 col_indices.size(),
1262 col_indices.data(),
1263 &values(i, 0),
1264 elide_zero_values);
1265 }
1266
1267
1268
1269 inline void
1270 MatrixBase::set(const size_type row,
1271 const std::vector<size_type> & col_indices,
1272 const std::vector<PetscScalar> &values,
1273 const bool elide_zero_values)
1274 {
1275 Assert(col_indices.size() == values.size(),
1276 ExcDimensionMismatch(col_indices.size(), values.size()));
1277
1278 set(row,
1279 col_indices.size(),
1280 col_indices.data(),
1281 values.data(),
1282 elide_zero_values);
1283 }
1284
1285
1286
1287 inline void
1288 MatrixBase::set(const size_type row,
1289 const size_type n_cols,
1290 const size_type * col_indices,
1291 const PetscScalar *values,
1292 const bool elide_zero_values)
1293 {
1294 prepare_action(VectorOperation::insert);
1295
1296 const PetscInt petsc_i = row;
1297 PetscInt const *col_index_ptr;
1298
1299 PetscScalar const *col_value_ptr;
1300 int n_columns;
1301
1302 // If we don't elide zeros, the pointers are already available...
1303 if (elide_zero_values == false)
1304 {
1305 col_index_ptr = reinterpret_cast<const PetscInt *>(col_indices);
1306 col_value_ptr = values;
1307 n_columns = n_cols;
1308 }
1309 else
1310 {
1311 // Otherwise, extract nonzero values in each row and get the
1312 // respective index.
1313 if (column_indices.size() < n_cols)
1314 {
1315 column_indices.resize(n_cols);
1316 column_values.resize(n_cols);
1317 }
1318
1319 n_columns = 0;
1320 for (size_type j = 0; j < n_cols; ++j)
1321 {
1322 const PetscScalar value = values[j];
1323 AssertIsFinite(value);
1324 if (value != PetscScalar())
1325 {
1326 column_indices[n_columns] = col_indices[j];
1327 column_values[n_columns] = value;
1328 n_columns++;
1329 }
1330 }
1331 AssertIndexRange(n_columns, n_cols + 1);
1332
1333 col_index_ptr = column_indices.data();
1334 col_value_ptr = column_values.data();
1335 }
1336
1337 const PetscErrorCode ierr = MatSetValues(matrix,
1338 1,
1339 &petsc_i,
1340 n_columns,
1341 col_index_ptr,
1342 col_value_ptr,
1343 INSERT_VALUES);
1344 AssertThrow(ierr == 0, ExcPETScError(ierr));
1345 }
1346
1347
1348
1349 inline void
1350 MatrixBase::add(const size_type i, const size_type j, const PetscScalar value)
1351 {
1352 AssertIsFinite(value);
1353
1354 if (value == PetscScalar())
1355 {
1356 // we have to check after using Insert/Add in any case to be
1357 // consistent with the MPI communication model, but we can save
1358 // some work if the addend is zero. However, these actions are done
1359 // in case we pass on to the other function.
1360 prepare_action(VectorOperation::add);
1361
1362 return;
1363 }
1364 else
1365 add(i, 1, &j, &value, false);
1366 }
1367
1368
1369
1370 inline void
1371 MatrixBase::add(const std::vector<size_type> & indices,
1373 const bool elide_zero_values)
1374 {
1375 Assert(indices.size() == values.m(),
1376 ExcDimensionMismatch(indices.size(), values.m()));
1377 Assert(values.m() == values.n(), ExcNotQuadratic());
1378
1379 for (size_type i = 0; i < indices.size(); ++i)
1380 add(indices[i],
1381 indices.size(),
1382 indices.data(),
1383 &values(i, 0),
1384 elide_zero_values);
1385 }
1386
1387
1388
1389 inline void
1390 MatrixBase::add(const std::vector<size_type> & row_indices,
1391 const std::vector<size_type> & col_indices,
1393 const bool elide_zero_values)
1394 {
1395 Assert(row_indices.size() == values.m(),
1396 ExcDimensionMismatch(row_indices.size(), values.m()));
1397 Assert(col_indices.size() == values.n(),
1398 ExcDimensionMismatch(col_indices.size(), values.n()));
1399
1400 for (size_type i = 0; i < row_indices.size(); ++i)
1401 add(row_indices[i],
1402 col_indices.size(),
1403 col_indices.data(),
1404 &values(i, 0),
1405 elide_zero_values);
1406 }
1407
1408
1409
1410 inline void
1411 MatrixBase::add(const size_type row,
1412 const std::vector<size_type> & col_indices,
1413 const std::vector<PetscScalar> &values,
1414 const bool elide_zero_values)
1415 {
1416 Assert(col_indices.size() == values.size(),
1417 ExcDimensionMismatch(col_indices.size(), values.size()));
1418
1419 add(row,
1420 col_indices.size(),
1421 col_indices.data(),
1422 values.data(),
1423 elide_zero_values);
1424 }
1425
1426
1427
1428 inline void
1429 MatrixBase::add(const size_type row,
1430 const size_type n_cols,
1431 const size_type * col_indices,
1432 const PetscScalar *values,
1433 const bool elide_zero_values,
1434 const bool /*col_indices_are_sorted*/)
1435 {
1436 (void)elide_zero_values;
1437
1438 prepare_action(VectorOperation::add);
1439
1440 const PetscInt petsc_i = row;
1441 PetscInt const *col_index_ptr;
1442
1443 PetscScalar const *col_value_ptr;
1444 int n_columns;
1445
1446 // If we don't elide zeros, the pointers are already available...
1447 if (elide_zero_values == false)
1448 {
1449 col_index_ptr = reinterpret_cast<const PetscInt *>(col_indices);
1450 col_value_ptr = values;
1451 n_columns = n_cols;
1452 }
1453 else
1454 {
1455 // Otherwise, extract nonzero values in each row and get the
1456 // respective index.
1457 if (column_indices.size() < n_cols)
1458 {
1459 column_indices.resize(n_cols);
1460 column_values.resize(n_cols);
1461 }
1462
1463 n_columns = 0;
1464 for (size_type j = 0; j < n_cols; ++j)
1465 {
1466 const PetscScalar value = values[j];
1467 AssertIsFinite(value);
1468 if (value != PetscScalar())
1469 {
1470 column_indices[n_columns] = col_indices[j];
1471 column_values[n_columns] = value;
1472 n_columns++;
1473 }
1474 }
1475 AssertIndexRange(n_columns, n_cols + 1);
1476
1477 col_index_ptr = column_indices.data();
1478 col_value_ptr = column_values.data();
1479 }
1480
1481 const PetscErrorCode ierr = MatSetValues(
1482 matrix, 1, &petsc_i, n_columns, col_index_ptr, col_value_ptr, ADD_VALUES);
1483 AssertThrow(ierr == 0, ExcPETScError(ierr));
1484 }
1485
1486
1487
1488 inline PetscScalar
1489 MatrixBase::operator()(const size_type i, const size_type j) const
1490 {
1491 return el(i, j);
1492 }
1493
1494
1495
1496 inline MatrixBase::const_iterator
1497 MatrixBase::begin() const
1498 {
1499 Assert(
1500 (in_local_range(0) && in_local_range(m() - 1)),
1501 ExcMessage(
1502 "begin() and end() can only be called on a processor owning the entire matrix. If this is a distributed matrix, use begin(row) and end(row) instead."));
1503
1504 // find the first non-empty row in order to make sure that the returned
1505 // iterator points to something useful
1506 size_type first_nonempty_row = 0;
1507 while ((first_nonempty_row < m()) && (row_length(first_nonempty_row) == 0))
1508 ++first_nonempty_row;
1509
1510 return const_iterator(this, first_nonempty_row, 0);
1511 }
1512
1513
1514 inline MatrixBase::const_iterator
1515 MatrixBase::end() const
1516 {
1517 Assert(
1518 (in_local_range(0) && in_local_range(m() - 1)),
1519 ExcMessage(
1520 "begin() and end() can only be called on a processor owning the entire matrix. If this is a distributed matrix, use begin(row) and end(row) instead."));
1521
1522 return const_iterator(this, m(), 0);
1523 }
1524
1525
1526 inline MatrixBase::const_iterator
1527 MatrixBase::begin(const size_type r) const
1528 {
1529 Assert(in_local_range(r),
1530 ExcIndexRange(r, local_range().first, local_range().second));
1531
1532 if (row_length(r) > 0)
1533 return const_iterator(this, r, 0);
1534 else
1535 return end(r);
1536 }
1537
1538
1539 inline MatrixBase::const_iterator
1540 MatrixBase::end(const size_type r) const
1541 {
1542 Assert(in_local_range(r),
1543 ExcIndexRange(r, local_range().first, local_range().second));
1544
1545 // place the iterator on the first entry past this line, or at the
1546 // end of the matrix
1547 //
1548 // in the parallel case, we need to put it on the first entry of
1549 // the first row after the locally owned range. this of course
1550 // doesn't exist, but we can nevertheless create such an
1551 // iterator. we need to check whether 'i' is past the locally
1552 // owned range of rows first, before we ask for the length of the
1553 // row since the latter query leads to an exception in case we ask
1554 // for a row that is not locally owned
1555 for (size_type i = r + 1; i < m(); ++i)
1556 if (i == local_range().second || (row_length(i) > 0))
1557 return const_iterator(this, i, 0);
1558
1559 // if there is no such line, then take the
1560 // end iterator of the matrix
1561 // we don't allow calling end() directly for distributed matrices so we need
1562 // to copy the code without the assertion.
1563 return {this, m(), 0};
1564 }
1565
1566
1567
1568 inline bool
1569 MatrixBase::in_local_range(const size_type index) const
1570 {
1571 PetscInt begin, end;
1572
1573 const PetscErrorCode ierr =
1574 MatGetOwnershipRange(static_cast<const Mat &>(matrix), &begin, &end);
1575 AssertThrow(ierr == 0, ExcPETScError(ierr));
1576
1577 return ((index >= static_cast<size_type>(begin)) &&
1578 (index < static_cast<size_type>(end)));
1579 }
1580
1581
1582
1583 inline void
1584 MatrixBase::prepare_action(const VectorOperation::values new_action)
1585 {
1586 if (last_action == VectorOperation::unknown)
1587 last_action = new_action;
1588
1589 Assert(last_action == new_action, ExcWrongMode(last_action, new_action));
1590 }
1591
1592
1593
1594 inline void
1595 MatrixBase::assert_is_compressed()
1596 {
1597 // compress() sets the last action to none, which allows us to check if
1598 // there are pending add/insert operations:
1600 ExcMessage("Error: missing compress() call."));
1601 }
1602
1603
1604
1605 inline void
1606 MatrixBase::prepare_add()
1607 {
1608 prepare_action(VectorOperation::add);
1609 }
1610
1611
1612
1613 inline void
1614 MatrixBase::prepare_set()
1615 {
1616 prepare_action(VectorOperation::insert);
1617 }
1618
1619# endif // DOXYGEN
1620} // namespace PETScWrappers
1621
1622
1624
1625
1626# endif // DEAL_II_WITH_PETSC
1627
1628#endif
1629/*---------------------------- petsc_matrix_base.h --------------------------*/
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
void add(const size_type i, const size_type j, const PetscScalar value)
void add(const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=true)
const_iterator end(const size_type r) const
size_type row_length(const size_type row) const
std::size_t memory_consumption() const
VectorOperation::values last_action
void vmult(VectorBase &dst, const VectorBase &src) const
void clear_rows(const std::vector< size_type > &rows, const PetscScalar new_diag_value=0)
std::vector< PetscScalar > column_values
PetscScalar diag_element(const size_type i) const
MatrixBase(const MatrixBase &)=delete
const_iterator begin() const
MatrixBase & operator/=(const PetscScalar factor)
void mmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
PetscBool is_symmetric(const double tolerance=1.e-12)
PetscReal frobenius_norm() const
virtual ~MatrixBase() override
const_iterator end() const
void Tvmult_add(VectorBase &dst, const VectorBase &src) const
void print(std::ostream &out, const bool alternative_output=false) const
const_iterator begin(const size_type r) const
void add(const size_type row, const std::vector< size_type > &col_indices, const std::vector< PetscScalar > &values, const bool elide_zero_values=true)
PetscScalar el(const size_type i, const size_type j) const
PetscScalar operator()(const size_type i, const size_type j) const
MatrixBase & operator=(const MatrixBase &)=delete
void set(const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=false)
void prepare_action(const VectorOperation::values new_action)
void add(const std::vector< size_type > &indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=true)
bool in_local_range(const size_type index) const
PetscBool is_hermitian(const double tolerance=1.e-12)
PetscScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
std::vector< PetscInt > column_indices
void Tvmult(VectorBase &dst, const VectorBase &src) const
MatrixBase & operator*=(const PetscScalar factor)
PetscScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
void Tmmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
void set(const size_type row, const std::vector< size_type > &col_indices, const std::vector< PetscScalar > &values, const bool elide_zero_values=false)
void add(const size_type row, const size_type n_cols, const size_type *col_indices, const PetscScalar *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
void vmult_add(VectorBase &dst, const VectorBase &src) const
std::pair< size_type, size_type > local_range() const
size_type n_nonzero_elements() const
void compress(const VectorOperation::values operation)
PetscScalar matrix_norm_square(const VectorBase &v) const
void set(const size_type row, const size_type n_cols, const size_type *col_indices, const PetscScalar *values, const bool elide_zero_values=false)
virtual const MPI_Comm & get_mpi_communicator() const =0
void set(const std::vector< size_type > &indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=false)
void set(const size_type i, const size_type j, const PetscScalar value)
void clear_row(const size_type row, const PetscScalar new_diag_value=0)
std::shared_ptr< const std::vector< PetscScalar > > value_cache
std::shared_ptr< const std::vector< size_type > > colnum_cache
Accessor(const MatrixBase *matrix, const size_type row, const size_type index)
bool operator!=(const const_iterator &) const
const_iterator(const MatrixBase *matrix, const size_type row, const size_type index)
bool operator==(const const_iterator &) const
bool operator<(const const_iterator &) const
std::enable_if< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typenameProductType< std::complex< T >, std::complex< U > >::type >::type operator*(const std::complex< T > &left, const std::complex< U > &right)
#define DEAL_II_DEPRECATED
Definition: config.h:162
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
Point< 2 > second
Definition: grid_out.cc:4588
Point< 2 > first
Definition: grid_out.cc:4587
__global__ void set(Number *val, const Number s, const size_type N)
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcAccessToNonlocalRow(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcIteratorPastEnd()
static ::ExceptionBase & ExcInvalidIndexWithinRow(int arg1, int arg2)
#define AssertIsFinite(number)
Definition: exceptions.h:1721
static ::ExceptionBase & ExcBeyondEndOfMatrix()
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcWrongMode(int arg1, int arg2)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:561
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
@ matrix
Contents is actually a matrix.
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)
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)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int global_dof_index
Definition: types.h:76
SynchronousIterators< Iterators > operator++(SynchronousIterators< Iterators > &a)
bool operator<(const SynchronousIterators< Iterators > &a, const SynchronousIterators< Iterators > &b)