16 #include <deal.II/lac/petsc_sparse_matrix.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_vector_base.h> 23 # include <deal.II/lac/sparsity_pattern.h> 24 # include <deal.II/lac/dynamic_sparsity_pattern.h> 26 DEAL_II_NAMESPACE_OPEN
33 const int m=0,
n=0, n_nonzero_per_row=0;
34 const PetscErrorCode ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,
m,
n,
35 n_nonzero_per_row,
nullptr, &
matrix);
44 const bool is_symmetric)
53 const std::vector<size_type> &row_lengths,
54 const bool is_symmetric)
61 template <
typename SparsityPatternType>
64 const bool preset_nonzero_locations)
66 do_reinit (sparsity_pattern, preset_nonzero_locations);
84 const bool is_symmetric)
99 const std::vector<size_type> &row_lengths,
100 const bool is_symmetric)
112 template <
typename SparsityPatternType>
115 reinit (
const SparsityPatternType &sparsity_pattern,
116 const bool preset_nonzero_locations)
123 do_reinit (sparsity_pattern, preset_nonzero_locations);
131 static MPI_Comm comm;
132 const PetscErrorCode ierr = PetscObjectGetComm((PetscObject)
matrix, &comm);
143 const bool is_symmetric)
148 const PetscErrorCode ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,
m,
n,
165 const std::vector<size_type> &row_lengths,
166 const bool is_symmetric)
168 Assert (row_lengths.size() ==
m,
179 const std::vector<PetscInt>
180 int_row_lengths (row_lengths.begin(), row_lengths.end());
182 const PetscErrorCode ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,
m,
n, 0,
183 int_row_lengths.data(), &
matrix);
195 template <
typename SparsityPatternType>
198 const bool preset_nonzero_locations)
200 std::vector<size_type> row_lengths (sparsity_pattern.n_rows());
201 for (
size_type i=0; i<sparsity_pattern.n_rows(); ++i)
202 row_lengths[i] = sparsity_pattern.row_length (i);
205 sparsity_pattern.n_cols(),
221 if (preset_nonzero_locations ==
true)
223 std::vector<PetscInt> row_entries;
224 std::vector<PetscScalar> row_values;
225 for (
size_type i=0; i<sparsity_pattern.n_rows(); ++i)
227 row_entries.resize (row_lengths[i]);
228 row_values.resize (row_lengths[i], 0.0);
229 for (
size_type j=0; j<row_lengths[i]; ++j)
230 row_entries[j] = sparsity_pattern.column_number (i,j);
232 const PetscInt int_row = i;
233 const PetscErrorCode ierr = MatSetValues (
matrix, 1, &int_row,
234 row_lengths[i], row_entries.data(),
235 row_values.data(), INSERT_VALUES);
249 const PetscErrorCode ierr = MatGetSize(
matrix, &
m, &
n);
259 const PetscErrorCode ierr = MatGetSize(
matrix, &
m, &
n);
312 DEAL_II_NAMESPACE_CLOSE
314 #endif // DEAL_II_WITH_PETSC PetscErrorCode destroy_matrix(Mat &matrix)
void Tmmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
#define AssertThrow(cond, exc)
void set_matrix_option(Mat &matrix, const MatOption option_name, const PetscBool option_value=PETSC_FALSE)
PetscBool is_symmetric(const double tolerance=1.e-12)
void compress(const VectorOperation::values operation)
MatrixBase & operator=(const MatrixBase &)=delete
SparseMatrix & operator=(const double d)
void reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
virtual const MPI_Comm & get_mpi_communicator() const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void mmult(MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
types::global_dof_index size_type
void set_keep_zero_rows(Mat &matrix)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void do_reinit(const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false)
void close_matrix(Mat &matrix)