16 #ifndef dealii__petsc_matrix_base_h 17 #define dealii__petsc_matrix_base_h 20 #include <deal.II/base/config.h> 22 #ifdef DEAL_II_WITH_PETSC 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> 29 # include <petscmat.h> 30 # include <deal.II/base/std_cxx11/shared_ptr.h> 35 DEAL_II_NAMESPACE_OPEN
46 namespace MatrixIterators
82 const size_type
index);
92 size_type
row()
const;
97 size_type
index()
const;
107 PetscScalar
value()
const;
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.");
156 std_cxx11::shared_ptr<const std::vector<PetscScalar> >
value_cache;
183 const size_type
index);
227 <<
"Attempt to access element " << arg2
228 <<
" of row " << arg1
229 <<
" which doesn't have that many elements.");
311 operator = (
const value_type d);
327 void set (
const size_type i,
329 const PetscScalar
value);
350 void set (
const std::vector<size_type> &indices,
352 const bool elide_zero_values =
false);
359 void set (
const std::vector<size_type> &row_indices,
360 const std::vector<size_type> &col_indices,
362 const bool elide_zero_values =
false);
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);
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);
412 void add (
const size_type i,
414 const PetscScalar
value);
435 void add (
const std::vector<size_type> &indices,
437 const bool elide_zero_values =
true);
444 void add (
const std::vector<size_type> &row_indices,
445 const std::vector<size_type> &col_indices,
447 const bool elide_zero_values =
true);
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);
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);
505 void clear_row (
const size_type row,
506 const PetscScalar new_diag_value = 0);
516 void clear_rows (
const std::vector<size_type> &rows,
517 const PetscScalar new_diag_value = 0);
530 void compress (
const VectorOperation::values operation);
543 PetscScalar operator () (
const size_type i,
544 const size_type j)
const;
553 PetscScalar el (
const size_type i,
554 const size_type j)
const;
565 PetscScalar diag_element (
const size_type i)
const;
570 size_type m ()
const;
575 size_type n ()
const;
585 size_type local_size ()
const;
595 std::pair<size_type, size_type>
596 local_range ()
const;
602 bool in_local_range (
const size_type
index)
const;
608 virtual const MPI_Comm &get_mpi_communicator ()
const = 0;
615 size_type n_nonzero_elements ()
const;
620 size_type row_length (
const size_type row)
const;
629 PetscReal l1_norm ()
const;
638 PetscReal linfty_norm ()
const;
644 PetscReal frobenius_norm ()
const;
666 PetscScalar matrix_norm_square (
const VectorBase &v)
const;
682 PetscScalar matrix_scalar_product (
const VectorBase &u,
686 #if DEAL_II_PETSC_VERSION_GTE(3,1,0) 691 PetscScalar trace ()
const;
697 MatrixBase &operator *= (
const PetscScalar factor);
702 MatrixBase &operator /= (
const PetscScalar factor);
709 const PetscScalar factor);
789 const_iterator begin ()
const;
794 const_iterator end ()
const;
804 const_iterator begin (
const size_type r)
const;
814 const_iterator end (
const size_type r)
const;
823 operator Mat ()
const;
834 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 839 is_symmetric (
const double tolerance = 1.e-12);
846 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 851 is_hermitian (
const double tolerance = 1.e-12);
860 void write_ascii (
const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
869 void print (std::ostream &out,
870 const bool alternative_output =
false)
const;
875 std::size_t memory_consumption()
const;
882 <<
"An error with error number " << arg1
883 <<
" occurred while calling a PETSc function");
894 <<
"You tried to do a " 899 <<
" operation but the matrix is currently in " 904 <<
" mode. You first have to call 'compress()'.");
923 void prepare_action(
const VectorOperation::values new_action);
930 void assert_is_compressed();
981 template <
class>
friend class ::BlockMatrixBase;
990 namespace MatrixIterators
997 const size_type
index)
999 matrix(const_cast<MatrixBase *>(matrix)),
1023 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1032 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1041 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1050 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
1058 const size_type row,
1059 const size_type index)
1097 const const_iterator old_state = *
this;
1134 return ! (*
this == other);
1141 operator < (
const const_iterator &other)
const 1163 const PetscScalar
value)
1167 set (i, 1, &j, &
value,
false);
1176 const bool elide_zero_values)
1178 Assert (indices.size() == values.
m(),
1179 ExcDimensionMismatch(indices.size(), values.
m()));
1180 Assert (values.
m() == values.
n(), ExcNotQuadratic());
1182 for (size_type i=0; i<indices.size(); ++i)
1183 set (indices[i], indices.size(), &indices[0], &values(i,0),
1192 const std::vector<size_type> &col_indices,
1194 const bool elide_zero_values)
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()));
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),
1211 const std::vector<size_type> &col_indices,
1212 const std::vector<PetscScalar> &values,
1213 const bool elide_zero_values)
1215 Assert (col_indices.size() == values.size(),
1216 ExcDimensionMismatch(col_indices.size(), values.size()));
1218 set (
row, col_indices.size(), &col_indices[0], &values[0],
1227 const size_type n_cols,
1228 const size_type *col_indices,
1229 const PetscScalar *values,
1230 const bool elide_zero_values)
1232 (void)elide_zero_values;
1234 prepare_action(VectorOperation::insert);
1236 const PetscInt petsc_i =
row;
1237 PetscInt *col_index_ptr;
1239 PetscScalar
const *col_value_ptr;
1243 #ifndef PETSC_USE_64BIT_INDICES 1244 if (elide_zero_values ==
false)
1246 col_index_ptr = (
int *)col_indices;
1247 col_value_ptr = values;
1255 if (column_indices.size() < n_cols)
1257 column_indices.resize(n_cols);
1258 column_values.resize(n_cols);
1262 for (size_type j=0; j<n_cols; ++j)
1264 const PetscScalar value = values[j];
1266 if (value != PetscScalar())
1268 column_indices[n_columns] = col_indices[j];
1269 column_values[n_columns] =
value;
1273 Assert(n_columns <= (
int)n_cols, ExcInternalError());
1275 col_index_ptr = &column_indices[0];
1276 col_value_ptr = &column_values[0];
1280 = MatSetValues (
matrix, 1, &petsc_i, n_columns, col_index_ptr,
1281 col_value_ptr, INSERT_VALUES);
1291 const PetscScalar value)
1296 if (value == PetscScalar())
1303 prepare_action(VectorOperation::add);
1308 add (i, 1, &j, &value,
false);
1317 const bool elide_zero_values)
1319 Assert (indices.size() == values.
m(),
1320 ExcDimensionMismatch(indices.size(), values.
m()));
1321 Assert (values.
m() == values.
n(), ExcNotQuadratic());
1323 for (size_type i=0; i<indices.size(); ++i)
1324 add (indices[i], indices.size(), &indices[0], &values(i,0),
1333 const std::vector<size_type> &col_indices,
1335 const bool elide_zero_values)
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()));
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),
1352 const std::vector<size_type> &col_indices,
1353 const std::vector<PetscScalar> &values,
1354 const bool elide_zero_values)
1356 Assert (col_indices.size() == values.size(),
1357 ExcDimensionMismatch(col_indices.size(), values.size()));
1359 add (row, col_indices.size(), &col_indices[0], &values[0],
1368 const size_type n_cols,
1369 const size_type *col_indices,
1370 const PetscScalar *values,
1371 const bool elide_zero_values,
1374 (void)elide_zero_values;
1376 prepare_action(VectorOperation::add);
1378 const PetscInt petsc_i =
row;
1379 PetscInt *col_index_ptr;
1381 PetscScalar
const *col_value_ptr;
1385 #ifndef PETSC_USE_64BIT_INDICES 1386 if (elide_zero_values ==
false)
1388 col_index_ptr = (
int *)col_indices;
1389 col_value_ptr = values;
1397 if (column_indices.size() < n_cols)
1399 column_indices.resize(n_cols);
1400 column_values.resize(n_cols);
1404 for (size_type j=0; j<n_cols; ++j)
1406 const PetscScalar value = values[j];
1408 if (value != PetscScalar())
1410 column_indices[n_columns] = col_indices[j];
1411 column_values[n_columns] =
value;
1415 Assert(n_columns <= (
int)n_cols, ExcInternalError());
1417 col_index_ptr = &column_indices[0];
1418 col_value_ptr = &column_values[0];
1422 = MatSetValues (
matrix, 1, &petsc_i, n_columns, col_index_ptr,
1423 col_value_ptr, ADD_VALUES);
1435 const size_type j)
const 1462 Assert (r < m(), ExcIndexRange(r, 0, m()));
1463 if (row_length(r) > 0)
1474 Assert (r < m(), ExcIndexRange(r, 0, m()));
1479 for (size_type i=r+1; i<m(); ++i)
1480 if (row_length(i) > 0)
1494 PetscInt begin, end;
1496 const int ierr = MatGetOwnershipRange (static_cast<const Mat &>(
matrix),
1500 return ((index >= static_cast<size_type>(begin)) &&
1501 (index < static_cast<size_type>(end)));
1510 if (last_action == VectorOperation::unknown)
1511 last_action = new_action;
1513 Assert (last_action == new_action, ExcWrongMode (last_action, new_action));
1524 AssertThrow (last_action == VectorOperation::unknown,
1525 ExcMessage(
"Error: missing compress() call."));
1534 prepare_action(VectorOperation::add);
1543 prepare_action(VectorOperation::insert);
1550 DEAL_II_NAMESPACE_CLOSE
1553 #endif // DEAL_II_WITH_PETSC types::global_dof_index size_type
bool operator==(const const_iterator &) const
void add(const size_type i, const size_type j, const PetscScalar value)
const_iterator begin() const
const Accessor & operator*() const
::ExceptionBase & ExcMessage(std::string arg1)
const Accessor * operator->() const
std_cxx11::shared_ptr< const std::vector< PetscScalar > > value_cache
size_type row_length(const size_type row) const
#define AssertThrow(cond, exc)
std::vector< PetscInt > column_indices
const_iterator(const MatrixBase *matrix, const size_type row, const size_type index)
types::global_dof_index size_type
std::vector< PetscScalar > column_values
MatrixIterators::const_iterator const_iterator
#define DeclException1(Exception1, type1, outsequence)
DeclException0(ExcBeyondEndOfMatrix)
unsigned int global_dof_index
Accessor(const MatrixBase *matrix, const size_type row, const size_type index)
bool operator<(const const_iterator &) const
#define Assert(cond, exc)
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
friend class const_iterator
VectorOperation::values last_action
void prepare_action(const VectorOperation::values new_action)
const_iterator & operator++()
bool in_local_range(const size_type index) const
PetscScalar value() 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)
void assert_is_compressed()