Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
petsc_matrix_base.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_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 
26 # include <deal.II/lac/exceptions.h>
27 # include <deal.II/lac/full_matrix.h>
31 
32 # include <petscmat.h>
33 
34 # include <cmath>
35 # include <memory>
36 # include <vector>
37 
39 
40 // Forward declarations
41 # ifndef DOXYGEN
42 template <typename Matrix>
43 class BlockMatrixBase;
44 # endif
45 
46 
47 namespace PETScWrappers
48 {
49  // forward declarations
50  class MatrixBase;
51 
52  namespace MatrixIterators
53  {
69  {
70  private:
74  class Accessor
75  {
76  public:
81 
86  Accessor(const MatrixBase *matrix,
87  const size_type row,
88  const size_type index);
89 
93  size_type
94  row() const;
95 
99  size_type
100  index() const;
101 
105  size_type
106  column() const;
107 
111  PetscScalar
112  value() const;
113 
122  int,
123  int,
124  int,
125  << "You tried to access row " << arg1
126  << " of a distributed matrix, but only rows " << arg2
127  << " through " << arg3
128  << " are stored locally and can be accessed.");
129 
130  private:
134  mutable MatrixBase *matrix;
135 
140 
145 
158  std::shared_ptr<const std::vector<size_type>> colnum_cache;
159 
163  std::shared_ptr<const std::vector<PetscScalar>> value_cache;
164 
170  void
172 
173  // Make enclosing class a friend.
174  friend class const_iterator;
175  };
176 
177  public:
182 
188  const size_type row,
189  const size_type index);
190 
195  operator++();
196 
201  operator++(int);
202 
206  const Accessor &operator*() const;
207 
211  const Accessor *operator->() const;
212 
217  bool
218  operator==(const const_iterator &) const;
222  bool
223  operator!=(const const_iterator &) const;
224 
230  bool
231  operator<(const const_iterator &) const;
232 
237  int,
238  int,
239  << "Attempt to access element " << arg2 << " of row "
240  << arg1 << " which doesn't have that many elements.");
241 
242  private:
247  };
248 
249  } // namespace MatrixIterators
250 
251 
284  class MatrixBase : public Subscriptor
285  {
286  public:
291 
296 
300  using value_type = PetscScalar;
301 
305  MatrixBase();
306 
312  MatrixBase(const MatrixBase &) = delete;
313 
319  MatrixBase &
320  operator=(const MatrixBase &) = delete;
321 
325  virtual ~MatrixBase() override;
326 
336  MatrixBase &
337  operator=(const value_type d);
342  void
343  clear();
344 
354  void
355  set(const size_type i, const size_type j, const PetscScalar value);
356 
376  void
377  set(const std::vector<size_type> & indices,
378  const FullMatrix<PetscScalar> &full_matrix,
379  const bool elide_zero_values = false);
380 
386  void
387  set(const std::vector<size_type> & row_indices,
388  const std::vector<size_type> & col_indices,
389  const FullMatrix<PetscScalar> &full_matrix,
390  const bool elide_zero_values = false);
391 
406  void
407  set(const size_type row,
408  const std::vector<size_type> & col_indices,
409  const std::vector<PetscScalar> &values,
410  const bool elide_zero_values = false);
411 
426  void
427  set(const size_type row,
428  const size_type n_cols,
429  const size_type * col_indices,
430  const PetscScalar *values,
431  const bool elide_zero_values = false);
432 
442  void
443  add(const size_type i, const size_type j, const PetscScalar value);
444 
464  void
465  add(const std::vector<size_type> & indices,
466  const FullMatrix<PetscScalar> &full_matrix,
467  const bool elide_zero_values = true);
468 
474  void
475  add(const std::vector<size_type> & row_indices,
476  const std::vector<size_type> & col_indices,
477  const FullMatrix<PetscScalar> &full_matrix,
478  const bool elide_zero_values = true);
479 
494  void
495  add(const size_type row,
496  const std::vector<size_type> & col_indices,
497  const std::vector<PetscScalar> &values,
498  const bool elide_zero_values = true);
499 
514  void
515  add(const size_type row,
516  const size_type n_cols,
517  const size_type * col_indices,
518  const PetscScalar *values,
519  const bool elide_zero_values = true,
520  const bool col_indices_are_sorted = false);
521 
538  void
539  clear_row(const size_type row, const PetscScalar new_diag_value = 0);
540 
549  void
550  clear_rows(const std::vector<size_type> &rows,
551  const PetscScalar new_diag_value = 0);
552 
564  void
565  compress(const VectorOperation::values operation);
566 
578  PetscScalar
579  operator()(const size_type i, const size_type j) const;
580 
588  PetscScalar
589  el(const size_type i, const size_type j) const;
590 
600  PetscScalar
601  diag_element(const size_type i) const;
602 
606  size_type
607  m() const;
608 
612  size_type
613  n() const;
614 
623  size_type
624  local_size() const;
625 
634  std::pair<size_type, size_type>
635  local_range() const;
636 
641  bool
642  in_local_range(const size_type index) const;
643 
648  virtual const MPI_Comm &
649  get_mpi_communicator() const = 0;
650 
656  size_type
657  n_nonzero_elements() const;
658 
662  size_type
663  row_length(const size_type row) const;
664 
672  PetscReal
673  l1_norm() const;
674 
682  PetscReal
683  linfty_norm() const;
684 
689  PetscReal
690  frobenius_norm() const;
691 
692 
712  PetscScalar
713  matrix_norm_square(const VectorBase &v) const;
714 
715 
729  PetscScalar
730  matrix_scalar_product(const VectorBase &u, const VectorBase &v) const;
731 
736  PetscScalar
737  trace() const;
738 
742  MatrixBase &
743  operator*=(const PetscScalar factor);
744 
748  MatrixBase &
749  operator/=(const PetscScalar factor);
750 
751 
756  MatrixBase &
757  add(const PetscScalar factor, const MatrixBase &other);
758 
759 
766  MatrixBase &
767  add(const MatrixBase &other, const PetscScalar factor);
768 
780  void
781  vmult(VectorBase &dst, const VectorBase &src) const;
782 
795  void
796  Tvmult(VectorBase &dst, const VectorBase &src) const;
797 
809  void
810  vmult_add(VectorBase &dst, const VectorBase &src) const;
811 
824  void
825  Tvmult_add(VectorBase &dst, const VectorBase &src) const;
826 
839  PetscScalar
840  residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const;
841 
848  begin() const;
849 
856  end() const;
857 
867  begin(const size_type r) const;
868 
878  end(const size_type r) const;
879 
887  operator Mat() const;
888 
894  Mat &
895  petsc_matrix();
896 
900  void
901  transpose();
902 
907  PetscBool
908  is_symmetric(const double tolerance = 1.e-12);
909 
915  PetscBool
916  is_hermitian(const double tolerance = 1.e-12);
917 
924  void
925  write_ascii(const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
926 
934  void
935  print(std::ostream &out, const bool alternative_output = false) const;
936 
940  std::size_t
941  memory_consumption() const;
942 
947  "You are attempting an operation on two matrices that "
948  "are the same object, but the operation requires that the "
949  "two objects are in fact different.");
950 
955  int,
956  int,
957  << "You tried to do a "
958  << (arg1 == 1 ? "'set'" : (arg1 == 2 ? "'add'" : "???"))
959  << " operation but the matrix is currently in "
960  << (arg2 == 1 ? "'set'" : (arg2 == 2 ? "'add'" : "???"))
961  << " mode. You first have to call 'compress()'.");
962 
963  protected:
968  Mat matrix;
969 
974 
980  void
981  prepare_action(const VectorOperation::values new_action);
982 
988  void
990 
1001  void
1002  prepare_add();
1008  void
1009  prepare_set();
1010 
1026  void
1027  mmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const;
1028 
1045  void
1046  Tmmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const;
1047 
1048  private:
1063  mutable std::vector<PetscInt> column_indices;
1064 
1073  mutable std::vector<PetscScalar> column_values;
1074 
1075 
1076  // To allow calling protected prepare_add() and prepare_set().
1077  template <class>
1078  friend class ::BlockMatrixBase;
1079  };
1080 
1081 
1082 
1083 # ifndef DOXYGEN
1084  // ---------------------- inline and template functions ---------------------
1085 
1086 
1087  namespace MatrixIterators
1088  {
1090  const size_type row,
1091  const size_type index)
1092  : matrix(const_cast<MatrixBase *>(matrix))
1093  , a_row(row)
1094  , a_index(index)
1095  {
1096  visit_present_row();
1097  }
1098 
1099 
1100 
1103  {
1104  Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1105  return a_row;
1106  }
1107 
1108 
1111  {
1112  Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1113  return (*colnum_cache)[a_index];
1114  }
1115 
1116 
1119  {
1120  Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1121  return a_index;
1122  }
1123 
1124 
1125  inline PetscScalar
1127  {
1128  Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
1129  return (*value_cache)[a_index];
1130  }
1131 
1132 
1133  inline const_iterator::const_iterator(const MatrixBase *matrix,
1134  const size_type row,
1135  const size_type index)
1136  : accessor(matrix, row, index)
1137  {}
1138 
1139 
1140 
1141  inline const_iterator &
1143  {
1144  Assert(accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
1145 
1146  ++accessor.a_index;
1147 
1148  // if at end of line: do one step, then cycle until we find a
1149  // row with a nonzero number of entries
1150  if (accessor.a_index >= accessor.colnum_cache->size())
1151  {
1152  accessor.a_index = 0;
1153  ++accessor.a_row;
1154 
1155  while ((accessor.a_row < accessor.matrix->m()) &&
1156  (accessor.a_row < accessor.matrix->local_range().second) &&
1157  (accessor.matrix->row_length(accessor.a_row) == 0))
1158  ++accessor.a_row;
1159 
1160  accessor.visit_present_row();
1161  }
1162  return *this;
1163  }
1164 
1165 
1166  inline const_iterator
1168  {
1169  const const_iterator old_state = *this;
1170  ++(*this);
1171  return old_state;
1172  }
1173 
1174 
1176  {
1177  return accessor;
1178  }
1179 
1180 
1181  inline const const_iterator::Accessor *const_iterator::operator->() const
1182  {
1183  return &accessor;
1184  }
1185 
1186 
1187  inline bool
1188  const_iterator::operator==(const const_iterator &other) const
1189  {
1190  return (accessor.a_row == other.accessor.a_row &&
1191  accessor.a_index == other.accessor.a_index);
1192  }
1193 
1194 
1195  inline bool
1196  const_iterator::operator!=(const const_iterator &other) const
1197  {
1198  return !(*this == other);
1199  }
1200 
1201 
1202  inline bool
1203  const_iterator::operator<(const const_iterator &other) const
1204  {
1205  return (accessor.row() < other.accessor.row() ||
1206  (accessor.row() == other.accessor.row() &&
1207  accessor.index() < other.accessor.index()));
1208  }
1209 
1210  } // namespace MatrixIterators
1211 
1212 
1213 
1214  // Inline the set() and add()
1215  // functions, since they will be
1216  // called frequently, and the
1217  // compiler can optimize away
1218  // some unnecessary loops when
1219  // the sizes are given at
1220  // compile time.
1221  inline void
1222  MatrixBase::set(const size_type i, const size_type j, const PetscScalar value)
1223  {
1225 
1226  set(i, 1, &j, &value, false);
1227  }
1228 
1229 
1230 
1231  inline void
1232  MatrixBase::set(const std::vector<size_type> & indices,
1233  const FullMatrix<PetscScalar> &values,
1234  const bool elide_zero_values)
1235  {
1236  Assert(indices.size() == values.m(),
1237  ExcDimensionMismatch(indices.size(), values.m()));
1238  Assert(values.m() == values.n(), ExcNotQuadratic());
1239 
1240  for (size_type i = 0; i < indices.size(); ++i)
1241  set(indices[i],
1242  indices.size(),
1243  indices.data(),
1244  &values(i, 0),
1245  elide_zero_values);
1246  }
1247 
1248 
1249 
1250  inline void
1251  MatrixBase::set(const std::vector<size_type> & row_indices,
1252  const std::vector<size_type> & col_indices,
1253  const FullMatrix<PetscScalar> &values,
1254  const bool elide_zero_values)
1255  {
1256  Assert(row_indices.size() == values.m(),
1257  ExcDimensionMismatch(row_indices.size(), values.m()));
1258  Assert(col_indices.size() == values.n(),
1259  ExcDimensionMismatch(col_indices.size(), values.n()));
1260 
1261  for (size_type i = 0; i < row_indices.size(); ++i)
1262  set(row_indices[i],
1263  col_indices.size(),
1264  col_indices.data(),
1265  &values(i, 0),
1266  elide_zero_values);
1267  }
1268 
1269 
1270 
1271  inline void
1272  MatrixBase::set(const size_type row,
1273  const std::vector<size_type> & col_indices,
1274  const std::vector<PetscScalar> &values,
1275  const bool elide_zero_values)
1276  {
1277  Assert(col_indices.size() == values.size(),
1278  ExcDimensionMismatch(col_indices.size(), values.size()));
1279 
1280  set(row,
1281  col_indices.size(),
1282  col_indices.data(),
1283  values.data(),
1284  elide_zero_values);
1285  }
1286 
1287 
1288 
1289  inline void
1290  MatrixBase::set(const size_type row,
1291  const size_type n_cols,
1292  const size_type * col_indices,
1293  const PetscScalar *values,
1294  const bool elide_zero_values)
1295  {
1296  prepare_action(VectorOperation::insert);
1297 
1298  const PetscInt petsc_i = row;
1299  PetscInt const *col_index_ptr;
1300 
1301  PetscScalar const *col_value_ptr;
1302  int n_columns;
1303 
1304  // If we don't elide zeros, the pointers are already available...
1305  if (elide_zero_values == false)
1306  {
1307  col_index_ptr = reinterpret_cast<const PetscInt *>(col_indices);
1308  col_value_ptr = values;
1309  n_columns = n_cols;
1310  }
1311  else
1312  {
1313  // Otherwise, extract nonzero values in each row and get the
1314  // respective index.
1315  if (column_indices.size() < n_cols)
1316  {
1317  column_indices.resize(n_cols);
1318  column_values.resize(n_cols);
1319  }
1320 
1321  n_columns = 0;
1322  for (size_type j = 0; j < n_cols; ++j)
1323  {
1324  const PetscScalar value = values[j];
1326  if (value != PetscScalar())
1327  {
1328  column_indices[n_columns] = col_indices[j];
1329  column_values[n_columns] = value;
1330  n_columns++;
1331  }
1332  }
1333  AssertIndexRange(n_columns, n_cols + 1);
1334 
1335  col_index_ptr = column_indices.data();
1336  col_value_ptr = column_values.data();
1337  }
1338 
1339  const PetscErrorCode ierr = MatSetValues(matrix,
1340  1,
1341  &petsc_i,
1342  n_columns,
1343  col_index_ptr,
1344  col_value_ptr,
1345  INSERT_VALUES);
1346  AssertThrow(ierr == 0, ExcPETScError(ierr));
1347  }
1348 
1349 
1350 
1351  inline void
1352  MatrixBase::add(const size_type i, const size_type j, const PetscScalar value)
1353  {
1355 
1356  if (value == PetscScalar())
1357  {
1358  // we have to check after using Insert/Add in any case to be
1359  // consistent with the MPI communication model, but we can save
1360  // some work if the addend is zero. However, these actions are done
1361  // in case we pass on to the other function.
1362  prepare_action(VectorOperation::add);
1363 
1364  return;
1365  }
1366  else
1367  add(i, 1, &j, &value, false);
1368  }
1369 
1370 
1371 
1372  inline void
1373  MatrixBase::add(const std::vector<size_type> & indices,
1374  const FullMatrix<PetscScalar> &values,
1375  const bool elide_zero_values)
1376  {
1377  Assert(indices.size() == values.m(),
1378  ExcDimensionMismatch(indices.size(), values.m()));
1379  Assert(values.m() == values.n(), ExcNotQuadratic());
1380 
1381  for (size_type i = 0; i < indices.size(); ++i)
1382  add(indices[i],
1383  indices.size(),
1384  indices.data(),
1385  &values(i, 0),
1386  elide_zero_values);
1387  }
1388 
1389 
1390 
1391  inline void
1392  MatrixBase::add(const std::vector<size_type> & row_indices,
1393  const std::vector<size_type> & col_indices,
1394  const FullMatrix<PetscScalar> &values,
1395  const bool elide_zero_values)
1396  {
1397  Assert(row_indices.size() == values.m(),
1398  ExcDimensionMismatch(row_indices.size(), values.m()));
1399  Assert(col_indices.size() == values.n(),
1400  ExcDimensionMismatch(col_indices.size(), values.n()));
1401 
1402  for (size_type i = 0; i < row_indices.size(); ++i)
1403  add(row_indices[i],
1404  col_indices.size(),
1405  col_indices.data(),
1406  &values(i, 0),
1407  elide_zero_values);
1408  }
1409 
1410 
1411 
1412  inline void
1413  MatrixBase::add(const size_type row,
1414  const std::vector<size_type> & col_indices,
1415  const std::vector<PetscScalar> &values,
1416  const bool elide_zero_values)
1417  {
1418  Assert(col_indices.size() == values.size(),
1419  ExcDimensionMismatch(col_indices.size(), values.size()));
1420 
1421  add(row,
1422  col_indices.size(),
1423  col_indices.data(),
1424  values.data(),
1425  elide_zero_values);
1426  }
1427 
1428 
1429 
1430  inline void
1431  MatrixBase::add(const size_type row,
1432  const size_type n_cols,
1433  const size_type * col_indices,
1434  const PetscScalar *values,
1435  const bool elide_zero_values,
1436  const bool /*col_indices_are_sorted*/)
1437  {
1438  (void)elide_zero_values;
1439 
1440  prepare_action(VectorOperation::add);
1441 
1442  const PetscInt petsc_i = row;
1443  PetscInt const *col_index_ptr;
1444 
1445  PetscScalar const *col_value_ptr;
1446  int n_columns;
1447 
1448  // If we don't elide zeros, the pointers are already available...
1449  if (elide_zero_values == false)
1450  {
1451  col_index_ptr = reinterpret_cast<const PetscInt *>(col_indices);
1452  col_value_ptr = values;
1453  n_columns = n_cols;
1454  }
1455  else
1456  {
1457  // Otherwise, extract nonzero values in each row and get the
1458  // respective index.
1459  if (column_indices.size() < n_cols)
1460  {
1461  column_indices.resize(n_cols);
1462  column_values.resize(n_cols);
1463  }
1464 
1465  n_columns = 0;
1466  for (size_type j = 0; j < n_cols; ++j)
1467  {
1468  const PetscScalar value = values[j];
1470  if (value != PetscScalar())
1471  {
1472  column_indices[n_columns] = col_indices[j];
1473  column_values[n_columns] = value;
1474  n_columns++;
1475  }
1476  }
1477  AssertIndexRange(n_columns, n_cols + 1);
1478 
1479  col_index_ptr = column_indices.data();
1480  col_value_ptr = column_values.data();
1481  }
1482 
1483  const PetscErrorCode ierr = MatSetValues(
1484  matrix, 1, &petsc_i, n_columns, col_index_ptr, col_value_ptr, ADD_VALUES);
1485  AssertThrow(ierr == 0, ExcPETScError(ierr));
1486  }
1487 
1488 
1489 
1490  inline PetscScalar
1491  MatrixBase::operator()(const size_type i, const size_type j) const
1492  {
1493  return el(i, j);
1494  }
1495 
1496 
1497 
1498  inline MatrixBase::const_iterator
1499  MatrixBase::begin() const
1500  {
1501  Assert(
1502  (in_local_range(0) && in_local_range(m() - 1)),
1503  ExcMessage(
1504  "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."));
1505 
1506  // find the first non-empty row in order to make sure that the returned
1507  // iterator points to something useful
1508  size_type first_nonempty_row = 0;
1509  while ((first_nonempty_row < m()) && (row_length(first_nonempty_row) == 0))
1510  ++first_nonempty_row;
1511 
1512  return const_iterator(this, first_nonempty_row, 0);
1513  }
1514 
1515 
1516  inline MatrixBase::const_iterator
1517  MatrixBase::end() const
1518  {
1519  Assert(
1520  (in_local_range(0) && in_local_range(m() - 1)),
1521  ExcMessage(
1522  "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."));
1523 
1524  return const_iterator(this, m(), 0);
1525  }
1526 
1527 
1528  inline MatrixBase::const_iterator
1529  MatrixBase::begin(const size_type r) const
1530  {
1531  Assert(in_local_range(r),
1532  ExcIndexRange(r, local_range().first, local_range().second));
1533 
1534  if (row_length(r) > 0)
1535  return const_iterator(this, r, 0);
1536  else
1537  return end(r);
1538  }
1539 
1540 
1541  inline MatrixBase::const_iterator
1542  MatrixBase::end(const size_type r) const
1543  {
1544  Assert(in_local_range(r),
1545  ExcIndexRange(r, local_range().first, local_range().second));
1546 
1547  // place the iterator on the first entry past this line, or at the
1548  // end of the matrix
1549  //
1550  // in the parallel case, we need to put it on the first entry of
1551  // the first row after the locally owned range. this of course
1552  // doesn't exist, but we can nevertheless create such an
1553  // iterator. we need to check whether 'i' is past the locally
1554  // owned range of rows first, before we ask for the length of the
1555  // row since the latter query leads to an exception in case we ask
1556  // for a row that is not locally owned
1557  for (size_type i = r + 1; i < m(); ++i)
1558  if (i == local_range().second || (row_length(i) > 0))
1559  return const_iterator(this, i, 0);
1560 
1561  // if there is no such line, then take the
1562  // end iterator of the matrix
1563  // we don't allow calling end() directly for distributed matrices so we need
1564  // to copy the code without the assertion.
1565  return {this, m(), 0};
1566  }
1567 
1568 
1569 
1570  inline bool
1571  MatrixBase::in_local_range(const size_type index) const
1572  {
1573  PetscInt begin, end;
1574 
1575  const PetscErrorCode ierr =
1576  MatGetOwnershipRange(static_cast<const Mat &>(matrix), &begin, &end);
1577  AssertThrow(ierr == 0, ExcPETScError(ierr));
1578 
1579  return ((index >= static_cast<size_type>(begin)) &&
1580  (index < static_cast<size_type>(end)));
1581  }
1582 
1583 
1584 
1585  inline void
1586  MatrixBase::prepare_action(const VectorOperation::values new_action)
1587  {
1588  if (last_action == VectorOperation::unknown)
1589  last_action = new_action;
1590 
1591  Assert(last_action == new_action, ExcWrongMode(last_action, new_action));
1592  }
1593 
1594 
1595 
1596  inline void
1597  MatrixBase::assert_is_compressed()
1598  {
1599  // compress() sets the last action to none, which allows us to check if
1600  // there are pending add/insert operations:
1601  AssertThrow(last_action == VectorOperation::unknown,
1602  ExcMessage("Error: missing compress() call."));
1603  }
1604 
1605 
1606 
1607  inline void
1608  MatrixBase::prepare_add()
1609  {
1610  prepare_action(VectorOperation::add);
1611  }
1612 
1613 
1614 
1615  inline void
1616  MatrixBase::prepare_set()
1617  {
1618  prepare_action(VectorOperation::insert);
1619  }
1620 
1621 # endif // DOXYGEN
1622 } // namespace PETScWrappers
1623 
1624 
1626 
1627 
1628 # endif // DEAL_II_WITH_PETSC
1629 
1630 #endif
1631 /*---------------------------- petsc_matrix_base.h --------------------------*/
PETScWrappers::MatrixBase::memory_consumption
std::size_t memory_consumption() const
Definition: petsc_matrix_base.cc:688
PETScWrappers::MatrixBase::column_indices
std::vector< PetscInt > column_indices
Definition: petsc_matrix_base.h:1063
DeclExceptionMsg
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:496
operator!=
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
Definition: aligned_vector.h:1192
PETScWrappers::MatrixIterators::const_iterator::size_type
types::global_dof_index size_type
Definition: petsc_matrix_base.h:181
PETScWrappers::MatrixIterators::const_iterator::Accessor::colnum_cache
std::shared_ptr< const std::vector< size_type > > colnum_cache
Definition: petsc_matrix_base.h:158
operator++
SynchronousIterators< Iterators > operator++(SynchronousIterators< Iterators > &a)
Definition: synchronous_iterator.h:248
operator<
bool operator<(const SynchronousIterators< Iterators > &a, const SynchronousIterators< Iterators > &b)
Definition: synchronous_iterator.h:113
PETScWrappers
Definition: petsc_block_sparse_matrix.h:36
petsc_compatibility.h
PETScWrappers::MatrixBase::Tmmult
void Tmmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
Definition: petsc_matrix_base.cc:572
PETScWrappers::MatrixBase::vmult_add
void vmult_add(VectorBase &dst, const VectorBase &src) const
Definition: petsc_matrix_base.cc:471
VectorOperation::insert
@ insert
Definition: vector_operation.h:49
PETScWrappers::MatrixBase::in_local_range
bool in_local_range(const size_type index) const
PETScWrappers::MatrixIterators::const_iterator::const_iterator
const_iterator(const MatrixBase *matrix, const size_type row, const size_type index)
petsc_vector_base.h
PETScWrappers::MatrixIterators::const_iterator::Accessor::ExcBeyondEndOfMatrix
static ::ExceptionBase & ExcBeyondEndOfMatrix()
operator==
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
Definition: aligned_vector.h:1170
TransposeTableIterators::Accessor
MatrixTableIterators::Accessor< TransposeTable< T >, Constness, MatrixTableIterators::Storage::column_major > Accessor
Definition: table.h:1907
PETScWrappers::MatrixBase::mmult
void mmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
Definition: petsc_matrix_base.cc:564
PETScWrappers::MatrixBase::matrix_scalar_product
PetscScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
Definition: petsc_matrix_base.cc:385
PETScWrappers::MatrixIterators::const_iterator::operator*
const Accessor & operator*() const
PETScWrappers::MatrixBase::matrix_norm_square
PetscScalar matrix_norm_square(const VectorBase &v) const
Definition: petsc_matrix_base.cc:376
FullMatrix::m
size_type m() const
LAPACKSupport::V
static const char V
Definition: lapack_support.h:175
PETScWrappers::MatrixBase::n_nonzero_elements
size_type n_nonzero_elements() const
Definition: petsc_matrix_base.cc:290
PETScWrappers::MatrixIterators::const_iterator::Accessor::size_type
types::global_dof_index size_type
Definition: petsc_matrix_base.h:80
PETScWrappers::FullMatrix
Definition: petsc_full_matrix.h:49
PETScWrappers::MatrixBase::clear
void clear()
Definition: petsc_matrix_base.cc:92
PETScWrappers::MatrixIterators::const_iterator::ExcInvalidIndexWithinRow
static ::ExceptionBase & ExcInvalidIndexWithinRow(int arg1, int arg2)
PETScWrappers::MatrixBase::local_size
size_type local_size() const
Definition: petsc_matrix_base.cc:263
PETScWrappers::MatrixBase::l1_norm
PetscReal l1_norm() const
Definition: petsc_matrix_base.cc:338
PETScWrappers::MatrixBase::set
void set(const size_type i, const size_type j, const PetscScalar value)
PETScWrappers::MatrixBase::n
size_type n() const
Definition: petsc_matrix_base.cc:250
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Physics::Elasticity::Kinematics::C
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
MPI_Comm
VectorOperation::add
@ add
Definition: vector_operation.h:53
Physics::Elasticity::Kinematics::e
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
PETScWrappers::MatrixBase::prepare_add
void prepare_add()
PETScWrappers::MatrixIterators::const_iterator::operator++
const_iterator & operator++()
FullMatrix::n
size_type n() const
exceptions.h
PETScWrappers::MatrixBase::print
void print(std::ostream &out, const bool alternative_output=false) const
Definition: petsc_matrix_base.cc:657
PETScWrappers::MatrixIterators::const_iterator::operator!=
bool operator!=(const const_iterator &) const
PETScWrappers::MatrixBase::is_symmetric
PetscBool is_symmetric(const double tolerance=1.e-12)
Definition: petsc_matrix_base.cc:620
PETScWrappers::MatrixBase::prepare_action
void prepare_action(const VectorOperation::values new_action)
second
Point< 2 > second
Definition: grid_out.cc:4353
PETScWrappers::MatrixIterators::const_iterator::Accessor::index
size_type index() const
PETScWrappers::MatrixIterators::const_iterator::Accessor::value_cache
std::shared_ptr< const std::vector< PetscScalar > > value_cache
Definition: petsc_matrix_base.h:163
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
PETScWrappers::MatrixIterators::const_iterator::Accessor::Accessor
Accessor(const MatrixBase *matrix, const size_type row, const size_type index)
PETScWrappers::MatrixBase::petsc_matrix
Mat & petsc_matrix()
Definition: petsc_matrix_base.cc:602
PETScWrappers::MatrixBase::operator*=
MatrixBase & operator*=(const PetscScalar factor)
Definition: petsc_matrix_base.cc:408
LinearAlgebra::CUDAWrappers::kernel::set
__global__ void set(Number *val, const Number s, const size_type N)
PETScWrappers::MatrixBase::vmult
void vmult(VectorBase &dst, const VectorBase &src) const
Definition: petsc_matrix_base.cc:449
PETScWrappers::MatrixBase::column_values
std::vector< PetscScalar > column_values
Definition: petsc_matrix_base.h:1073
PETScWrappers::MatrixBase::frobenius_norm
PetscReal frobenius_norm() const
Definition: petsc_matrix_base.cc:364
LACExceptions::ExcNotQuadratic
static ::ExceptionBase & ExcNotQuadratic()
PETScWrappers::MatrixBase::diag_element
PetscScalar diag_element(const size_type i) const
Definition: petsc_matrix_base.cc:181
Subscriptor
Definition: subscriptor.h:62
PETScWrappers::MatrixBase::row_length
size_type row_length(const size_type row) const
Definition: petsc_matrix_base.cc:302
PETScWrappers::MatrixBase::trace
PetscScalar trace() const
Definition: petsc_matrix_base.cc:395
StandardExceptions::ExcIndexRange
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
PETScWrappers::MatrixIterators::const_iterator::Accessor::value
PetscScalar value() const
PETScWrappers::MatrixIterators::const_iterator::Accessor::column
size_type column() const
AssertIsFinite
#define AssertIsFinite(number)
Definition: exceptions.h:1681
PETScWrappers::MatrixBase::Tvmult
void Tvmult(VectorBase &dst, const VectorBase &src) const
Definition: petsc_matrix_base.cc:460
vector_operation.h
BlockMatrixBase
Definition: affine_constraints.h:1903
PETScWrappers::MatrixBase::write_ascii
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
Definition: petsc_matrix_base.cc:642
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
types::global_dof_index
unsigned int global_dof_index
Definition: types.h:76
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
PETScWrappers::MatrixIterators::const_iterator::operator==
bool operator==(const const_iterator &) const
PETScWrappers::MatrixIterators::const_iterator::Accessor
Definition: petsc_matrix_base.h:74
PETScWrappers::MatrixIterators::const_iterator::Accessor::visit_present_row
void visit_present_row()
subscriptor.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
PETScWrappers::MatrixBase::local_range
std::pair< size_type, size_type > local_range() const
Definition: petsc_matrix_base.cc:276
PETScWrappers::MatrixBase
Definition: petsc_matrix_base.h:284
PETScWrappers::MatrixBase::get_mpi_communicator
virtual const MPI_Comm & get_mpi_communicator() const =0
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
PETScWrappers::MatrixBase::clear_rows
void clear_rows(const std::vector< size_type > &rows, const PetscScalar new_diag_value=0)
Definition: petsc_matrix_base.cc:136
PETScWrappers::MatrixBase::add
void add(const size_type i, const size_type j, const PetscScalar value)
DeclException3
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:564
PETScWrappers::MatrixBase::operator/=
MatrixBase & operator/=(const PetscScalar factor)
Definition: petsc_matrix_base.cc:419
PETScWrappers::MatrixBase::is_hermitian
PetscBool is_hermitian(const double tolerance=1.e-12)
Definition: petsc_matrix_base.cc:630
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
DEAL_II_DEPRECATED
#define DEAL_II_DEPRECATED
Definition: config.h:98
VectorOperation::values
values
Definition: vector_operation.h:40
PETScWrappers::MatrixBase::compress
void compress(const VectorOperation::values operation)
Definition: petsc_matrix_base.cc:193
PETScWrappers::MatrixBase::assert_is_compressed
void assert_is_compressed()
PETScWrappers::MatrixIterators::const_iterator::accessor
Accessor accessor
Definition: petsc_matrix_base.h:246
VectorOperation::unknown
@ unknown
Definition: vector_operation.h:45
PETScWrappers::MatrixIterators::const_iterator
Definition: petsc_matrix_base.h:68
PETScWrappers::MatrixBase::matrix
Mat matrix
Definition: petsc_matrix_base.h:968
unsigned int
value
static const bool value
Definition: dof_tools_constraints.cc:433
PETScWrappers::MatrixBase::MatrixBase
MatrixBase()
Definition: petsc_matrix_base.cc:77
LinearAlgebra::CUDAWrappers::kernel::size_type
types::global_dof_index size_type
Definition: cuda_kernels.h:45
PETScWrappers::MatrixBase::transpose
void transpose()
Definition: petsc_matrix_base.cc:608
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
PETScWrappers::MatrixBase::clear_row
void clear_row(const size_type row, const PetscScalar new_diag_value=0)
Definition: petsc_matrix_base.cc:127
PETScWrappers::MatrixIterators::const_iterator::Accessor::a_row
size_type a_row
Definition: petsc_matrix_base.h:139
PETScWrappers::MatrixBase::value_type
PetscScalar value_type
Definition: petsc_matrix_base.h:300
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
PETScWrappers::MatrixBase::residual
PetscScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
Definition: petsc_matrix_base.cc:580
LACExceptions::ExcPETScError
Definition: exceptions.h:56
PETScWrappers::MatrixBase::Tvmult_add
void Tvmult_add(VectorBase &dst, const VectorBase &src) const
Definition: petsc_matrix_base.cc:482
PETScWrappers::MatrixBase::end
const_iterator end() const
PETScWrappers::MatrixBase::ExcSourceEqualsDestination
static ::ExceptionBase & ExcSourceEqualsDestination()
PETScWrappers::MatrixBase::begin
const_iterator begin() const
PETScWrappers::MatrixBase::m
size_type m() const
Definition: petsc_matrix_base.cc:237
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
PETScWrappers::MatrixIterators::const_iterator::operator<
bool operator<(const const_iterator &) const
PETScWrappers::MatrixBase::operator=
MatrixBase & operator=(const MatrixBase &)=delete
PETScWrappers::MatrixBase::el
PetscScalar el(const size_type i, const size_type j) const
Definition: petsc_matrix_base.cc:165
config.h
FullMatrix
Definition: full_matrix.h:71
PETScWrappers::MatrixBase::prepare_set
void prepare_set()
PETScWrappers::MatrixBase::ExcWrongMode
static ::ExceptionBase & ExcWrongMode(int arg1, int arg2)
PETScWrappers::MatrixBase::~MatrixBase
virtual ~MatrixBase() override
Definition: petsc_matrix_base.cc:84
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
first
Point< 2 > first
Definition: grid_out.cc:4352
PETScWrappers::MatrixIterators::const_iterator::Accessor::a_index
size_type a_index
Definition: petsc_matrix_base.h:144
PETScWrappers::MatrixIterators::const_iterator::Accessor::ExcAccessToNonlocalRow
static ::ExceptionBase & ExcAccessToNonlocalRow(int arg1, int arg2, int arg3)
PETScWrappers::VectorBase
Definition: petsc_vector_base.h:240
PETScWrappers::MatrixBase::last_action
VectorOperation::values last_action
Definition: petsc_matrix_base.h:973
PETScWrappers::MatrixIterators::const_iterator::Accessor::row
size_type row() const
PETScWrappers::MatrixBase::linfty_norm
PetscReal linfty_norm() const
Definition: petsc_matrix_base.cc:351
DeclException2
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:541
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
full_matrix.h
PETScWrappers::MatrixIterators::const_iterator::Accessor::matrix
MatrixBase * matrix
Definition: petsc_matrix_base.h:134
StandardExceptions::ExcIteratorPastEnd
static ::ExceptionBase & ExcIteratorPastEnd()
operator*
std::enable_if< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typename ProductType< std::complex< T >, std::complex< U > >::type >::type operator*(const std::complex< T > &left, const std::complex< U > &right)
Definition: complex_overloads.h:43
PETScWrappers::MatrixBase::operator()
PetscScalar operator()(const size_type i, const size_type j) const
PETScWrappers::MatrixIterators::const_iterator::operator->
const Accessor * operator->() const
TrilinosWrappers::SparseMatrixIterators::ExcBeyondEndOfMatrix
static ::ExceptionBase & ExcBeyondEndOfMatrix()