16 #ifndef dealii_petsc_block_vector_h 17 #define dealii_petsc_block_vector_h 22 #ifdef DEAL_II_WITH_PETSC 98 const MPI_Comm & communicator,
115 BlockVector(
const std::vector<size_type> &block_sizes,
116 const MPI_Comm & communicator,
117 const std::vector<size_type> &local_elements);
123 explicit BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
124 const MPI_Comm &communicator = MPI_COMM_WORLD);
129 BlockVector(
const std::vector<IndexSet> ¶llel_partitioning,
130 const std::vector<IndexSet> &ghost_indices,
131 const MPI_Comm & communicator);
164 const MPI_Comm & communicator,
167 const bool omit_zeroing_entries =
false);
190 reinit(
const std::vector<size_type> &block_sizes,
191 const MPI_Comm & communicator,
192 const std::vector<size_type> &local_sizes,
193 const bool omit_zeroing_entries =
false);
217 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
218 const MPI_Comm & communicator);
224 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
225 const std::vector<IndexSet> &ghost_entries,
226 const MPI_Comm & communicator);
235 reinit(
const unsigned int num_blocks);
274 print(std::ostream & out,
275 const unsigned int precision = 3,
276 const bool scientific =
true,
277 const bool across =
true)
const;
294 const MPI_Comm & communicator,
298 reinit(n_blocks, communicator, block_size, local_size);
304 const std::vector<size_type> &block_sizes,
305 const MPI_Comm & communicator,
306 const std::vector<size_type> &local_elements)
308 reinit(block_sizes, communicator, local_elements,
false);
318 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
323 const std::vector<IndexSet> ¶llel_partitioning,
324 const MPI_Comm & communicator)
326 reinit(parallel_partitioning, communicator);
330 const std::vector<IndexSet> ¶llel_partitioning,
331 const std::vector<IndexSet> &ghost_indices,
332 const MPI_Comm & communicator)
334 reinit(parallel_partitioning, ghost_indices, communicator);
367 const MPI_Comm & communicator,
370 const bool omit_zeroing_entries)
372 reinit(std::vector<size_type>(n_blocks, block_size),
374 std::vector<size_type>(n_blocks, local_size),
375 omit_zeroing_entries);
382 const MPI_Comm & communicator,
383 const std::vector<size_type> &local_sizes,
384 const bool omit_zeroing_entries)
390 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
394 omit_zeroing_entries);
405 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
406 block(i).reinit(v.
block(i), omit_zeroing_entries);
411 const MPI_Comm & communicator)
413 std::vector<size_type> sizes(parallel_partitioning.size());
414 for (
unsigned int i = 0; i < parallel_partitioning.size(); ++i)
415 sizes[i] = parallel_partitioning[i].
size();
421 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
422 block(i).reinit(parallel_partitioning[i], communicator);
427 const std::vector<IndexSet> &ghost_entries,
428 const MPI_Comm & communicator)
430 std::vector<types::global_dof_index> sizes(parallel_partitioning.size());
431 for (
unsigned int i = 0; i < parallel_partitioning.size(); ++i)
432 sizes[i] = parallel_partitioning[i].
size();
438 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
439 block(i).reinit(parallel_partitioning[i],
446 inline const MPI_Comm &
449 return block(0).get_mpi_communicator();
455 bool ghosted =
block(0).has_ghost_elements();
457 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
476 const unsigned int precision,
477 const bool scientific,
478 const bool across)
const 480 for (
unsigned int i = 0; i < this->
n_blocks(); ++i)
483 out <<
'C' << i <<
':';
485 out <<
"Component " << i << std::endl;
486 this->
components[i].print(out, precision, scientific, across);
512 namespace LinearOperatorImplementation
525 template <
typename Matrix>
531 v.
reinit(matrix.locally_owned_range_indices(),
532 matrix.get_mpi_communicator());
535 template <
typename Matrix>
541 v.
reinit(matrix.locally_owned_domain_indices(),
542 matrix.get_mpi_communicator());
562 #endif // DEAL_II_WITH_PETSC void reinit(const unsigned int n_blocks, const MPI_Comm &communicator, const size_type block_size, const size_type local_size, const bool omit_zeroing_entries=false)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
void swap(BlockVector &v)
Contents is actually a matrix.
BaseClass::reference reference
const MPI_Comm & get_mpi_communicator() const
BaseClass::size_type size_type
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
BaseClass::const_reference const_reference
types::global_dof_index size_type
~BlockVector() override=default
::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > const_iterator
static ::ExceptionBase & ExcNonMatchingBlockVectors()
BlockIndices block_indices
void swap(BlockVector &u, BlockVector &v)
BaseClass::pointer pointer
::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > iterator
BlockVector< double > BlockVector
const BlockIndices & get_block_indices() const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
#define DeclException0(Exception0)
#define DEAL_II_NAMESPACE_CLOSE
static void reinit_range_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
std::vector< Vector > components
BaseClass::BlockType BlockType
unsigned int n_blocks() const
static void reinit_domain_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
typename BlockType::value_type value_type
#define DEAL_II_NAMESPACE_OPEN
typename BlockType::const_reference const_reference
bool has_ghost_elements() const
BaseClass::const_pointer const_pointer
BlockVector & operator=(const value_type s)
BaseClass::value_type value_type
BlockType & block(const unsigned int i)
BlockVectorBase & operator=(const value_type s)
typename BlockType::reference reference
const value_type * const_pointer
static ::ExceptionBase & ExcInternalError()