Reference documentation for deal.II version 8.4.2
petsc_matrix_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2016 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/full_matrix.h>
26 # include <deal.II/lac/exceptions.h>
27 # include <deal.II/lac/vector.h>
28 
29 # include <petscmat.h>
30 # include <deal.II/base/std_cxx11/shared_ptr.h>
31 
32 # include <vector>
33 # include <cmath>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 template <typename Matrix> class BlockMatrixBase;
38 
39 
40 namespace PETScWrappers
41 {
42  // forward declarations
43  class VectorBase;
44  class MatrixBase;
45 
46  namespace MatrixIterators
47  {
63  {
64  private:
68  class Accessor
69  {
70  public:
75 
80  Accessor (const MatrixBase *matrix,
81  const size_type row,
82  const size_type index);
83 
87  Accessor (const Accessor &a);
88 
92  size_type row() const;
93 
97  size_type index() const;
98 
102  size_type column() const;
103 
107  PetscScalar value() const;
108 
112  DeclException0 (ExcBeyondEndOfMatrix);
116  DeclException3 (ExcAccessToNonlocalRow,
117  int, int, int,
118  << "You tried to access row " << arg1
119  << " of a distributed matrix, but only rows "
120  << arg2 << " through " << arg3
121  << " are stored locally and can be accessed.");
122 
123  private:
127  mutable MatrixBase *matrix;
128 
132  size_type a_row;
133 
137  size_type a_index;
138 
151  std_cxx11::shared_ptr<const std::vector<size_type> > colnum_cache;
152 
156  std_cxx11::shared_ptr<const std::vector<PetscScalar> > value_cache;
157 
163  void visit_present_row ();
164 
168  friend class const_iterator;
169  };
170 
171  public:
176 
182  const size_type row,
183  const size_type index);
184 
189 
194 
198  const Accessor &operator* () const;
199 
203  const Accessor *operator-> () const;
204 
209  bool operator == (const const_iterator &) const;
213  bool operator != (const const_iterator &) const;
214 
220  bool operator < (const const_iterator &) const;
221 
225  DeclException2 (ExcInvalidIndexWithinRow,
226  int, int,
227  << "Attempt to access element " << arg2
228  << " of row " << arg1
229  << " which doesn't have that many elements.");
230 
231  private:
236  };
237 
238  }
239 
240 
273  class MatrixBase : public Subscriptor
274  {
275  public:
280 
285 
289  typedef PetscScalar value_type;
290 
294  MatrixBase ();
295 
299  virtual ~MatrixBase ();
300 
310  MatrixBase &
311  operator = (const value_type d);
316  void clear ();
317 
327  void set (const size_type i,
328  const size_type j,
329  const PetscScalar value);
330 
350  void set (const std::vector<size_type> &indices,
351  const FullMatrix<PetscScalar> &full_matrix,
352  const bool elide_zero_values = false);
353 
359  void set (const std::vector<size_type> &row_indices,
360  const std::vector<size_type> &col_indices,
361  const FullMatrix<PetscScalar> &full_matrix,
362  const bool elide_zero_values = false);
363 
378  void set (const size_type row,
379  const std::vector<size_type > &col_indices,
380  const std::vector<PetscScalar> &values,
381  const bool elide_zero_values = false);
382 
397  void set (const size_type row,
398  const size_type n_cols,
399  const size_type *col_indices,
400  const PetscScalar *values,
401  const bool elide_zero_values = false);
402 
412  void add (const size_type i,
413  const size_type j,
414  const PetscScalar value);
415 
435  void add (const std::vector<size_type> &indices,
436  const FullMatrix<PetscScalar> &full_matrix,
437  const bool elide_zero_values = true);
438 
444  void add (const std::vector<size_type> &row_indices,
445  const std::vector<size_type> &col_indices,
446  const FullMatrix<PetscScalar> &full_matrix,
447  const bool elide_zero_values = true);
448 
463  void add (const size_type row,
464  const std::vector<size_type> &col_indices,
465  const std::vector<PetscScalar> &values,
466  const bool elide_zero_values = true);
467 
482  void add (const size_type row,
483  const size_type n_cols,
484  const size_type *col_indices,
485  const PetscScalar *values,
486  const bool elide_zero_values = true,
487  const bool col_indices_are_sorted = false);
488 
505  void clear_row (const size_type row,
506  const PetscScalar new_diag_value = 0);
507 
516  void clear_rows (const std::vector<size_type> &rows,
517  const PetscScalar new_diag_value = 0);
518 
530  void compress (const VectorOperation::values operation);
531 
543  PetscScalar operator () (const size_type i,
544  const size_type j) const;
545 
553  PetscScalar el (const size_type i,
554  const size_type j) const;
555 
565  PetscScalar diag_element (const size_type i) const;
566 
570  size_type m () const;
571 
575  size_type n () const;
576 
585  size_type local_size () const;
586 
595  std::pair<size_type, size_type>
596  local_range () const;
597 
602  bool in_local_range (const size_type index) const;
603 
608  virtual const MPI_Comm &get_mpi_communicator () const = 0;
609 
615  size_type n_nonzero_elements () const;
616 
620  size_type row_length (const size_type row) const;
621 
629  PetscReal l1_norm () const;
630 
638  PetscReal linfty_norm () const;
639 
644  PetscReal frobenius_norm () const;
645 
646 
666  PetscScalar matrix_norm_square (const VectorBase &v) const;
667 
668 
682  PetscScalar matrix_scalar_product (const VectorBase &u,
683  const VectorBase &v) const;
684 
685 
686 #if DEAL_II_PETSC_VERSION_GTE(3,1,0)
687 
691  PetscScalar trace () const;
692 #endif
693 
697  MatrixBase &operator *= (const PetscScalar factor);
698 
702  MatrixBase &operator /= (const PetscScalar factor);
703 
708  MatrixBase &add (const MatrixBase &other,
709  const PetscScalar factor);
710 
722  void vmult (VectorBase &dst,
723  const VectorBase &src) const;
724 
737  void Tvmult (VectorBase &dst,
738  const VectorBase &src) const;
739 
751  void vmult_add (VectorBase &dst,
752  const VectorBase &src) const;
753 
766  void Tvmult_add (VectorBase &dst,
767  const VectorBase &src) const;
768 
769 
782  PetscScalar residual (VectorBase &dst,
783  const VectorBase &x,
784  const VectorBase &b) const;
785 
789  const_iterator begin () const;
790 
794  const_iterator end () const;
795 
804  const_iterator begin (const size_type r) const;
805 
814  const_iterator end (const size_type r) const;
815 
823  operator Mat () const;
824 
828  void transpose ();
829 
834 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
835  PetscTruth
836 #else
837  PetscBool
838 #endif
839  is_symmetric (const double tolerance = 1.e-12);
840 
846 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
847  PetscTruth
848 #else
849  PetscBool
850 #endif
851  is_hermitian (const double tolerance = 1.e-12);
852 
860  void write_ascii (const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
861 
869  void print (std::ostream &out,
870  const bool alternative_output = false) const;
871 
875  std::size_t memory_consumption() const;
876 
880  DeclException1 (ExcPETScError,
881  int,
882  << "An error with error number " << arg1
883  << " occurred while calling a PETSc function");
887  DeclException0 (ExcSourceEqualsDestination);
888 
892  DeclException2 (ExcWrongMode,
893  int, int,
894  << "You tried to do a "
895  << (arg1 == 1 ?
896  "'set'" :
897  (arg1 == 2 ?
898  "'add'" : "???"))
899  << " operation but the matrix is currently in "
900  << (arg2 == 1 ?
901  "'set'" :
902  (arg2 == 2 ?
903  "'add'" : "???"))
904  << " mode. You first have to call 'compress()'.");
905 
906  protected:
911  Mat matrix;
912 
916  VectorOperation::values last_action;
917 
923  void prepare_action(const VectorOperation::values new_action);
924 
930  void assert_is_compressed();
931 
942  void prepare_add();
948  void prepare_set();
949 
950 
951 
952  private:
953 
957  MatrixBase(const MatrixBase &);
961  MatrixBase &operator=(const MatrixBase &);
962 
968  std::vector<PetscInt> column_indices;
969 
975  std::vector<PetscScalar> column_values;
976 
977 
981  template <class> friend class ::BlockMatrixBase;
982  };
983 
984 
985 
986 #ifndef DOXYGEN
987 // -------------------------- inline and template functions ----------------------
988 
989 
990  namespace MatrixIterators
991  {
992 
993  inline
995  Accessor (const MatrixBase *matrix,
996  const size_type row,
997  const size_type index)
998  :
999  matrix(const_cast<MatrixBase *>(matrix)),
1000  a_row(row),
1001  a_index(index)
1002  {
1003  visit_present_row ();
1004  }
1005 
1006 
1007  inline
1009  Accessor (const Accessor &a)
1010  :
1011  matrix(a.matrix),
1012  a_row(a.a_row),
1013  a_index(a.a_index),
1014  colnum_cache (a.colnum_cache),
1015  value_cache (a.value_cache)
1016  {}
1017 
1018 
1019  inline
1022  {
1023  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1024  return a_row;
1025  }
1026 
1027 
1028  inline
1031  {
1032  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1033  return (*colnum_cache)[a_index];
1034  }
1035 
1036 
1037  inline
1040  {
1041  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1042  return a_index;
1043  }
1044 
1045 
1046  inline
1047  PetscScalar
1049  {
1050  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1051  return (*value_cache)[a_index];
1052  }
1053 
1054 
1055  inline
1057  const_iterator(const MatrixBase *matrix,
1058  const size_type row,
1059  const size_type index)
1060  :
1061  accessor(matrix, row, index)
1062  {}
1063 
1064 
1065 
1066  inline
1067  const_iterator &
1069  {
1070  Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
1071 
1072  ++accessor.a_index;
1073 
1074  // if at end of line: do one step, then
1075  // cycle until we find a row with a
1076  // nonzero number of entries
1077  if (accessor.a_index >= accessor.colnum_cache->size())
1078  {
1079  accessor.a_index = 0;
1080  ++accessor.a_row;
1081 
1082  while ((accessor.a_row < accessor.matrix->m())
1083  &&
1085  ++accessor.a_row;
1086 
1088  }
1089  return *this;
1090  }
1091 
1092 
1093  inline
1094  const_iterator
1096  {
1097  const const_iterator old_state = *this;
1098  ++(*this);
1099  return old_state;
1100  }
1101 
1102 
1103  inline
1104  const const_iterator::Accessor &
1106  {
1107  return accessor;
1108  }
1109 
1110 
1111  inline
1112  const const_iterator::Accessor *
1114  {
1115  return &accessor;
1116  }
1117 
1118 
1119  inline
1120  bool
1122  operator == (const const_iterator &other) const
1123  {
1124  return (accessor.a_row == other.accessor.a_row &&
1125  accessor.a_index == other.accessor.a_index);
1126  }
1127 
1128 
1129  inline
1130  bool
1132  operator != (const const_iterator &other) const
1133  {
1134  return ! (*this == other);
1135  }
1136 
1137 
1138  inline
1139  bool
1141  operator < (const const_iterator &other) const
1142  {
1143  return (accessor.row() < other.accessor.row() ||
1144  (accessor.row() == other.accessor.row() &&
1145  accessor.index() < other.accessor.index()));
1146  }
1147 
1148  }
1149 
1150 
1151 
1152  // Inline the set() and add()
1153  // functions, since they will be
1154  // called frequently, and the
1155  // compiler can optimize away
1156  // some unnecessary loops when
1157  // the sizes are given at
1158  // compile time.
1159  inline
1160  void
1161  MatrixBase::set (const size_type i,
1162  const size_type j,
1163  const PetscScalar value)
1164  {
1165  AssertIsFinite(value);
1166 
1167  set (i, 1, &j, &value, false);
1168  }
1169 
1170 
1171 
1172  inline
1173  void
1174  MatrixBase::set (const std::vector<size_type> &indices,
1175  const FullMatrix<PetscScalar> &values,
1176  const bool elide_zero_values)
1177  {
1178  Assert (indices.size() == values.m(),
1179  ExcDimensionMismatch(indices.size(), values.m()));
1180  Assert (values.m() == values.n(), ExcNotQuadratic());
1181 
1182  for (size_type i=0; i<indices.size(); ++i)
1183  set (indices[i], indices.size(), &indices[0], &values(i,0),
1184  elide_zero_values);
1185  }
1186 
1187 
1188 
1189  inline
1190  void
1191  MatrixBase::set (const std::vector<size_type> &row_indices,
1192  const std::vector<size_type> &col_indices,
1193  const FullMatrix<PetscScalar> &values,
1194  const bool elide_zero_values)
1195  {
1196  Assert (row_indices.size() == values.m(),
1197  ExcDimensionMismatch(row_indices.size(), values.m()));
1198  Assert (col_indices.size() == values.n(),
1199  ExcDimensionMismatch(col_indices.size(), values.n()));
1200 
1201  for (size_type i=0; i<row_indices.size(); ++i)
1202  set (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1203  elide_zero_values);
1204  }
1205 
1206 
1207 
1208  inline
1209  void
1210  MatrixBase::set (const size_type row,
1211  const std::vector<size_type> &col_indices,
1212  const std::vector<PetscScalar> &values,
1213  const bool elide_zero_values)
1214  {
1215  Assert (col_indices.size() == values.size(),
1216  ExcDimensionMismatch(col_indices.size(), values.size()));
1217 
1218  set (row, col_indices.size(), &col_indices[0], &values[0],
1219  elide_zero_values);
1220  }
1221 
1222 
1223 
1224  inline
1225  void
1226  MatrixBase::set (const size_type row,
1227  const size_type n_cols,
1228  const size_type *col_indices,
1229  const PetscScalar *values,
1230  const bool elide_zero_values)
1231  {
1232  (void)elide_zero_values;
1233 
1234  prepare_action(VectorOperation::insert);
1235 
1236  const PetscInt petsc_i = row;
1237  PetscInt *col_index_ptr;
1238 
1239  PetscScalar const *col_value_ptr;
1240  int n_columns;
1241 
1242  // If we don't elide zeros, the pointers are already available...
1243 #ifndef PETSC_USE_64BIT_INDICES
1244  if (elide_zero_values == false)
1245  {
1246  col_index_ptr = (int *)col_indices;
1247  col_value_ptr = values;
1248  n_columns = n_cols;
1249  }
1250  else
1251 #endif
1252  {
1253  // Otherwise, extract nonzero values in each row and get the
1254  // respective index.
1255  if (column_indices.size() < n_cols)
1256  {
1257  column_indices.resize(n_cols);
1258  column_values.resize(n_cols);
1259  }
1260 
1261  n_columns = 0;
1262  for (size_type j=0; j<n_cols; ++j)
1263  {
1264  const PetscScalar value = values[j];
1265  AssertIsFinite(value);
1266  if (value != PetscScalar())
1267  {
1268  column_indices[n_columns] = col_indices[j];
1269  column_values[n_columns] = value;
1270  n_columns++;
1271  }
1272  }
1273  Assert(n_columns <= (int)n_cols, ExcInternalError());
1274 
1275  col_index_ptr = &column_indices[0];
1276  col_value_ptr = &column_values[0];
1277  }
1278 
1279  const int ierr
1280  = MatSetValues (matrix, 1, &petsc_i, n_columns, col_index_ptr,
1281  col_value_ptr, INSERT_VALUES);
1282  AssertThrow (ierr == 0, ExcPETScError(ierr));
1283  }
1284 
1285 
1286 
1287  inline
1288  void
1289  MatrixBase::add (const size_type i,
1290  const size_type j,
1291  const PetscScalar value)
1292  {
1293 
1294  AssertIsFinite(value);
1295 
1296  if (value == PetscScalar())
1297  {
1298  // we have to check after using Insert/Add in any case to be
1299  // consistent with the MPI communication model (see the comments in
1300  // the documentation of TrilinosWrappers::Vector), but we can save
1301  // some work if the addend is zero. However, these actions are done
1302  // in case we pass on to the other function.
1303  prepare_action(VectorOperation::add);
1304 
1305  return;
1306  }
1307  else
1308  add (i, 1, &j, &value, false);
1309  }
1310 
1311 
1312 
1313  inline
1314  void
1315  MatrixBase::add (const std::vector<size_type> &indices,
1316  const FullMatrix<PetscScalar> &values,
1317  const bool elide_zero_values)
1318  {
1319  Assert (indices.size() == values.m(),
1320  ExcDimensionMismatch(indices.size(), values.m()));
1321  Assert (values.m() == values.n(), ExcNotQuadratic());
1322 
1323  for (size_type i=0; i<indices.size(); ++i)
1324  add (indices[i], indices.size(), &indices[0], &values(i,0),
1325  elide_zero_values);
1326  }
1327 
1328 
1329 
1330  inline
1331  void
1332  MatrixBase::add (const std::vector<size_type> &row_indices,
1333  const std::vector<size_type> &col_indices,
1334  const FullMatrix<PetscScalar> &values,
1335  const bool elide_zero_values)
1336  {
1337  Assert (row_indices.size() == values.m(),
1338  ExcDimensionMismatch(row_indices.size(), values.m()));
1339  Assert (col_indices.size() == values.n(),
1340  ExcDimensionMismatch(col_indices.size(), values.n()));
1341 
1342  for (size_type i=0; i<row_indices.size(); ++i)
1343  add (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1344  elide_zero_values);
1345  }
1346 
1347 
1348 
1349  inline
1350  void
1351  MatrixBase::add (const size_type row,
1352  const std::vector<size_type> &col_indices,
1353  const std::vector<PetscScalar> &values,
1354  const bool elide_zero_values)
1355  {
1356  Assert (col_indices.size() == values.size(),
1357  ExcDimensionMismatch(col_indices.size(), values.size()));
1358 
1359  add (row, col_indices.size(), &col_indices[0], &values[0],
1360  elide_zero_values);
1361  }
1362 
1363 
1364 
1365  inline
1366  void
1367  MatrixBase::add (const size_type row,
1368  const size_type n_cols,
1369  const size_type *col_indices,
1370  const PetscScalar *values,
1371  const bool elide_zero_values,
1372  const bool /*col_indices_are_sorted*/)
1373  {
1374  (void)elide_zero_values;
1375 
1376  prepare_action(VectorOperation::add);
1377 
1378  const PetscInt petsc_i = row;
1379  PetscInt *col_index_ptr;
1380 
1381  PetscScalar const *col_value_ptr;
1382  int n_columns;
1383 
1384  // If we don't elide zeros, the pointers are already available...
1385 #ifndef PETSC_USE_64BIT_INDICES
1386  if (elide_zero_values == false)
1387  {
1388  col_index_ptr = (int *)col_indices;
1389  col_value_ptr = values;
1390  n_columns = n_cols;
1391  }
1392  else
1393 #endif
1394  {
1395  // Otherwise, extract nonzero values in each row and get the
1396  // respective index.
1397  if (column_indices.size() < n_cols)
1398  {
1399  column_indices.resize(n_cols);
1400  column_values.resize(n_cols);
1401  }
1402 
1403  n_columns = 0;
1404  for (size_type j=0; j<n_cols; ++j)
1405  {
1406  const PetscScalar value = values[j];
1407  AssertIsFinite(value);
1408  if (value != PetscScalar())
1409  {
1410  column_indices[n_columns] = col_indices[j];
1411  column_values[n_columns] = value;
1412  n_columns++;
1413  }
1414  }
1415  Assert(n_columns <= (int)n_cols, ExcInternalError());
1416 
1417  col_index_ptr = &column_indices[0];
1418  col_value_ptr = &column_values[0];
1419  }
1420 
1421  const int ierr
1422  = MatSetValues (matrix, 1, &petsc_i, n_columns, col_index_ptr,
1423  col_value_ptr, ADD_VALUES);
1424  AssertThrow (ierr == 0, ExcPETScError(ierr));
1425  }
1426 
1427 
1428 
1429 
1430 
1431 
1432  inline
1433  PetscScalar
1434  MatrixBase::operator() (const size_type i,
1435  const size_type j) const
1436  {
1437  return el(i,j);
1438  }
1439 
1440 
1441 
1442  inline
1444  MatrixBase::begin() const
1445  {
1446  return const_iterator(this, 0, 0);
1447  }
1448 
1449 
1450  inline
1452  MatrixBase::end() const
1453  {
1454  return const_iterator(this, m(), 0);
1455  }
1456 
1457 
1458  inline
1460  MatrixBase::begin(const size_type r) const
1461  {
1462  Assert (r < m(), ExcIndexRange(r, 0, m()));
1463  if (row_length(r) > 0)
1464  return const_iterator(this, r, 0);
1465  else
1466  return end (r);
1467  }
1468 
1469 
1470  inline
1472  MatrixBase::end(const size_type r) const
1473  {
1474  Assert (r < m(), ExcIndexRange(r, 0, m()));
1475 
1476  // place the iterator on the first entry
1477  // past this line, or at the end of the
1478  // matrix
1479  for (size_type i=r+1; i<m(); ++i)
1480  if (row_length(i) > 0)
1481  return const_iterator(this, i, 0);
1482 
1483  // if there is no such line, then take the
1484  // end iterator of the matrix
1485  return end();
1486  }
1487 
1488 
1489 
1490  inline
1491  bool
1492  MatrixBase::in_local_range (const size_type index) const
1493  {
1494  PetscInt begin, end;
1495 
1496  const int ierr = MatGetOwnershipRange (static_cast<const Mat &>(matrix),
1497  &begin, &end);
1498  AssertThrow (ierr == 0, ExcPETScError(ierr));
1499 
1500  return ((index >= static_cast<size_type>(begin)) &&
1501  (index < static_cast<size_type>(end)));
1502  }
1503 
1504 
1505 
1506  inline
1507  void
1508  MatrixBase::prepare_action(const VectorOperation::values new_action)
1509  {
1510  if (last_action == VectorOperation::unknown)
1511  last_action = new_action;
1512 
1513  Assert (last_action == new_action, ExcWrongMode (last_action, new_action));
1514  }
1515 
1516 
1517 
1518  inline
1519  void
1521  {
1522  // compress() sets the last action to none, which allows us to check if there
1523  // are pending add/insert operations:
1524  AssertThrow (last_action == VectorOperation::unknown,
1525  ExcMessage("Error: missing compress() call."));
1526  }
1527 
1528 
1529 
1530  inline
1531  void
1533  {
1534  prepare_action(VectorOperation::add);
1535  }
1536 
1537 
1538 
1539  inline
1540  void
1542  {
1543  prepare_action(VectorOperation::insert);
1544  }
1545 
1546 #endif // DOXYGEN
1547 }
1548 
1549 
1550 DEAL_II_NAMESPACE_CLOSE
1551 
1552 
1553 #endif // DEAL_II_WITH_PETSC
1554 
1555 
1556 /*---------------------------- petsc_matrix_base.h ---------------------------*/
1557 
1558 #endif
1559 /*---------------------------- petsc_matrix_base.h ---------------------------*/
bool operator==(const const_iterator &) const
void add(const size_type i, const size_type j, const PetscScalar value)
const_iterator begin() const
::ExceptionBase & ExcMessage(std::string arg1)
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:358
std::vector< PetscInt > column_indices
const_iterator(const MatrixBase *matrix, const size_type row, const size_type index)
std::vector< PetscScalar > column_values
MatrixIterators::const_iterator const_iterator
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:542
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:294
types::global_dof_index size_type
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
VectorOperation::values last_action
void prepare_action(const VectorOperation::values new_action)
bool in_local_range(const size_type index) const
void set(const size_type i, const size_type j, const PetscScalar value)
DeclException3(ExcAccessToNonlocalRow, int, int, int,<< "You tried to access row "<< arg1<< " of a distributed matrix, but only rows "<< arg2<< " through "<< arg3<< " are stored locally and can be accessed.")
DeclException2(ExcInvalidIndexWithinRow, int, int,<< "Attempt to access element "<< arg2<< " of row "<< arg1<< " which doesn't have that many elements.")
bool operator!=(const const_iterator &) const
#define AssertIsFinite(number)
Definition: exceptions.h:1096