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> 29 DEAL_II_NAMESPACE_OPEN
78 typedef BaseClass::pointer pointer;
79 typedef BaseClass::const_pointer const_pointer;
80 typedef BaseClass::reference reference;
81 typedef BaseClass::const_reference const_reference;
82 typedef BaseClass::size_type size_type;
98 const MPI_Comm &communicator,
99 const size_type block_size,
100 const size_type local_size);
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);
181 const MPI_Comm &communicator,
182 const size_type block_size,
183 const size_type local_size,
184 const bool omit_zeroing_entries =
false);
206 void reinit (
const std::vector<size_type> &block_sizes,
207 const MPI_Comm &communicator,
208 const std::vector<size_type> &local_sizes,
209 const bool omit_zeroing_entries=
false);
226 const bool omit_zeroing_entries=
false);
232 void reinit (
const std::vector<IndexSet> ¶llel_partitioning,
233 const MPI_Comm &communicator);
238 void reinit (
const std::vector<IndexSet> ¶llel_partitioning,
239 const std::vector<IndexSet> &ghost_entries,
240 const MPI_Comm &communicator);
248 void reinit (
const unsigned int num_blocks);
283 void print (std::ostream &out,
284 const unsigned int precision = 3,
285 const bool scientific =
true,
286 const bool across =
true)
const;
311 const MPI_Comm &communicator,
312 const size_type block_size,
313 const size_type local_size)
315 reinit (n_blocks, communicator, block_size, local_size);
322 const MPI_Comm &communicator,
323 const std::vector<size_type> &local_elements)
325 reinit (block_sizes, communicator, local_elements,
false);
337 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
343 const MPI_Comm &communicator)
345 reinit(parallel_partitioning, communicator);
350 const std::vector<IndexSet> &ghost_indices,
351 const MPI_Comm &communicator)
353 reinit(parallel_partitioning, ghost_indices, communicator);
376 for (size_type i=0; i<this->
n_blocks(); ++i)
392 const MPI_Comm &communicator,
393 const size_type block_size,
394 const size_type local_size,
395 const bool omit_zeroing_entries)
397 reinit(std::vector<size_type>(n_blocks, block_size),
399 std::vector<size_type>(n_blocks, local_size),
400 omit_zeroing_entries);
408 const MPI_Comm &communicator,
409 const std::vector<size_type> &local_sizes,
410 const bool omit_zeroing_entries)
416 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
418 local_sizes[i], omit_zeroing_entries);
425 const bool omit_zeroing_entries)
431 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
438 const MPI_Comm &communicator)
440 std::vector<size_type> sizes(parallel_partitioning.size());
441 for (
unsigned int i=0; i<parallel_partitioning.size(); ++i)
442 sizes[i] = parallel_partitioning[i].
size();
448 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
449 block(i).
reinit(parallel_partitioning[i], communicator);
455 const std::vector<IndexSet> &ghost_entries,
456 const MPI_Comm &communicator)
458 std::vector<types::global_dof_index> sizes(parallel_partitioning.size());
459 for (
unsigned int i=0; i<parallel_partitioning.size(); ++i)
460 sizes[i] = parallel_partitioning[i].
size();
466 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
467 block(i).
reinit(parallel_partitioning[i], ghost_entries[i], communicator);
476 return block(0).get_mpi_communicator();
483 bool ghosted=
block(0).has_ghost_elements();
485 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
486 Assert(
block(i).has_ghost_elements()==ghosted, ExcInternalError());
506 const unsigned int precision,
507 const bool scientific,
508 const bool across)
const 510 for (
unsigned int i=0; i<this->
n_blocks(); ++i)
513 out <<
'C' << i <<
':';
515 out <<
"Component " << i << std::endl;
516 this->
components[i].print(out, precision, scientific, across);
541 DEAL_II_NAMESPACE_CLOSE
543 #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)
DeclException0(ExcIteratorRangeDoesNotMatchVectorSize)
BaseClass::value_type value_type
BlockIndices block_indices
const BlockIndices & get_block_indices() const
#define Assert(cond, exc)
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
std::vector< Vector > components
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)