Reference documentation for deal.II version 8.5.1
petsc_matrix_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__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 
24 # include <deal.II/base/subscriptor.h>
25 # include <deal.II/lac/exceptions.h>
26 # include <deal.II/lac/full_matrix.h>
27 # include <deal.II/lac/petsc_compatibility.h>
28 # include <deal.II/lac/vector.h>
29 
30 # include <petscmat.h>
31 # include <deal.II/base/std_cxx11/shared_ptr.h>
32 
33 # include <vector>
34 # include <cmath>
35 
36 DEAL_II_NAMESPACE_OPEN
37 
38 template <typename Matrix> class BlockMatrixBase;
39 
40 
41 namespace PETScWrappers
42 {
43  // forward declarations
44  class VectorBase;
45  class MatrixBase;
46 
47  namespace MatrixIterators
48  {
64  {
65  private:
69  class Accessor
70  {
71  public:
76 
81  Accessor (const MatrixBase *matrix,
82  const size_type row,
83  const size_type index);
84 
88  size_type row() const;
89 
93  size_type index() const;
94 
98  size_type column() const;
99 
103  PetscScalar value() const;
104 
113  int, int, int,
114  << "You tried to access row " << arg1
115  << " of a distributed matrix, but only rows "
116  << arg2 << " through " << arg3
117  << " are stored locally and can be accessed.");
118 
119  private:
123  mutable MatrixBase *matrix;
124 
129 
134 
147  std_cxx11::shared_ptr<const std::vector<size_type> > colnum_cache;
148 
152  std_cxx11::shared_ptr<const std::vector<PetscScalar> > value_cache;
153 
159  void visit_present_row ();
160 
164  friend class const_iterator;
165  };
166 
167  public:
172 
177  const_iterator (const MatrixBase *matrix,
178  const size_type row,
179  const size_type index);
180 
185 
190 
194  const Accessor &operator* () const;
195 
199  const Accessor *operator-> () const;
200 
205  bool operator == (const const_iterator &) const;
209  bool operator != (const const_iterator &) const;
210 
216  bool operator < (const const_iterator &) const;
217 
222  int, int,
223  << "Attempt to access element " << arg2
224  << " of row " << arg1
225  << " which doesn't have that many elements.");
226 
227  private:
232  };
233 
234  }
235 
236 
269  class MatrixBase : public Subscriptor
270  {
271  public:
276 
281 
285  typedef PetscScalar value_type;
286 
290  MatrixBase ();
291 
295  virtual ~MatrixBase ();
296 
306  MatrixBase &
307  operator = (const value_type d);
312  void clear ();
313 
323  void set (const size_type i,
324  const size_type j,
325  const PetscScalar value);
326 
346  void set (const std::vector<size_type> &indices,
347  const FullMatrix<PetscScalar> &full_matrix,
348  const bool elide_zero_values = false);
349 
355  void set (const std::vector<size_type> &row_indices,
356  const std::vector<size_type> &col_indices,
357  const FullMatrix<PetscScalar> &full_matrix,
358  const bool elide_zero_values = false);
359 
374  void set (const size_type row,
375  const std::vector<size_type > &col_indices,
376  const std::vector<PetscScalar> &values,
377  const bool elide_zero_values = false);
378 
393  void set (const size_type row,
394  const size_type n_cols,
395  const size_type *col_indices,
396  const PetscScalar *values,
397  const bool elide_zero_values = false);
398 
408  void add (const size_type i,
409  const size_type j,
410  const PetscScalar value);
411 
431  void add (const std::vector<size_type> &indices,
432  const FullMatrix<PetscScalar> &full_matrix,
433  const bool elide_zero_values = true);
434 
440  void add (const std::vector<size_type> &row_indices,
441  const std::vector<size_type> &col_indices,
442  const FullMatrix<PetscScalar> &full_matrix,
443  const bool elide_zero_values = true);
444 
459  void add (const size_type row,
460  const std::vector<size_type> &col_indices,
461  const std::vector<PetscScalar> &values,
462  const bool elide_zero_values = true);
463 
478  void add (const size_type row,
479  const size_type n_cols,
480  const size_type *col_indices,
481  const PetscScalar *values,
482  const bool elide_zero_values = true,
483  const bool col_indices_are_sorted = false);
484 
501  void clear_row (const size_type row,
502  const PetscScalar new_diag_value = 0);
503 
512  void clear_rows (const std::vector<size_type> &rows,
513  const PetscScalar new_diag_value = 0);
514 
526  void compress (const VectorOperation::values operation);
527 
539  PetscScalar operator () (const size_type i,
540  const size_type j) const;
541 
549  PetscScalar el (const size_type i,
550  const size_type j) const;
551 
561  PetscScalar diag_element (const size_type i) const;
562 
566  size_type m () const;
567 
571  size_type n () const;
572 
581  size_type local_size () const;
582 
591  std::pair<size_type, size_type>
592  local_range () const;
593 
598  bool in_local_range (const size_type index) const;
599 
604  virtual const MPI_Comm &get_mpi_communicator () const = 0;
605 
611  size_type n_nonzero_elements () const;
612 
616  size_type row_length (const size_type row) const;
617 
625  PetscReal l1_norm () const;
626 
634  PetscReal linfty_norm () const;
635 
640  PetscReal frobenius_norm () const;
641 
642 
662  PetscScalar matrix_norm_square (const VectorBase &v) const;
663 
664 
678  PetscScalar matrix_scalar_product (const VectorBase &u,
679  const VectorBase &v) const;
680 
681 
682 #if DEAL_II_PETSC_VERSION_GTE(3,1,0)
683 
687  PetscScalar trace () const;
688 #endif
689 
693  MatrixBase &operator *= (const PetscScalar factor);
694 
698  MatrixBase &operator /= (const PetscScalar factor);
699 
700 
705  MatrixBase &add (const PetscScalar factor,
706  const MatrixBase &other);
707 
708 
714  MatrixBase &add (const MatrixBase &other,
715  const PetscScalar factor) DEAL_II_DEPRECATED;
716 
728  void vmult (VectorBase &dst,
729  const VectorBase &src) const;
730 
743  void Tvmult (VectorBase &dst,
744  const VectorBase &src) const;
745 
757  void vmult_add (VectorBase &dst,
758  const VectorBase &src) const;
759 
772  void Tvmult_add (VectorBase &dst,
773  const VectorBase &src) const;
774 
775 
788  PetscScalar residual (VectorBase &dst,
789  const VectorBase &x,
790  const VectorBase &b) const;
791 
795  const_iterator begin () const;
796 
800  const_iterator end () const;
801 
810  const_iterator begin (const size_type r) const;
811 
820  const_iterator end (const size_type r) const;
821 
829  operator Mat () const;
830 
834  void transpose ();
835 
840  PetscBooleanType
841  is_symmetric (const double tolerance = 1.e-12);
842 
848  PetscBooleanType
849  is_hermitian (const double tolerance = 1.e-12);
850 
857  void write_ascii (const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
858 
866  void print (std::ostream &out,
867  const bool alternative_output = false) const;
868 
872  std::size_t memory_consumption() const;
873 
878 
883  int, int,
884  << "You tried to do a "
885  << (arg1 == 1 ?
886  "'set'" :
887  (arg1 == 2 ?
888  "'add'" : "???"))
889  << " operation but the matrix is currently in "
890  << (arg2 == 1 ?
891  "'set'" :
892  (arg2 == 2 ?
893  "'add'" : "???"))
894  << " mode. You first have to call 'compress()'.");
895 
896  protected:
901  Mat matrix;
902 
907 
913  void prepare_action(const VectorOperation::values new_action);
914 
920  void assert_is_compressed();
921 
932  void prepare_add();
938  void prepare_set();
939 
940 
941 
942  private:
943 
947  MatrixBase(const MatrixBase &);
951  MatrixBase &operator=(const MatrixBase &);
952 
958  std::vector<PetscInt> column_indices;
959 
965  std::vector<PetscScalar> column_values;
966 
967 
971  template <class> friend class ::BlockMatrixBase;
972  };
973 
974 
975 
976 #ifndef DOXYGEN
977 // -------------------------- inline and template functions ----------------------
978 
979 
980  namespace MatrixIterators
981  {
982 
983  inline
985  Accessor (const MatrixBase *matrix,
986  const size_type row,
987  const size_type index)
988  :
989  matrix(const_cast<MatrixBase *>(matrix)),
990  a_row(row),
991  a_index(index)
992  {
993  visit_present_row ();
994  }
995 
996 
997 
998  inline
1001  {
1002  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1003  return a_row;
1004  }
1005 
1006 
1007  inline
1010  {
1011  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1012  return (*colnum_cache)[a_index];
1013  }
1014 
1015 
1016  inline
1019  {
1020  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1021  return a_index;
1022  }
1023 
1024 
1025  inline
1026  PetscScalar
1028  {
1029  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1030  return (*value_cache)[a_index];
1031  }
1032 
1033 
1034  inline
1036  const_iterator(const MatrixBase *matrix,
1037  const size_type row,
1038  const size_type index)
1039  :
1040  accessor(matrix, row, index)
1041  {}
1042 
1043 
1044 
1045  inline
1046  const_iterator &
1047  const_iterator::operator++ ()
1048  {
1049  Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
1050 
1051  ++accessor.a_index;
1052 
1053  // if at end of line: do one step, then
1054  // cycle until we find a row with a
1055  // nonzero number of entries
1056  if (accessor.a_index >= accessor.colnum_cache->size())
1057  {
1058  accessor.a_index = 0;
1059  ++accessor.a_row;
1060 
1061  while ((accessor.a_row < accessor.matrix->m())
1062  &&
1063  (accessor.matrix->row_length(accessor.a_row) == 0))
1064  ++accessor.a_row;
1065 
1066  accessor.visit_present_row();
1067  }
1068  return *this;
1069  }
1070 
1071 
1072  inline
1073  const_iterator
1074  const_iterator::operator++ (int)
1075  {
1076  const const_iterator old_state = *this;
1077  ++(*this);
1078  return old_state;
1079  }
1080 
1081 
1082  inline
1083  const const_iterator::Accessor &
1084  const_iterator::operator* () const
1085  {
1086  return accessor;
1087  }
1088 
1089 
1090  inline
1091  const const_iterator::Accessor *
1092  const_iterator::operator-> () const
1093  {
1094  return &accessor;
1095  }
1096 
1097 
1098  inline
1099  bool
1100  const_iterator::
1101  operator == (const const_iterator &other) const
1102  {
1103  return (accessor.a_row == other.accessor.a_row &&
1104  accessor.a_index == other.accessor.a_index);
1105  }
1106 
1107 
1108  inline
1109  bool
1110  const_iterator::
1111  operator != (const const_iterator &other) const
1112  {
1113  return ! (*this == other);
1114  }
1115 
1116 
1117  inline
1118  bool
1119  const_iterator::
1120  operator < (const const_iterator &other) const
1121  {
1122  return (accessor.row() < other.accessor.row() ||
1123  (accessor.row() == other.accessor.row() &&
1124  accessor.index() < other.accessor.index()));
1125  }
1126 
1127  }
1128 
1129 
1130 
1131  // Inline the set() and add()
1132  // functions, since they will be
1133  // called frequently, and the
1134  // compiler can optimize away
1135  // some unnecessary loops when
1136  // the sizes are given at
1137  // compile time.
1138  inline
1139  void
1140  MatrixBase::set (const size_type i,
1141  const size_type j,
1142  const PetscScalar value)
1143  {
1144  AssertIsFinite(value);
1145 
1146  set (i, 1, &j, &value, false);
1147  }
1148 
1149 
1150 
1151  inline
1152  void
1153  MatrixBase::set (const std::vector<size_type> &indices,
1154  const FullMatrix<PetscScalar> &values,
1155  const bool elide_zero_values)
1156  {
1157  Assert (indices.size() == values.m(),
1158  ExcDimensionMismatch(indices.size(), values.m()));
1159  Assert (values.m() == values.n(), ExcNotQuadratic());
1160 
1161  for (size_type i=0; i<indices.size(); ++i)
1162  set (indices[i], indices.size(), &indices[0], &values(i,0),
1163  elide_zero_values);
1164  }
1165 
1166 
1167 
1168  inline
1169  void
1170  MatrixBase::set (const std::vector<size_type> &row_indices,
1171  const std::vector<size_type> &col_indices,
1172  const FullMatrix<PetscScalar> &values,
1173  const bool elide_zero_values)
1174  {
1175  Assert (row_indices.size() == values.m(),
1176  ExcDimensionMismatch(row_indices.size(), values.m()));
1177  Assert (col_indices.size() == values.n(),
1178  ExcDimensionMismatch(col_indices.size(), values.n()));
1179 
1180  for (size_type i=0; i<row_indices.size(); ++i)
1181  set (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1182  elide_zero_values);
1183  }
1184 
1185 
1186 
1187  inline
1188  void
1189  MatrixBase::set (const size_type row,
1190  const std::vector<size_type> &col_indices,
1191  const std::vector<PetscScalar> &values,
1192  const bool elide_zero_values)
1193  {
1194  Assert (col_indices.size() == values.size(),
1195  ExcDimensionMismatch(col_indices.size(), values.size()));
1196 
1197  set (row, col_indices.size(), &col_indices[0], &values[0],
1198  elide_zero_values);
1199  }
1200 
1201 
1202 
1203  inline
1204  void
1205  MatrixBase::set (const size_type row,
1206  const size_type n_cols,
1207  const size_type *col_indices,
1208  const PetscScalar *values,
1209  const bool elide_zero_values)
1210  {
1211  (void)elide_zero_values;
1212 
1213  prepare_action(VectorOperation::insert);
1214 
1215  const PetscInt petsc_i = row;
1216  PetscInt *col_index_ptr;
1217 
1218  PetscScalar const *col_value_ptr;
1219  int n_columns;
1220 
1221  // If we don't elide zeros, the pointers are already available...
1222 #ifndef PETSC_USE_64BIT_INDICES
1223  if (elide_zero_values == false)
1224  {
1225  col_index_ptr = (int *)col_indices;
1226  col_value_ptr = values;
1227  n_columns = n_cols;
1228  }
1229  else
1230 #endif
1231  {
1232  // Otherwise, extract nonzero values in each row and get the
1233  // respective index.
1234  if (column_indices.size() < n_cols)
1235  {
1236  column_indices.resize(n_cols);
1237  column_values.resize(n_cols);
1238  }
1239 
1240  n_columns = 0;
1241  for (size_type j=0; j<n_cols; ++j)
1242  {
1243  const PetscScalar value = values[j];
1244  AssertIsFinite(value);
1245  if (value != PetscScalar())
1246  {
1247  column_indices[n_columns] = col_indices[j];
1248  column_values[n_columns] = value;
1249  n_columns++;
1250  }
1251  }
1252  Assert(n_columns <= (int)n_cols, ExcInternalError());
1253 
1254  col_index_ptr = &column_indices[0];
1255  col_value_ptr = &column_values[0];
1256  }
1257 
1258  const PetscErrorCode ierr = MatSetValues (matrix, 1, &petsc_i, n_columns,
1259  col_index_ptr,
1260  col_value_ptr, INSERT_VALUES);
1261  AssertThrow (ierr == 0, ExcPETScError(ierr));
1262  }
1263 
1264 
1265 
1266  inline
1267  void
1268  MatrixBase::add (const size_type i,
1269  const size_type j,
1270  const PetscScalar value)
1271  {
1272 
1273  AssertIsFinite(value);
1274 
1275  if (value == PetscScalar())
1276  {
1277  // we have to check after using Insert/Add in any case to be
1278  // consistent with the MPI communication model (see the comments in
1279  // the documentation of TrilinosWrappers::Vector), but we can save
1280  // some work if the addend is zero. However, these actions are done
1281  // in case we pass on to the other function.
1282  prepare_action(VectorOperation::add);
1283 
1284  return;
1285  }
1286  else
1287  add (i, 1, &j, &value, false);
1288  }
1289 
1290 
1291 
1292  inline
1293  void
1294  MatrixBase::add (const std::vector<size_type> &indices,
1295  const FullMatrix<PetscScalar> &values,
1296  const bool elide_zero_values)
1297  {
1298  Assert (indices.size() == values.m(),
1299  ExcDimensionMismatch(indices.size(), values.m()));
1300  Assert (values.m() == values.n(), ExcNotQuadratic());
1301 
1302  for (size_type i=0; i<indices.size(); ++i)
1303  add (indices[i], indices.size(), &indices[0], &values(i,0),
1304  elide_zero_values);
1305  }
1306 
1307 
1308 
1309  inline
1310  void
1311  MatrixBase::add (const std::vector<size_type> &row_indices,
1312  const std::vector<size_type> &col_indices,
1313  const FullMatrix<PetscScalar> &values,
1314  const bool elide_zero_values)
1315  {
1316  Assert (row_indices.size() == values.m(),
1317  ExcDimensionMismatch(row_indices.size(), values.m()));
1318  Assert (col_indices.size() == values.n(),
1319  ExcDimensionMismatch(col_indices.size(), values.n()));
1320 
1321  for (size_type i=0; i<row_indices.size(); ++i)
1322  add (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1323  elide_zero_values);
1324  }
1325 
1326 
1327 
1328  inline
1329  void
1330  MatrixBase::add (const size_type row,
1331  const std::vector<size_type> &col_indices,
1332  const std::vector<PetscScalar> &values,
1333  const bool elide_zero_values)
1334  {
1335  Assert (col_indices.size() == values.size(),
1336  ExcDimensionMismatch(col_indices.size(), values.size()));
1337 
1338  add (row, col_indices.size(), &col_indices[0], &values[0],
1339  elide_zero_values);
1340  }
1341 
1342 
1343 
1344  inline
1345  void
1346  MatrixBase::add (const size_type row,
1347  const size_type n_cols,
1348  const size_type *col_indices,
1349  const PetscScalar *values,
1350  const bool elide_zero_values,
1351  const bool /*col_indices_are_sorted*/)
1352  {
1353  (void)elide_zero_values;
1354 
1355  prepare_action(VectorOperation::add);
1356 
1357  const PetscInt petsc_i = row;
1358  PetscInt *col_index_ptr;
1359 
1360  PetscScalar const *col_value_ptr;
1361  int n_columns;
1362 
1363  // If we don't elide zeros, the pointers are already available...
1364 #ifndef PETSC_USE_64BIT_INDICES
1365  if (elide_zero_values == false)
1366  {
1367  col_index_ptr = (int *)col_indices;
1368  col_value_ptr = values;
1369  n_columns = n_cols;
1370  }
1371  else
1372 #endif
1373  {
1374  // Otherwise, extract nonzero values in each row and get the
1375  // respective index.
1376  if (column_indices.size() < n_cols)
1377  {
1378  column_indices.resize(n_cols);
1379  column_values.resize(n_cols);
1380  }
1381 
1382  n_columns = 0;
1383  for (size_type j=0; j<n_cols; ++j)
1384  {
1385  const PetscScalar value = values[j];
1386  AssertIsFinite(value);
1387  if (value != PetscScalar())
1388  {
1389  column_indices[n_columns] = col_indices[j];
1390  column_values[n_columns] = value;
1391  n_columns++;
1392  }
1393  }
1394  Assert(n_columns <= (int)n_cols, ExcInternalError());
1395 
1396  col_index_ptr = &column_indices[0];
1397  col_value_ptr = &column_values[0];
1398  }
1399 
1400  const PetscErrorCode ierr = MatSetValues (matrix, 1, &petsc_i, n_columns,
1401  col_index_ptr, col_value_ptr,
1402  ADD_VALUES);
1403  AssertThrow (ierr == 0, ExcPETScError(ierr));
1404  }
1405 
1406 
1407 
1408 
1409 
1410 
1411  inline
1412  PetscScalar
1413  MatrixBase::operator() (const size_type i,
1414  const size_type j) const
1415  {
1416  return el(i,j);
1417  }
1418 
1419 
1420 
1421  inline
1422  MatrixBase::const_iterator
1423  MatrixBase::begin() const
1424  {
1425  return const_iterator(this, 0, 0);
1426  }
1427 
1428 
1429  inline
1430  MatrixBase::const_iterator
1431  MatrixBase::end() const
1432  {
1433  return const_iterator(this, m(), 0);
1434  }
1435 
1436 
1437  inline
1438  MatrixBase::const_iterator
1439  MatrixBase::begin(const size_type r) const
1440  {
1441  Assert (r < m(), ExcIndexRange(r, 0, m()));
1442  if (row_length(r) > 0)
1443  return const_iterator(this, r, 0);
1444  else
1445  return end (r);
1446  }
1447 
1448 
1449  inline
1450  MatrixBase::const_iterator
1451  MatrixBase::end(const size_type r) const
1452  {
1453  Assert (r < m(), ExcIndexRange(r, 0, m()));
1454 
1455  // place the iterator on the first entry
1456  // past this line, or at the end of the
1457  // matrix
1458  for (size_type i=r+1; i<m(); ++i)
1459  if (row_length(i) > 0)
1460  return const_iterator(this, i, 0);
1461 
1462  // if there is no such line, then take the
1463  // end iterator of the matrix
1464  return end();
1465  }
1466 
1467 
1468 
1469  inline
1470  bool
1471  MatrixBase::in_local_range (const size_type index) const
1472  {
1473  PetscInt begin, end;
1474 
1475  const PetscErrorCode ierr = MatGetOwnershipRange (static_cast<const Mat &>(matrix),
1476  &begin, &end);
1477  AssertThrow (ierr == 0, ExcPETScError(ierr));
1478 
1479  return ((index >= static_cast<size_type>(begin)) &&
1480  (index < static_cast<size_type>(end)));
1481  }
1482 
1483 
1484 
1485  inline
1486  void
1487  MatrixBase::prepare_action(const VectorOperation::values new_action)
1488  {
1489  if (last_action == VectorOperation::unknown)
1490  last_action = new_action;
1491 
1492  Assert (last_action == new_action, ExcWrongMode (last_action, new_action));
1493  }
1494 
1495 
1496 
1497  inline
1498  void
1499  MatrixBase::assert_is_compressed ()
1500  {
1501  // compress() sets the last action to none, which allows us to check if there
1502  // are pending add/insert operations:
1503  AssertThrow (last_action == VectorOperation::unknown,
1504  ExcMessage("Error: missing compress() call."));
1505  }
1506 
1507 
1508 
1509  inline
1510  void
1511  MatrixBase::prepare_add()
1512  {
1513  prepare_action(VectorOperation::add);
1514  }
1515 
1516 
1517 
1518  inline
1519  void
1520  MatrixBase::prepare_set()
1521  {
1522  prepare_action(VectorOperation::insert);
1523  }
1524 
1525 #endif // DOXYGEN
1526 }
1527 
1528 
1529 DEAL_II_NAMESPACE_CLOSE
1530 
1531 
1532 #endif // DEAL_II_WITH_PETSC
1533 
1534 
1535 /*---------------------------- petsc_matrix_base.h ---------------------------*/
1536 
1537 #endif
1538 /*---------------------------- petsc_matrix_base.h ---------------------------*/
size_type m() const
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
bool operator==(const const_iterator &) const
PetscScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
static ::ExceptionBase & ExcBeyondEndOfMatrix()
void add(const size_type i, const size_type j, const PetscScalar value)
void vmult(VectorBase &dst, const VectorBase &src) const
const_iterator begin() const
PetscReal frobenius_norm() const
static ::ExceptionBase & ExcInvalidIndexWithinRow(int arg1, int arg2)
std_cxx11::shared_ptr< const std::vector< PetscScalar > > value_cache
size_type row_length(const size_type row) const
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void clear_rows(const std::vector< size_type > &rows, const PetscScalar new_diag_value=0)
std::vector< PetscInt > column_indices
const_iterator(const MatrixBase *matrix, const size_type row, const size_type index)
void compress(const VectorOperation::values operation)
std::vector< PetscScalar > column_values
MatrixIterators::const_iterator const_iterator
size_type n_nonzero_elements() const
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type n() const
unsigned int global_dof_index
Definition: types.h:88
Accessor(const MatrixBase *matrix, const size_type row, const size_type index)
bool operator<(const const_iterator &) const
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
MatrixBase & operator=(const value_type d)
#define DeclException0(Exception0)
Definition: exceptions.h:541
void vmult_add(VectorBase &dst, const VectorBase &src) const
PetscScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
types::global_dof_index size_type
PetscScalar diag_element(const size_type i) const
void Tvmult(VectorBase &dst, const VectorBase &src) const
PetscScalar el(const size_type i, const size_type j) const
void Tvmult_add(VectorBase &dst, const VectorBase &src) const
std_cxx11::shared_ptr< const std::vector< size_type > > colnum_cache
const_iterator end() const
PetscScalar operator()(const size_type i, const size_type j) const
static ::ExceptionBase & ExcIteratorPastEnd()
void print(std::ostream &out, const bool alternative_output=false) const
static ::ExceptionBase & ExcNotQuadratic()
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
VectorOperation::values last_action
void prepare_action(const VectorOperation::values new_action)
PetscBooleanType is_symmetric(const double tolerance=1.e-12)
PetscBooleanType is_hermitian(const double tolerance=1.e-12)
bool in_local_range(const size_type index) const
MatrixBase & operator*=(const PetscScalar factor)
std::pair< size_type, size_type > local_range() const
static ::ExceptionBase & ExcAccessToNonlocalRow(int arg1, int arg2, int arg3)
PetscScalar matrix_norm_square(const VectorBase &v) const
MatrixBase & operator/=(const PetscScalar factor)
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcWrongMode(int arg1, int arg2)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:588
bool operator!=(const const_iterator &) const
virtual const MPI_Comm & get_mpi_communicator() const =0
#define AssertIsFinite(number)
Definition: exceptions.h:1197
void clear_row(const size_type row, const PetscScalar new_diag_value=0)
std::size_t memory_consumption() const
static ::ExceptionBase & ExcInternalError()