15#ifndef dealii_petsc_block_vector_h
16#define dealii_petsc_block_vector_h
21#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,
143 template <
size_t num_blocks>
144 explicit BlockVector(
const std::array<Vec, num_blocks> &);
193 const bool omit_zeroing_entries =
false);
216 reinit(
const std::vector<size_type> &block_sizes,
218 const std::vector<size_type> &locally_owned_sizes,
219 const bool omit_zeroing_entries =
false);
243 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
250 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
251 const std::vector<IndexSet> &ghost_entries,
263 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
265 const bool make_ghosted =
true);
294 reinit(
const unsigned int num_blocks);
315 operator const Vec &()
const;
349 print(std::ostream &out,
350 const unsigned int precision = 3,
351 const bool scientific =
true,
352 const bool across =
true)
const;
385 , petsc_nest_vector(nullptr)
402 const std::vector<size_type> &block_sizes,
404 const std::vector<size_type> &local_elements)
407 reinit(block_sizes, communicator, local_elements,
false);
418 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
427 const std::vector<IndexSet> ¶llel_partitioning,
431 reinit(parallel_partitioning, communicator);
437 const std::vector<IndexSet> ¶llel_partitioning,
438 const std::vector<IndexSet> &ghost_indices,
442 reinit(parallel_partitioning, ghost_indices, communicator);
455 template <
size_t num_blocks>
462 for (
auto i = 0; i < num_blocks; ++i)
490 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
505 const bool omit_zeroing_entries)
510 omit_zeroing_entries);
518 const std::vector<size_type> &locally_owned_sizes,
519 const bool omit_zeroing_entries)
524 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
527 locally_owned_sizes[i],
528 omit_zeroing_entries);
541 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
558 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
569 const std::vector<IndexSet> &ghost_entries,
579 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
592 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
594 const bool make_ghosted)
601 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
623 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
644 const unsigned int precision,
645 const bool scientific,
646 const bool across)
const
648 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
651 out <<
'C' << i <<
':';
653 out <<
"Component " << i << std::endl;
654 this->
components[i].print(out, precision, scientific, across);
678 namespace LinearOperatorImplementation
691 template <
typename Matrix>
697 v.
reinit(matrix.locally_owned_range_indices(),
698 matrix.get_mpi_communicator());
701 template <
typename Matrix>
707 v.
reinit(matrix.locally_owned_domain_indices(),
708 matrix.get_mpi_communicator());
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > iterator
unsigned int n_blocks() const
::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > const_iterator
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)
std::vector< Vector > components
typename BlockType::reference reference
BlockIndices block_indices
BlockType & block(const unsigned int i)
std::size_t locally_owned_size() const
const BlockIndices & get_block_indices() const
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
bool has_ghost_elements() const
void swap(BlockVector &v)
BaseClass::const_reference const_reference
void swap(BlockVector &u, BlockVector &v)
MPI_Comm get_mpi_communicator() const
BaseClass::pointer pointer
BlockVector & operator=(const value_type s)
BaseClass::reference reference
void compress(VectorOperation::values operation)
BaseClass::const_pointer const_pointer
bool has_ghost_elements() const
static void reinit_range_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
static void reinit_domain_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcNonMatchingBlockVectors()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()