19#ifdef DEAL_II_WITH_PETSC
35 ierr = MatCreate(
comm, &dummy);
37 ierr = MatSetSizes(dummy, lr, lc, gr, gc);
39 ierr = MatSetType(dummy, MATAIJ);
41 ierr = MatSeqAIJSetPreallocation(dummy, 0,
nullptr);
43 ierr = MatMPIAIJSetPreallocation(dummy, 0,
nullptr, 0,
nullptr);
45 ierr = MatSetUp(dummy);
47 ierr = MatSetOption(dummy, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);
49 ierr = MatAssemblyBegin(dummy, MAT_FINAL_ASSEMBLY);
51 ierr = MatAssemblyEnd(dummy, MAT_FINAL_ASSEMBLY);
94 this->
sub_objects.reinit(n_block_rows, n_block_columns);
112 const std::vector<IndexSet> & cols,
123 std::vector<types::global_dof_index> row_sizes;
128 std::vector<types::global_dof_index> col_sizes;
142 p->
reinit(rows[r], cols[c], bdsp.
block(r, c), com);
154 reinit(sizes, sizes, bdsp, com);
169 std::vector<size_type> row_sizes(
m,
size_type(-1));
170 std::vector<size_type> col_sizes(
n,
size_type(-1));
171 std::vector<size_type> row_local_sizes(
m,
size_type(-1));
172 std::vector<size_type> col_local_sizes(
n,
size_type(-1));
183 row_local_sizes[r] = this->
sub_objects[r][c]->local_size();
198 "When passing empty sub-blocks of a block matrix, you need to make "
199 "sure that at least one block in each block row and block column is "
200 "non-empty. However, block row " +
202 " is completely empty "
203 "and so it is not possible to determine how many rows it should have."));
207 "When passing empty sub-blocks of a block matrix, you need to make "
208 "sure that at least one block in each block row and block column is "
209 "non-empty. However, block column " +
211 " is completely empty "
212 "and so it is not possible to determine how many columns it should have."));
213 Mat dummy = ::create_dummy_mat(
215 static_cast<PetscInt
>(row_local_sizes[r]),
216 static_cast<PetscInt
>(row_sizes[r]),
217 static_cast<PetscInt
>(col_local_sizes[c]),
218 static_cast<PetscInt
>(col_sizes[c]));
223 ierr = MatDestroy(&dummy);
250 std::vector<Mat> psub_objects(
m *
n);
251 for (
unsigned int r = 0; r <
m; r++)
252 for (
unsigned int c = 0; c <
n; c++)
255 psub_objects[r *
n + c] = this->
sub_objects[r][c]->petsc_matrix();
257 ierr = MatCreateNest(
276 std::vector<IndexSet>
279 std::vector<IndexSet> index_sets;
282 index_sets.push_back(this->block(0, i).locally_owned_domain_indices());
289 std::vector<IndexSet>
292 std::vector<IndexSet> index_sets;
295 index_sets.push_back(this->block(i, 0).locally_owned_range_indices());
305 std::uint64_t n_nonzero = 0;
321 BlockSparseMatrix::operator
const Mat &()
const
323 return petsc_nest_matrix;
340 PetscInt nr = 1, nc = 1;
342 PetscErrorCode ierr =
343 PetscObjectTypeCompare(
reinterpret_cast<PetscObject
>(A),
347 std::vector<Mat> mats;
348 bool need_empty_matrices =
false;
351 ierr = MatNestGetSize(A, &nr, &nc);
353 for (PetscInt i = 0; i < nr; ++i)
355 for (PetscInt j = 0; j < nc; ++j)
358 ierr = MatNestGetSubMat(A, i, j, &sA);
361 need_empty_matrices =
true;
370 std::vector<size_type> r_block_sizes(nr, 0);
371 std::vector<size_type> c_block_sizes(nc, 0);
375 for (PetscInt i = 0; i < nr; ++i)
377 for (PetscInt j = 0; j < nc; ++j)
379 if (mats[i * nc + j])
385 if (need_empty_matrices)
389 if (need_empty_matrices || !isnest)
395 ierr = PetscObjectReference(
reinterpret_cast<PetscObject
>(A));
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
BlockIndices column_block_indices
unsigned int n_block_rows() const
void compress(VectorOperation::values operation)
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
unsigned int n_block_cols() const
BlockIndices row_block_indices
BlockType & block(const unsigned int row, const unsigned int column)
size_type n_block_rows() const
SparsityPatternType & block(const size_type row, const size_type column)
size_type n_block_cols() const
~BlockSparseMatrix() override
BaseClass::size_type size_type
void reinit(const size_type n_block_rows, const size_type n_block_columns)
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
MPI_Comm get_mpi_communicator() const
void create_empty_matrices_if_needed()
std::uint64_t n_nonzero_elements() const
BaseClass::BlockType BlockType
std::vector< IndexSet > locally_owned_domain_indices() const
std::vector< IndexSet > locally_owned_range_indices() const
void compress(VectorOperation::values operation)
std::size_t n_nonzero_elements() const
virtual void reinit(const SparsityPattern &sparsity)
Subscriptor & operator=(const Subscriptor &)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertNothrow(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void petsc_increment_state_counter(Vec v)