15#ifndef dealii_trilinos_parallel_block_vector_h
16#define dealii_trilinos_parallel_block_vector_h
21#ifdef DEAL_II_WITH_TRILINOS
34template <
typename Number>
110 explicit BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
111 const MPI_Comm communicator = MPI_COMM_WORLD);
118 BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
119 const std::vector<IndexSet> &ghost_values,
121 const bool vector_writable =
false);
184 template <
typename Number>
186 operator=(const ::BlockVector<Number> &v);
197 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
198 const MPI_Comm communicator = MPI_COMM_WORLD,
199 const bool omit_zeroing_entries =
false);
219 reinit(
const std::vector<IndexSet> &partitioning,
220 const std::vector<IndexSet> &ghost_values,
221 const MPI_Comm communicator = MPI_COMM_WORLD,
222 const bool vector_writable =
false);
236 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
238 const bool make_ghosted =
true,
239 const bool vector_writable =
false);
321 print(std::ostream &out,
322 const unsigned int precision = 3,
323 const bool scientific =
true,
324 const bool across =
true)
const;
341 const std::vector<IndexSet> ¶llel_partitioning,
344 reinit(parallel_partitioning, communicator,
false);
350 const std::vector<IndexSet> ¶llel_partitioning,
351 const std::vector<IndexSet> &ghost_values,
353 const bool vector_writable)
355 reinit(parallel_partitioning,
376 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
393 template <
typename Number>
402 if (this->
n_blocks() != v.n_blocks())
406 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
419 bool ghosted =
block(0).has_ghost_elements();
421 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
432 std::swap(this->components, v.components);
434 ::swap(this->block_indices, v.block_indices);
461 namespace LinearOperatorImplementation
474 template <
typename Matrix>
478 bool omit_zeroing_entries)
480 v.
reinit(matrix.locally_owned_range_indices(),
481 matrix.get_mpi_communicator(),
482 omit_zeroing_entries);
485 template <
typename Matrix>
489 bool omit_zeroing_entries)
491 v.
reinit(matrix.locally_owned_domain_indices(),
492 matrix.get_mpi_communicator(),
493 omit_zeroing_entries);
::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
std::vector< MPI::Vector > components
typename BlockType::reference reference
BlockIndices block_indices
BlockType & block(const unsigned int i)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
bool has_ghost_elements() const
BaseClass::BlockType BlockType
BaseClass::const_pointer const_pointer
BaseClass::size_type size_type
void swap(BlockVector &u, BlockVector &v) noexcept
BaseClass::const_reference const_reference
void swap(BlockVector &v) noexcept
BlockVector & operator=(const value_type s)
void reinit(const std::vector< IndexSet > ¶llel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
void import_nonlocal_data_for_fe(const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v)
BaseClass::pointer pointer
BaseClass::value_type value_type
BaseClass::reference reference
~BlockVector() override=default
static void reinit_range_vector(const Matrix &matrix, TrilinosWrappers::MPI::BlockVector &v, bool omit_zeroing_entries)
static void reinit_domain_vector(const Matrix &matrix, TrilinosWrappers::MPI::BlockVector &v, bool omit_zeroing_entries)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
#define Assert(cond, exc)
static ::ExceptionBase & ExcNonMatchingBlockVectors()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
void swap(BlockVector &u, BlockVector &v) noexcept