16 #include <deal.II/lac/petsc_matrix_base.h> 18 #ifdef DEAL_II_WITH_PETSC 20 # include <deal.II/lac/exceptions.h> 21 # include <deal.II/lac/petsc_compatibility.h> 22 # include <deal.II/lac/petsc_full_matrix.h> 23 # include <deal.II/lac/petsc_sparse_matrix.h> 24 # include <deal.II/lac/petsc_parallel_sparse_matrix.h> 25 # include <deal.II/lac/petsc_vector.h> 27 DEAL_II_NAMESPACE_OPEN
31 namespace MatrixIterators
41 if (this->a_row == matrix->m())
43 colnum_cache.reset ();
51 const PetscInt *colnums;
52 const PetscScalar *values;
54 PetscErrorCode ierr = MatGetRow(*matrix, this->a_row, &ncols, &colnums,
65 colnum_cache.reset (
new std::vector<size_type> (colnums, colnums+ncols));
66 value_cache.reset (
new std::vector<PetscScalar> (values, values+ncols));
69 ierr = MatRestoreRow(*matrix, this->a_row, &ncols, &colnums, &values);
102 const int m=0,
n=0, n_nonzero_per_row=0;
103 const PetscErrorCode ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,
m,
n,
119 const PetscErrorCode ierr = MatZeroEntries (
matrix);
129 const PetscScalar new_diag_value)
131 std::vector<size_type> rows (1, row);
139 const PetscScalar new_diag_value)
145 const std::vector<PetscInt> petsc_rows (rows.begin(), rows.end());
152 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 154 &petsc_rows[0], &index_set);
157 &petsc_rows[0], PETSC_COPY_VALUES, &index_set);
160 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 161 const PetscErrorCode ierr = MatZeroRowsIS(
matrix, index_set, new_diag_value);
163 const PetscErrorCode ierr = MatZeroRowsIS(
matrix, index_set, new_diag_value,
164 PETSC_NULL, PETSC_NULL);
168 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 169 ISDestroy (index_set);
171 ISDestroy (&index_set);
181 PetscInt petsc_i = i, petsc_j = j;
185 const PetscErrorCode ierr = MatGetValues (
matrix, 1, &petsc_i, 1, &petsc_j,
211 #ifdef DEAL_II_WITH_MPI 215 int all_int_last_action;
217 const int ierr = MPI_Allreduce
218 (&my_int_last_action, &all_int_last_action, 1, MPI_INT, MPI_BOR,
223 ExcMessage(
"Error: not all processors agree on the last " 224 "VectorOperation before this compress() call."));
231 ExcMessage(
"Missing compress() or calling with wrong VectorOperation argument."));
234 PetscErrorCode ierr = MatAssemblyBegin (
matrix,MAT_FINAL_ASSEMBLY);
237 ierr = MatAssemblyEnd (
matrix,MAT_FINAL_ASSEMBLY);
248 PetscInt n_rows, n_cols;
250 const PetscErrorCode ierr = MatGetSize (
matrix, &n_rows, &n_cols);
261 PetscInt n_rows, n_cols;
263 const PetscErrorCode ierr = MatGetSize (
matrix, &n_rows, &n_cols);
274 PetscInt n_rows, n_cols;
276 const PetscErrorCode ierr = MatGetLocalSize (
matrix, &n_rows, &n_cols);
284 std::pair<MatrixBase::size_type, MatrixBase::size_type>
289 const PetscErrorCode ierr = MatGetOwnershipRange (static_cast<const Mat &>(
matrix),
302 const PetscErrorCode ierr = MatGetInfo (
matrix, MAT_GLOBAL_SUM, &mat_info);
305 return static_cast<size_type>(mat_info.nz_used);
325 const PetscInt *colnums;
326 const PetscScalar *values;
330 PetscErrorCode ierr = MatGetRow(*
this, row, &ncols, &colnums, &values);
338 const PetscInt ncols_saved = ncols;
339 ierr = MatRestoreRow(*
this, row, &ncols, &colnums, &values);
351 const PetscErrorCode ierr = MatNorm (
matrix, NORM_1, &result);
364 const PetscErrorCode ierr = MatNorm (
matrix, NORM_INFINITY, &result);
377 const PetscErrorCode ierr = MatNorm (
matrix, NORM_FROBENIUS, &result);
403 #if DEAL_II_PETSC_VERSION_GTE(3,1,0) 405 MatrixBase::trace ()
const 409 const PetscErrorCode ierr = MatGetTrace (
matrix, &result);
421 const PetscErrorCode ierr = MatScale (
matrix, a);
432 const PetscScalar factor = 1./a;
433 const PetscErrorCode ierr = MatScale (
matrix, factor);
444 const PetscErrorCode ierr = MatAXPY (
matrix, factor,
445 other, DIFFERENT_NONZERO_PATTERN);
455 const PetscScalar factor)
457 return add(factor, other);
467 const PetscErrorCode ierr = MatMult (
matrix, src, dst);
479 const PetscErrorCode ierr = MatMultTranspose (
matrix, src, dst);
491 const PetscErrorCode ierr = MatMultAdd (
matrix, src, dst, dst);
503 const PetscErrorCode ierr = MatMultTransposeAdd (
matrix, src, dst, dst);
525 MatrixBase::operator Mat ()
const 533 const PetscErrorCode ierr = MatTranspose(
matrix, MAT_REUSE_MATRIX, &
matrix);
540 PetscBooleanType truth;
542 const PetscErrorCode ierr = MatIsSymmetric (
matrix, tolerance, &truth);
550 PetscBooleanType truth;
553 const PetscErrorCode ierr = MatIsHermitian (
matrix, tolerance, &truth);
565 PetscErrorCode ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
570 ierr = MatView (
matrix, PETSC_VIEWER_STDOUT_WORLD);
578 std::pair<MatrixBase::size_type, MatrixBase::size_type>
582 const PetscInt *colnums;
583 const PetscScalar *values;
586 for (row = loc_range.first; row < loc_range.second; ++row)
588 PetscErrorCode ierr = MatGetRow(*
this, row, &ncols, &colnums, &values);
591 for (PetscInt col = 0; col < ncols; ++col)
593 out <<
"(" << row <<
"," << colnums[col] <<
") " << values[col] << std::endl;
596 ierr = MatRestoreRow(*
this, row, &ncols, &colnums, &values);
609 const PetscErrorCode ierr = MatGetInfo(
matrix, MAT_LOCAL, &info);
612 return sizeof(*this) +
static_cast<size_type>(info.memory);
617 DEAL_II_NAMESPACE_CLOSE
619 #endif // DEAL_II_WITH_PETSC PetscErrorCode destroy_matrix(Mat &matrix)
PetscScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
void add(const size_type i, const size_type j, const PetscScalar value)
static ::ExceptionBase & ExcIO()
void vmult(VectorBase &dst, const VectorBase &src) const
const_iterator begin() const
real_type l2_norm() const
PetscReal frobenius_norm() const
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
size_type row_length(const size_type row) const
#define AssertThrow(cond, exc)
void clear_rows(const std::vector< size_type > &rows, const PetscScalar new_diag_value=0)
PetscReal linfty_norm() const
void compress(const VectorOperation::values operation)
size_type n_nonzero_elements() const
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type local_size() const
#define Assert(cond, exc)
MatrixBase & operator=(const value_type d)
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
const_iterator end() const
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
PetscBooleanType is_symmetric(const double tolerance=1.e-12)
PetscBooleanType is_hermitian(const double tolerance=1.e-12)
#define AssertThrowMPI(error_code)
MatrixBase & operator*=(const PetscScalar factor)
std::pair< size_type, size_type > local_range() const
PetscScalar matrix_norm_square(const VectorBase &v) const
MatrixBase & operator/=(const PetscScalar factor)
PetscReal l1_norm() const
static ::ExceptionBase & ExcSourceEqualsDestination()
virtual const MPI_Comm & get_mpi_communicator() const =0
void clear_row(const size_type row, const PetscScalar new_diag_value=0)
std::size_t memory_consumption() const
void assert_is_compressed()
static ::ExceptionBase & ExcInternalError()