16 #ifndef dealii__petsc_parallel_block_vector_h 17 #define dealii__petsc_parallel_block_vector_h 20 #include <deal.II/base/config.h> 22 #ifdef DEAL_II_WITH_PETSC 24 # include <deal.II/lac/petsc_parallel_vector.h> 25 # include <deal.II/lac/block_indices.h> 26 # include <deal.II/lac/block_vector_base.h> 27 # include <deal.II/lac/exceptions.h> 28 # include <deal.II/lac/vector_type_traits.h> 30 DEAL_II_NAMESPACE_OPEN
79 typedef BaseClass::pointer pointer;
80 typedef BaseClass::const_pointer const_pointer;
81 typedef BaseClass::reference reference;
82 typedef BaseClass::const_reference const_reference;
83 typedef BaseClass::size_type size_type;
99 const MPI_Comm &communicator,
100 const size_type block_size,
101 const size_type local_size);
116 BlockVector (
const std::vector<size_type> &block_sizes,
117 const MPI_Comm &communicator,
118 const std::vector<size_type> &local_elements);
124 explicit BlockVector (
const std::vector<IndexSet> ¶llel_partitioning,
125 const MPI_Comm &communicator = MPI_COMM_WORLD);
130 BlockVector (
const std::vector<IndexSet> ¶llel_partitioning,
131 const std::vector<IndexSet> &ghost_indices,
132 const MPI_Comm &communicator);
182 const MPI_Comm &communicator,
183 const size_type block_size,
184 const size_type local_size,
185 const bool omit_zeroing_entries =
false);
207 void reinit (
const std::vector<size_type> &block_sizes,
208 const MPI_Comm &communicator,
209 const std::vector<size_type> &local_sizes,
210 const bool omit_zeroing_entries=
false);
227 const bool omit_zeroing_entries=
false);
233 void reinit (
const std::vector<IndexSet> ¶llel_partitioning,
234 const MPI_Comm &communicator);
239 void reinit (
const std::vector<IndexSet> ¶llel_partitioning,
240 const std::vector<IndexSet> &ghost_entries,
241 const MPI_Comm &communicator);
249 void reinit (
const unsigned int num_blocks);
284 void print (std::ostream &out,
285 const unsigned int precision = 3,
286 const bool scientific =
true,
287 const bool across =
true)
const;
312 const MPI_Comm &communicator,
313 const size_type block_size,
314 const size_type local_size)
323 const MPI_Comm &communicator,
324 const std::vector<size_type> &local_elements)
326 reinit (block_sizes, communicator, local_elements,
false);
338 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
344 const MPI_Comm &communicator)
346 reinit(parallel_partitioning, communicator);
351 const std::vector<IndexSet> &ghost_indices,
352 const MPI_Comm &communicator)
354 reinit(parallel_partitioning, ghost_indices, communicator);
377 for (size_type i=0; i<this->
n_blocks(); ++i)
393 const MPI_Comm &communicator,
394 const size_type block_size,
395 const size_type local_size,
396 const bool omit_zeroing_entries)
400 std::vector<size_type>(
n_blocks, local_size),
401 omit_zeroing_entries);
409 const MPI_Comm &communicator,
410 const std::vector<size_type> &local_sizes,
411 const bool omit_zeroing_entries)
417 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
419 local_sizes[i], omit_zeroing_entries);
426 const bool omit_zeroing_entries)
432 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
439 const MPI_Comm &communicator)
441 std::vector<size_type> sizes(parallel_partitioning.size());
442 for (
unsigned int i=0; i<parallel_partitioning.size(); ++i)
443 sizes[i] = parallel_partitioning[i].
size();
449 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
450 block(i).
reinit(parallel_partitioning[i], communicator);
456 const std::vector<IndexSet> &ghost_entries,
457 const MPI_Comm &communicator)
459 std::vector<types::global_dof_index> sizes(parallel_partitioning.size());
460 for (
unsigned int i=0; i<parallel_partitioning.size(); ++i)
461 sizes[i] = parallel_partitioning[i].
size();
467 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
468 block(i).
reinit(parallel_partitioning[i], ghost_entries[i], communicator);
477 return block(0).get_mpi_communicator();
484 bool ghosted=
block(0).has_ghost_elements();
486 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
507 const unsigned int precision,
508 const bool scientific,
509 const bool across)
const 511 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
514 out <<
'C' << i <<
':';
516 out <<
"Component " << i << std::endl;
517 this->
components[i].print(out, precision, scientific, across);
546 template <
typename>
class ReinitHelper;
556 template <
typename Matrix>
562 v.
reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator());
565 template <
typename Matrix>
571 v.
reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator());
590 DEAL_II_NAMESPACE_CLOSE
592 #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)
const MPI_Comm & get_mpi_communicator() const
void swap(BlockVector &u, BlockVector &v)
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
static ::ExceptionBase & ExcNonMatchingBlockVectors()
BaseClass::value_type value_type
BlockIndices block_indices
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)
std::vector< Vector > components
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
BaseClass::BlockType BlockType
unsigned int n_blocks() const
BlockVectorBase< Vector > BaseClass
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
bool has_ghost_elements() const
BlockVector & operator=(const value_type s)
BlockType & block(const unsigned int i)
BlockVectorBase & operator=(const value_type s)
static ::ExceptionBase & ExcInternalError()