16#ifndef dealii_petsc_block_vector_h
17#define dealii_petsc_block_vector_h
22#ifdef DEAL_II_WITH_PETSC
114 BlockVector(
const std::vector<size_type> &block_sizes,
116 const std::vector<size_type> &local_elements);
122 explicit BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
123 const MPI_Comm &communicator = MPI_COMM_WORLD);
128 BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
129 const std::vector<IndexSet> &ghost_indices,
166 const bool omit_zeroing_entries =
false);
189 reinit(
const std::vector<size_type> &block_sizes,
191 const std::vector<size_type> &locally_owned_sizes,
192 const bool omit_zeroing_entries =
false);
216 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
223 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
224 const std::vector<IndexSet> &ghost_entries,
234 reinit(
const unsigned int num_blocks);
273 print(std::ostream & out,
274 const unsigned int precision = 3,
275 const bool scientific =
true,
276 const bool across =
true)
const;
303 const std::vector<size_type> &block_sizes,
305 const std::vector<size_type> &local_elements)
307 reinit(block_sizes, communicator, local_elements,
false);
317 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
322 const std::vector<IndexSet> ¶llel_partitioning,
325 reinit(parallel_partitioning, communicator);
329 const std::vector<IndexSet> ¶llel_partitioning,
330 const std::vector<IndexSet> &ghost_indices,
333 reinit(parallel_partitioning, ghost_indices, communicator);
369 const bool omit_zeroing_entries)
374 omit_zeroing_entries);
382 const std::vector<size_type> &locally_owned_sizes,
383 const bool omit_zeroing_entries)
386 if (this->
components.size() != this->n_blocks())
389 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
392 locally_owned_sizes[i],
393 omit_zeroing_entries);
401 if (this->
components.size() != this->n_blocks())
404 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
412 std::vector<size_type> sizes(parallel_partitioning.size());
413 for (
unsigned int i = 0; i < parallel_partitioning.size(); ++i)
414 sizes[i] = parallel_partitioning[i].
size();
417 if (this->
components.size() != this->n_blocks())
420 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
421 block(i).
reinit(parallel_partitioning[i], communicator);
426 const std::vector<IndexSet> &ghost_entries,
429 std::vector<types::global_dof_index> sizes(parallel_partitioning.size());
430 for (
unsigned int i = 0; i < parallel_partitioning.size(); ++i)
431 sizes[i] = parallel_partitioning[i].
size();
434 if (this->
components.size() != this->n_blocks())
437 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
448 return block(0).get_mpi_communicator();
456 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
475 const unsigned int precision,
476 const bool scientific,
477 const bool across)
const
479 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
482 out <<
'C' << i <<
':';
484 out <<
"Component " << i << std::endl;
485 this->
components[i].print(out, precision, scientific, across);
510 namespace LinearOperatorImplementation
523 template <
typename Matrix>
529 v.
reinit(matrix.locally_owned_range_indices(),
530 matrix.get_mpi_communicator());
533 template <
typename Matrix>
539 v.
reinit(matrix.locally_owned_domain_indices(),
540 matrix.get_mpi_communicator());
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcNonMatchingBlockVectors()
::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > iterator
unsigned int n_blocks() const
const value_type * const_pointer
typename BlockType::const_reference const_reference
types::global_dof_index size_type
typename BlockType::value_type value_type
BlockVectorBase & operator=(const value_type s)
#define Assert(cond, exc)
std::vector< Vector > components
::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > const_iterator
typename BlockType::reference reference
BlockIndices block_indices
static ::ExceptionBase & ExcInternalError()
BlockType & block(const unsigned int i)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::size_t locally_owned_size() const
const BlockIndices & get_block_indices() const
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
BaseClass::value_type value_type
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
BaseClass::size_type size_type
static void reinit_range_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
bool has_ghost_elements() const
void reinit(const unsigned int n_blocks, const MPI_Comm &communicator, const size_type block_size, const size_type locally_owned_size, const bool omit_zeroing_entries=false)
void swap(BlockVector &v)
BaseClass::const_reference const_reference
void swap(BlockVector &u, BlockVector &v)
BaseClass::pointer pointer
~BlockVector() override=default
BlockVector & operator=(const value_type s)
BaseClass::BlockType BlockType
BaseClass::reference reference
BaseClass::const_pointer const_pointer
const MPI_Comm & get_mpi_communicator() const
static void reinit_domain_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
bool has_ghost_elements() const
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)