Reference documentation for deal.II version 8.5.1
trilinos_sparsity_pattern.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 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__trilinos_sparsity_pattern_h
17 #define dealii__trilinos_sparsity_pattern_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
24 # include <deal.II/base/subscriptor.h>
25 # include <deal.II/base/index_set.h>
26 # include <deal.II/lac/exceptions.h>
27 
28 # include <vector>
29 # include <cmath>
30 # include <memory>
31 
32 # include <deal.II/base/std_cxx11/shared_ptr.h>
33 
34 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
35 # include <Epetra_FECrsGraph.h>
36 # include <Epetra_Map.h>
37 # ifdef DEAL_II_WITH_MPI
38 # include <Epetra_MpiComm.h>
39 # include "mpi.h"
40 # else
41 # include "Epetra_SerialComm.h"
42 # endif
43 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
44 
45 
46 DEAL_II_NAMESPACE_OPEN
47 
48 // forward declarations
49 class SparsityPattern;
51 
52 namespace TrilinosWrappers
53 {
54  // forward declarations
55  class SparsityPattern;
56  class SparseMatrix;
57 
58  namespace SparsityPatternIterators
59  {
60  // forward declaration
61  class Iterator;
62 
76  class Accessor
77  {
78  public:
82  typedef ::types::global_dof_index size_type;
83 
88  const size_type row,
89  const size_type index);
90 
94  size_type row() const;
95 
99  size_type index() const;
100 
104  size_type column() const;
105 
110 
116  << "You tried to access row " << arg1
117  << " of a distributed sparsity pattern, "
118  << " but only rows " << arg2 << " through " << arg3
119  << " are stored locally and can be accessed.");
120 
121  private:
126 
131 
136 
149  std_cxx11::shared_ptr<const std::vector<size_type> > colnum_cache;
150 
156  void visit_present_row ();
157 
161  friend class Iterator;
162  };
163 
169  class Iterator
170  {
171  public:
175  typedef ::types::global_dof_index size_type;
176 
181  Iterator (const SparsityPattern *sparsity_pattern,
182  const size_type row,
183  const size_type index);
184 
188  Iterator (const Iterator &i);
189 
194 
198  Iterator operator++ (int);
199 
203  const Accessor &operator* () const;
204 
208  const Accessor *operator-> () const;
209 
214  bool operator == (const Iterator &) const;
215 
219  bool operator != (const Iterator &) const;
220 
226  bool operator < (const Iterator &) const;
227 
233  << "Attempt to access element " << arg2
234  << " of row " << arg1
235  << " which doesn't have that many elements.");
236 
237  private:
242 
244  };
245 
246  }
247 
248 
268  {
269  public:
270 
274  typedef ::types::global_dof_index size_type;
275 
280 
288  SparsityPattern ();
289 
304  SparsityPattern (const size_type m,
305  const size_type n,
306  const size_type n_entries_per_row = 0);
307 
316  SparsityPattern (const size_type m,
317  const size_type n,
318  const std::vector<size_type> &n_entries_per_row);
319 
324  SparsityPattern (const SparsityPattern &input_sparsity_pattern);
325 
329  virtual ~SparsityPattern ();
330 
342  void
343  reinit (const size_type m,
344  const size_type n,
345  const size_type n_entries_per_row = 0);
346 
355  void
356  reinit (const size_type m,
357  const size_type n,
358  const std::vector<size_type> &n_entries_per_row);
359 
364  void
365  copy_from (const SparsityPattern &input_sparsity_pattern);
366 
372  template<typename SparsityPatternType>
373  void
374  copy_from (const SparsityPatternType &nontrilinos_sparsity_pattern);
375 
382  SparsityPattern &operator = (const SparsityPattern &input_sparsity_pattern);
383 
391  void clear ();
392 
402  void compress ();
404 
408 
423  SparsityPattern (const Epetra_Map &parallel_partitioning,
424  const size_type n_entries_per_row = 0) DEAL_II_DEPRECATED;
425 
438  SparsityPattern (const Epetra_Map &parallel_partitioning,
439  const std::vector<size_type> &n_entries_per_row) DEAL_II_DEPRECATED;
440 
459  SparsityPattern (const Epetra_Map &row_parallel_partitioning,
460  const Epetra_Map &col_parallel_partitioning,
461  const size_type n_entries_per_row = 0) DEAL_II_DEPRECATED;
462 
476  SparsityPattern (const Epetra_Map &row_parallel_partitioning,
477  const Epetra_Map &col_parallel_partitioning,
478  const std::vector<size_type> &n_entries_per_row) DEAL_II_DEPRECATED;
479 
496  void
497  reinit (const Epetra_Map &parallel_partitioning,
498  const size_type n_entries_per_row = 0) DEAL_II_DEPRECATED;
499 
512  void
513  reinit (const Epetra_Map &parallel_partitioning,
514  const std::vector<size_type> &n_entries_per_row) DEAL_II_DEPRECATED;
515 
534  void
535  reinit (const Epetra_Map &row_parallel_partitioning,
536  const Epetra_Map &col_parallel_partitioning,
537  const size_type n_entries_per_row = 0) DEAL_II_DEPRECATED;
538 
552  void
553  reinit (const Epetra_Map &row_parallel_partitioning,
554  const Epetra_Map &col_parallel_partitioning,
555  const std::vector<size_type> &n_entries_per_row) DEAL_II_DEPRECATED;
556 
567  template<typename SparsityPatternType>
568  void
569  reinit (const Epetra_Map &row_parallel_partitioning,
570  const Epetra_Map &col_parallel_partitioning,
571  const SparsityPatternType &nontrilinos_sparsity_pattern,
572  const bool exchange_data = false) DEAL_II_DEPRECATED;
573 
584  template<typename SparsityPatternType>
585  void
586  reinit (const Epetra_Map &parallel_partitioning,
587  const SparsityPatternType &nontrilinos_sparsity_pattern,
588  const bool exchange_data = false) DEAL_II_DEPRECATED;
590 
594 
606  SparsityPattern (const IndexSet &parallel_partitioning,
607  const MPI_Comm &communicator = MPI_COMM_WORLD,
608  const size_type n_entries_per_row = 0);
609 
620  SparsityPattern (const IndexSet &parallel_partitioning,
621  const MPI_Comm &communicator,
622  const std::vector<size_type> &n_entries_per_row);
623 
638  SparsityPattern (const IndexSet &row_parallel_partitioning,
639  const IndexSet &col_parallel_partitioning,
640  const MPI_Comm &communicator = MPI_COMM_WORLD,
641  const size_type n_entries_per_row = 0);
642 
654  SparsityPattern (const IndexSet &row_parallel_partitioning,
655  const IndexSet &col_parallel_partitioning,
656  const MPI_Comm &communicator,
657  const std::vector<size_type> &n_entries_per_row);
658 
685  SparsityPattern (const IndexSet &row_parallel_partitioning,
686  const IndexSet &col_parallel_partitioning,
687  const IndexSet &writable_rows,
688  const MPI_Comm &communicator = MPI_COMM_WORLD,
689  const size_type n_entries_per_row = 0);
690 
706  void
707  reinit (const IndexSet &parallel_partitioning,
708  const MPI_Comm &communicator = MPI_COMM_WORLD,
709  const size_type n_entries_per_row = 0);
710 
721  void
722  reinit (const IndexSet &parallel_partitioning,
723  const MPI_Comm &communicator,
724  const std::vector<size_type> &n_entries_per_row);
725 
742  void
743  reinit (const IndexSet &row_parallel_partitioning,
744  const IndexSet &col_parallel_partitioning,
745  const MPI_Comm &communicator = MPI_COMM_WORLD,
746  const size_type n_entries_per_row = 0);
747 
773  void
774  reinit (const IndexSet &row_parallel_partitioning,
775  const IndexSet &col_parallel_partitioning,
776  const IndexSet &writeable_rows,
777  const MPI_Comm &communicator = MPI_COMM_WORLD,
778  const size_type n_entries_per_row = 0);
779 
784  void
785  reinit (const IndexSet &row_parallel_partitioning,
786  const IndexSet &col_parallel_partitioning,
787  const MPI_Comm &communicator,
788  const std::vector<size_type> &n_entries_per_row);
789 
799  template<typename SparsityPatternType>
800  void
801  reinit (const IndexSet &row_parallel_partitioning,
802  const IndexSet &col_parallel_partitioning,
803  const SparsityPatternType &nontrilinos_sparsity_pattern,
804  const MPI_Comm &communicator = MPI_COMM_WORLD,
805  const bool exchange_data = false);
806 
815  template<typename SparsityPatternType>
816  void
817  reinit (const IndexSet &parallel_partitioning,
818  const SparsityPatternType &nontrilinos_sparsity_pattern,
819  const MPI_Comm &communicator = MPI_COMM_WORLD,
820  const bool exchange_data = false);
822 
826 
831  bool is_compressed () const;
832 
836  unsigned int max_entries_per_row () const;
837 
841  size_type n_rows () const;
842 
846  size_type n_cols () const;
847 
857  unsigned int local_size () const;
858 
867  std::pair<size_type, size_type>
868  local_range () const;
869 
874  bool in_local_range (const size_type index) const;
875 
879  size_type n_nonzero_elements () const;
880 
884  size_type row_length (const size_type row) const;
885 
892  size_type bandwidth () const;
893 
898  bool empty () const;
899 
904  bool exists (const size_type i,
905  const size_type j) const;
906 
911  std::size_t memory_consumption () const;
912 
914 
921  void add (const size_type i,
922  const size_type j);
923 
924 
928  template <typename ForwardIterator>
929  void add_entries (const size_type row,
930  ForwardIterator begin,
931  ForwardIterator end,
932  const bool indices_are_sorted = false);
934 
938 
943  const Epetra_FECrsGraph &trilinos_sparsity_pattern () const;
944 
953  const Epetra_Map &domain_partitioner () const DEAL_II_DEPRECATED;
954 
963  const Epetra_Map &range_partitioner () const DEAL_II_DEPRECATED;
964 
972  const Epetra_Map &row_partitioner () const DEAL_II_DEPRECATED;
973 
983  const Epetra_Map &col_partitioner () const DEAL_II_DEPRECATED;
984 
990  const Epetra_Comm &trilinos_communicator () const DEAL_II_DEPRECATED;
991 
995  MPI_Comm get_mpi_communicator () const;
997 
1002 
1009 
1016 
1018 
1023 
1027  const_iterator begin () const;
1028 
1032  const_iterator end () const;
1033 
1042  const_iterator begin (const size_type r) const;
1043 
1052  const_iterator end (const size_type r) const;
1053 
1055 
1059 
1065  void write_ascii ();
1066 
1074  void print (std::ostream &out,
1075  const bool write_extended_trilinos_info = false) const;
1076 
1091  void print_gnuplot (std::ostream &out) const;
1092 
1094 
1102  int,
1103  << "An error with error number " << arg1
1104  << " occurred while calling a Trilinos function");
1105 
1111  << "The entry with index <" << arg1 << ',' << arg2
1112  << "> does not exist.");
1113 
1118 
1124  << "You tried to access element (" << arg1
1125  << "/" << arg2 << ")"
1126  << " of a distributed matrix, but only rows "
1127  << arg3 << " through " << arg4
1128  << " are stored locally and can be accessed.");
1129 
1135  << "You tried to access element (" << arg1
1136  << "/" << arg2 << ")"
1137  << " of a sparse matrix, but it appears to not"
1138  << " exist in the Trilinos sparsity pattern.");
1140  private:
1141 
1146  std_cxx11::shared_ptr<Epetra_Map> column_space_map;
1147 
1153  std_cxx11::shared_ptr<Epetra_FECrsGraph> graph;
1154 
1161  std_cxx11::shared_ptr<Epetra_CrsGraph> nonlocal_graph;
1162 
1163  friend class TrilinosWrappers::SparseMatrix;
1164  friend class SparsityPatternIterators::Accessor;
1165  friend class SparsityPatternIterators::Iterator;
1166  };
1167 
1168 
1169 
1170 // -------------------------- inline and template functions ----------------------
1171 
1172 
1173 #ifndef DOXYGEN
1174 
1175  namespace SparsityPatternIterators
1176  {
1177 
1178  inline
1179  Accessor::Accessor (const SparsityPattern *sp,
1180  const size_type row,
1181  const size_type index)
1182  :
1183  sparsity_pattern(const_cast<SparsityPattern *>(sp)),
1184  a_row(row),
1185  a_index(index)
1186  {
1187  visit_present_row ();
1188  }
1189 
1190 
1191 
1192  inline
1193  Accessor::size_type
1194  Accessor::row() const
1195  {
1196  Assert (a_row < sparsity_pattern->n_rows(), ExcBeyondEndOfSparsityPattern());
1197  return a_row;
1198  }
1199 
1200 
1201 
1202  inline
1203  Accessor::size_type
1204  Accessor::column() const
1205  {
1206  Assert (a_row < sparsity_pattern->n_rows(), ExcBeyondEndOfSparsityPattern());
1207  return (*colnum_cache)[a_index];
1208  }
1209 
1210 
1211 
1212  inline
1213  Accessor::size_type
1214  Accessor::index() const
1215  {
1216  Assert (a_row < sparsity_pattern->n_rows(), ExcBeyondEndOfSparsityPattern());
1217  return a_index;
1218  }
1219 
1220 
1221 
1222  inline
1223  Iterator::Iterator(const SparsityPattern *sp,
1224  const size_type row,
1225  const size_type index)
1226  :
1227  accessor(sp, row, index)
1228  {}
1229 
1230 
1231  inline
1232  Iterator::Iterator(const Iterator &i)
1233  :
1234  accessor(i.accessor)
1235  {}
1236 
1237 
1238 
1239  inline
1240  Iterator &
1241  Iterator::operator++ ()
1242  {
1243  Assert (accessor.a_row < accessor.sparsity_pattern->n_rows(),
1244  ExcIteratorPastEnd());
1245 
1246  ++accessor.a_index;
1247 
1248  // If at end of line: do one
1249  // step, then cycle until we
1250  // find a row with a nonzero
1251  // number of entries.
1252  if (accessor.a_index >= accessor.colnum_cache->size())
1253  {
1254  accessor.a_index = 0;
1255  ++accessor.a_row;
1256 
1257  while ((accessor.a_row < accessor.sparsity_pattern->n_rows())
1258  &&
1259  (accessor.sparsity_pattern->row_length(accessor.a_row) == 0))
1260  ++accessor.a_row;
1261 
1262  accessor.visit_present_row();
1263  }
1264  return *this;
1265  }
1266 
1267 
1268 
1269  inline
1270  Iterator
1271  Iterator::operator++ (int)
1272  {
1273  const Iterator old_state = *this;
1274  ++(*this);
1275  return old_state;
1276  }
1277 
1278 
1279 
1280  inline
1281  const Accessor &
1282  Iterator::operator* () const
1283  {
1284  return accessor;
1285  }
1286 
1287 
1288 
1289  inline
1290  const Accessor *
1291  Iterator::operator-> () const
1292  {
1293  return &accessor;
1294  }
1295 
1296 
1297 
1298  inline
1299  bool
1300  Iterator::operator == (const Iterator &other) const
1301  {
1302  return (accessor.a_row == other.accessor.a_row &&
1303  accessor.a_index == other.accessor.a_index);
1304  }
1305 
1306 
1307 
1308  inline
1309  bool
1310  Iterator::operator != (const Iterator &other) const
1311  {
1312  return ! (*this == other);
1313  }
1314 
1315 
1316 
1317  inline
1318  bool
1319  Iterator::operator < (const Iterator &other) const
1320  {
1321  return (accessor.row() < other.accessor.row() ||
1322  (accessor.row() == other.accessor.row() &&
1323  accessor.index() < other.accessor.index()));
1324  }
1325 
1326  }
1327 
1328 
1329 
1330  inline
1332  SparsityPattern::begin() const
1333  {
1334  return const_iterator(this, 0, 0);
1335  }
1336 
1337 
1338 
1339  inline
1341  SparsityPattern::end() const
1342  {
1343  return const_iterator(this, n_rows(), 0);
1344  }
1345 
1346 
1347 
1348  inline
1350  SparsityPattern::begin(const size_type r) const
1351  {
1352  Assert (r < n_rows(), ExcIndexRangeType<size_type>(r, 0, n_rows()));
1353  if (row_length(r) > 0)
1354  return const_iterator(this, r, 0);
1355  else
1356  return end (r);
1357  }
1358 
1359 
1360 
1361  inline
1363  SparsityPattern::end(const size_type r) const
1364  {
1365  Assert (r < n_rows(), ExcIndexRangeType<size_type>(r, 0, n_rows()));
1366 
1367  // place the iterator on the first entry
1368  // past this line, or at the end of the
1369  // matrix
1370  for (size_type i=r+1; i<n_rows(); ++i)
1371  if (row_length(i) > 0)
1372  return const_iterator(this, i, 0);
1373 
1374  // if there is no such line, then take the
1375  // end iterator of the matrix
1376  return end();
1377  }
1378 
1379 
1380 
1381  inline
1382  bool
1383  SparsityPattern::in_local_range (const size_type index) const
1384  {
1385  TrilinosWrappers::types::int_type begin, end;
1386 #ifndef DEAL_II_WITH_64BIT_INDICES
1387  begin = graph->RowMap().MinMyGID();
1388  end = graph->RowMap().MaxMyGID()+1;
1389 #else
1390  begin = graph->RowMap().MinMyGID64();
1391  end = graph->RowMap().MaxMyGID64()+1;
1392 #endif
1393 
1394  return ((index >= static_cast<size_type>(begin)) &&
1395  (index < static_cast<size_type>(end)));
1396  }
1397 
1398 
1399 
1400  inline
1401  bool
1403  {
1404  return graph->Filled();
1405  }
1406 
1407 
1408 
1409  inline
1410  bool
1411  SparsityPattern::empty () const
1412  {
1413  return ((n_rows() == 0) && (n_cols() == 0));
1414  }
1415 
1416 
1417 
1418  inline
1419  void
1421  const size_type j)
1422  {
1423  add_entries (i, &j, &j+1);
1424  }
1425 
1426 
1427 
1428  template <typename ForwardIterator>
1429  inline
1430  void
1432  ForwardIterator begin,
1433  ForwardIterator end,
1434  const bool /*indices_are_sorted*/)
1435  {
1436  if (begin == end)
1437  return;
1438 
1439  // verify that the size of the data type Trilinos expects matches that the
1440  // iterator points to. we allow for some slippage between signed and
1441  // unsigned and only compare that they are both either 32 or 64 bit. to
1442  // write this test properly, not that we cannot compare the size of
1443  // '*begin' because 'begin' may be an iterator and '*begin' may be an
1444  // accessor class. consequently, we need to somehow get an actual value
1445  // from it which we can by evaluating an expression such as when
1446  // multiplying the value produced by 2
1447  Assert (sizeof(TrilinosWrappers::types::int_type) ==
1448  sizeof((*begin)*2),
1449  ExcNotImplemented());
1450 
1451  TrilinosWrappers::types::int_type *col_index_ptr =
1452  (TrilinosWrappers::types::int_type *)(&*begin);
1453  const int n_cols = static_cast<int>(end - begin);
1454 
1455  int ierr;
1456  if ( graph->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row)) != -1)
1457  ierr = graph->InsertGlobalIndices (row, n_cols, col_index_ptr);
1458  else if (nonlocal_graph.get() != 0)
1459  {
1460  // this is the case when we have explicitly set the off-processor rows
1461  // and want to create a separate matrix object for them (to retain
1462  // thread-safety)
1463  Assert (nonlocal_graph->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row)) != -1,
1464  ExcMessage("Attempted to write into off-processor matrix row "
1465  "that has not be specified as being writable upon "
1466  "initialization"));
1467  ierr = nonlocal_graph->InsertGlobalIndices (row, n_cols, col_index_ptr);
1468  }
1469  else
1470  ierr = graph->InsertGlobalIndices
1471  (1, (TrilinosWrappers::types::int_type *)&row, n_cols, col_index_ptr);
1472 
1473  AssertThrow (ierr >= 0, ExcTrilinosError(ierr));
1474  }
1475 
1476 
1477 
1478  inline
1479  const Epetra_FECrsGraph &
1481  {
1482  return *graph;
1483  }
1484 
1485 
1486 
1487  inline
1488  IndexSet
1490  {
1491  return IndexSet(static_cast<const Epetra_Map &>(graph->DomainMap()));
1492  }
1493 
1494 
1495 
1496  inline
1497  IndexSet
1499  {
1500  return IndexSet(static_cast<const Epetra_Map &>(graph->RangeMap()));
1501  }
1502 
1503 #endif // DOXYGEN
1504 }
1505 
1506 
1507 DEAL_II_NAMESPACE_CLOSE
1508 
1509 
1510 #endif // DEAL_II_WITH_TRILINOS
1511 
1512 
1513 /*-------------------- trilinos_sparsity_pattern.h --------------------*/
1514 
1515 #endif
1516 /*-------------------- trilinos_sparsity_pattern.h --------------------*/
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
static ::ExceptionBase & ExcBeyondEndOfSparsityPattern()
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
size_type row_length(const size_type row) const
IndexSet locally_owned_range_indices() const
const Epetra_FECrsGraph & trilinos_sparsity_pattern() const
static ::ExceptionBase & ExcSourceEqualsDestination()
std_cxx11::shared_ptr< Epetra_FECrsGraph > graph
std_cxx11::shared_ptr< Epetra_CrsGraph > nonlocal_graph
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
void add(const size_type i, const size_type j)
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false)
const_iterator begin() const
bool exists(const size_type i, const size_type j) const
static ::ExceptionBase & ExcTrilinosError(int arg1)
std_cxx11::shared_ptr< Epetra_Map > column_space_map
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
std_cxx11::shared_ptr< const std::vector< size_type > > colnum_cache
#define Assert(cond, exc)
Definition: exceptions.h:313
void copy_from(const SparsityPattern &input_sparsity_pattern)
static ::ExceptionBase & ExcInvalidIndexWithinRow(size_type arg1, size_type arg2)
#define DeclException0(Exception0)
Definition: exceptions.h:541
const Epetra_Comm & trilinos_communicator() const 1
const Epetra_Map & domain_partitioner() const 1
void print_gnuplot(std::ostream &out) const
const_iterator end() const
SparsityPatternIterators::Iterator const_iterator
static ::ExceptionBase & ExcIteratorPastEnd()
const Epetra_Map & range_partitioner() const 1
SparsityPattern & operator=(const SparsityPattern &input_sparsity_pattern)
Iterator(const SparsityPattern *sparsity_pattern, const size_type row, const size_type index)
static ::ExceptionBase & ExcAccessToNonlocalRow(size_type arg1, size_type arg2, size_type arg3)
const Epetra_Map & col_partitioner() const 1
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:600
static ::ExceptionBase & ExcNotImplemented()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:588
bool in_local_range(const size_type index) const
const Epetra_Map & row_partitioner() const 1
void reinit(const size_type m, const size_type n, const size_type n_entries_per_row=0)
IndexSet locally_owned_domain_indices() const
std::pair< size_type, size_type > local_range() const
Accessor(const SparsityPattern *sparsity_pattern, const size_type row, const size_type index)