Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
petsc_block_vector.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_petsc_block_vector_h
17 #define dealii_petsc_block_vector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_PETSC
23 
26 # include <deal.II/lac/exceptions.h>
29 
31 
32 
33 namespace PETScWrappers
34 {
35  // forward declaration
36  class BlockVector;
37 
38  namespace MPI
39  {
61  class BlockVector : public BlockVectorBase<Vector>
62  {
63  public:
68 
73 
85 
89  BlockVector() = default;
90 
97  explicit BlockVector(const unsigned int n_blocks,
98  const MPI_Comm & communicator,
99  const size_type block_size,
100  const size_type local_size);
101 
106  BlockVector(const BlockVector &V);
107 
115  BlockVector(const std::vector<size_type> &block_sizes,
116  const MPI_Comm & communicator,
117  const std::vector<size_type> &local_elements);
118 
123  explicit BlockVector(const std::vector<IndexSet> &parallel_partitioning,
124  const MPI_Comm &communicator = MPI_COMM_WORLD);
125 
129  BlockVector(const std::vector<IndexSet> &parallel_partitioning,
130  const std::vector<IndexSet> &ghost_indices,
131  const MPI_Comm & communicator);
132 
133 
134 
138  ~BlockVector() override = default;
139 
144  BlockVector &
145  operator=(const value_type s);
146 
150  BlockVector &
151  operator=(const BlockVector &V);
152 
162  void
163  reinit(const unsigned int n_blocks,
164  const MPI_Comm & communicator,
165  const size_type block_size,
166  const size_type local_size,
167  const bool omit_zeroing_entries = false);
168 
189  void
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);
194 
209  void
210  reinit(const BlockVector &V, const bool omit_zeroing_entries = false);
211 
216  void
217  reinit(const std::vector<IndexSet> &parallel_partitioning,
218  const MPI_Comm & communicator);
219 
223  void
224  reinit(const std::vector<IndexSet> &parallel_partitioning,
225  const std::vector<IndexSet> &ghost_entries,
226  const MPI_Comm & communicator);
227 
234  void
235  reinit(const unsigned int num_blocks);
236 
240  bool
241  has_ghost_elements() const;
242 
247  const MPI_Comm &
248  get_mpi_communicator() const;
249 
267  void
268  swap(BlockVector &v);
269 
273  void
274  print(std::ostream & out,
275  const unsigned int precision = 3,
276  const bool scientific = true,
277  const bool across = true) const;
278 
287  };
288 
291  /*--------------------- Inline functions --------------------------------*/
292 
293  inline BlockVector::BlockVector(const unsigned int n_blocks,
294  const MPI_Comm & communicator,
295  const size_type block_size,
296  const size_type local_size)
297  {
298  reinit(n_blocks, communicator, block_size, local_size);
299  }
300 
301 
302 
304  const std::vector<size_type> &block_sizes,
305  const MPI_Comm & communicator,
306  const std::vector<size_type> &local_elements)
307  {
308  reinit(block_sizes, communicator, local_elements, false);
309  }
310 
311 
314  {
315  this->components.resize(v.n_blocks());
316  this->block_indices = v.block_indices;
317 
318  for (unsigned int i = 0; i < this->n_blocks(); ++i)
319  this->components[i] = v.components[i];
320  }
321 
323  const std::vector<IndexSet> &parallel_partitioning,
324  const MPI_Comm & communicator)
325  {
326  reinit(parallel_partitioning, communicator);
327  }
328 
330  const std::vector<IndexSet> &parallel_partitioning,
331  const std::vector<IndexSet> &ghost_indices,
332  const MPI_Comm & communicator)
333  {
334  reinit(parallel_partitioning, ghost_indices, communicator);
335  }
336 
337  inline BlockVector &
339  {
341  return *this;
342  }
343 
344  inline BlockVector &
346  {
347  // we only allow assignment to vectors with the same number of blocks
348  // or to an empty BlockVector
349  Assert(n_blocks() == 0 || n_blocks() == v.n_blocks(),
351 
352  if (this->n_blocks() != v.n_blocks())
353  reinit(v.n_blocks());
354 
355  for (size_type i = 0; i < this->n_blocks(); ++i)
356  this->components[i] = v.block(i);
357 
358  collect_sizes();
359 
360  return *this;
361  }
362 
363 
364 
365  inline void
366  BlockVector::reinit(const unsigned int n_blocks,
367  const MPI_Comm & communicator,
368  const size_type block_size,
369  const size_type local_size,
370  const bool omit_zeroing_entries)
371  {
372  reinit(std::vector<size_type>(n_blocks, block_size),
373  communicator,
374  std::vector<size_type>(n_blocks, local_size),
375  omit_zeroing_entries);
376  }
377 
378 
379 
380  inline void
381  BlockVector::reinit(const std::vector<size_type> &block_sizes,
382  const MPI_Comm & communicator,
383  const std::vector<size_type> &local_sizes,
384  const bool omit_zeroing_entries)
385  {
386  this->block_indices.reinit(block_sizes);
387  if (this->components.size() != this->n_blocks())
388  this->components.resize(this->n_blocks());
389 
390  for (unsigned int i = 0; i < this->n_blocks(); ++i)
391  this->components[i].reinit(communicator,
392  block_sizes[i],
393  local_sizes[i],
394  omit_zeroing_entries);
395  }
396 
397 
398  inline void
399  BlockVector::reinit(const BlockVector &v, const bool omit_zeroing_entries)
400  {
401  this->block_indices = v.get_block_indices();
402  if (this->components.size() != this->n_blocks())
403  this->components.resize(this->n_blocks());
404 
405  for (unsigned int i = 0; i < this->n_blocks(); ++i)
406  block(i).reinit(v.block(i), omit_zeroing_entries);
407  }
408 
409  inline void
410  BlockVector::reinit(const std::vector<IndexSet> &parallel_partitioning,
411  const MPI_Comm & communicator)
412  {
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();
416 
417  this->block_indices.reinit(sizes);
418  if (this->components.size() != this->n_blocks())
419  this->components.resize(this->n_blocks());
420 
421  for (unsigned int i = 0; i < this->n_blocks(); ++i)
422  block(i).reinit(parallel_partitioning[i], communicator);
423  }
424 
425  inline void
426  BlockVector::reinit(const std::vector<IndexSet> &parallel_partitioning,
427  const std::vector<IndexSet> &ghost_entries,
428  const MPI_Comm & communicator)
429  {
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();
433 
434  this->block_indices.reinit(sizes);
435  if (this->components.size() != this->n_blocks())
436  this->components.resize(this->n_blocks());
437 
438  for (unsigned int i = 0; i < this->n_blocks(); ++i)
439  block(i).reinit(parallel_partitioning[i],
440  ghost_entries[i],
441  communicator);
442  }
443 
444 
445 
446  inline const MPI_Comm &
448  {
449  return block(0).get_mpi_communicator();
450  }
451 
452  inline bool
454  {
455  bool ghosted = block(0).has_ghost_elements();
456 # ifdef DEBUG
457  for (unsigned int i = 0; i < this->n_blocks(); ++i)
458  Assert(block(i).has_ghost_elements() == ghosted, ExcInternalError());
459 # endif
460  return ghosted;
461  }
462 
463 
464  inline void
466  {
467  std::swap(this->components, v.components);
468 
470  }
471 
472 
473 
474  inline void
475  BlockVector::print(std::ostream & out,
476  const unsigned int precision,
477  const bool scientific,
478  const bool across) const
479  {
480  for (unsigned int i = 0; i < this->n_blocks(); ++i)
481  {
482  if (across)
483  out << 'C' << i << ':';
484  else
485  out << "Component " << i << std::endl;
486  this->components[i].print(out, precision, scientific, across);
487  }
488  }
489 
490 
491 
500  inline void
502  {
503  u.swap(v);
504  }
505 
506  } // namespace MPI
507 
508 } // namespace PETScWrappers
509 
510 namespace internal
511 {
512  namespace LinearOperatorImplementation
513  {
514  template <typename>
515  class ReinitHelper;
516 
521  template <>
523  {
524  public:
525  template <typename Matrix>
526  static void
527  reinit_range_vector(const Matrix & matrix,
529  bool /*omit_zeroing_entries*/)
530  {
531  v.reinit(matrix.locally_owned_range_indices(),
532  matrix.get_mpi_communicator());
533  }
534 
535  template <typename Matrix>
536  static void
537  reinit_domain_vector(const Matrix & matrix,
539  bool /*omit_zeroing_entries*/)
540  {
541  v.reinit(matrix.locally_owned_domain_indices(),
542  matrix.get_mpi_communicator());
543  }
544  };
545 
546  } // namespace LinearOperatorImplementation
547 } /* namespace internal */
548 
549 
555 template <>
557 {};
558 
559 
561 
562 #endif // DEAL_II_WITH_PETSC
563 
564 #endif
PETScWrappers::MPI::BlockVector::~BlockVector
~BlockVector() override=default
BlockVectorBase< Vector >::value_type
typename BlockType::value_type value_type
Definition: block_vector_base.h:468
PETScWrappers::MPI::BlockVector::reinit
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)
Definition: petsc_block_vector.h:366
CUDAWrappers::block_size
constexpr int block_size
Definition: cuda_size.h:29
PETScWrappers
Definition: petsc_block_sparse_matrix.h:36
BlockVector
Definition: block_vector.h:71
block_vector_base.h
BlockVectorBase< Vector >::size
std::size_t size() const
PETScWrappers::MPI::BlockVector::operator=
BlockVector & operator=(const value_type s)
Definition: petsc_block_vector.h:338
LAPACKSupport::V
static const char V
Definition: lapack_support.h:175
PETScWrappers::MPI::swap
void swap(BlockVector &u, BlockVector &v)
Definition: petsc_block_vector.h:501
block_indices.h
MPI_Comm
exceptions.h
internal::LinearOperatorImplementation::ReinitHelper< PETScWrappers::MPI::BlockVector >::reinit_domain_vector
static void reinit_domain_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
Definition: petsc_block_vector.h:537
internal::LinearOperatorImplementation::ReinitHelper< PETScWrappers::MPI::BlockVector >::reinit_range_vector
static void reinit_range_vector(const Matrix &matrix, PETScWrappers::MPI::BlockVector &v, bool)
Definition: petsc_block_vector.h:527
PETScWrappers::MPI::BlockVector::ExcIteratorRangeDoesNotMatchVectorSize
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
BlockVectorBase< Vector >::reference
typename BlockType::reference reference
Definition: block_vector_base.h:475
BlockVectorBase< Vector >::const_reference
typename BlockType::const_reference const_reference
Definition: block_vector_base.h:476
PETScWrappers::MPI::BlockVector::size_type
BaseClass::size_type size_type
Definition: petsc_block_vector.h:82
PETScWrappers::MPI::BlockVector::BlockVector
BlockVector()=default
MatrixFreeOperators::BlockHelper::n_blocks
std::enable_if< IsBlockVector< VectorType >::value, unsigned int >::type n_blocks(const VectorType &vector)
Definition: operators.h:49
PETScWrappers::MPI::BlockVector::get_mpi_communicator
const MPI_Comm & get_mpi_communicator() const
Definition: petsc_block_vector.h:447
BlockVectorBase< Vector >::collect_sizes
void collect_sizes()
PETScWrappers::MPI::BlockVector::ExcNonMatchingBlockVectors
static ::ExceptionBase & ExcNonMatchingBlockVectors()
PETScWrappers::MPI::BlockVector::swap
void swap(BlockVector &v)
Definition: petsc_block_vector.h:465
PETScWrappers::MPI::BlockVector::const_reference
BaseClass::const_reference const_reference
Definition: petsc_block_vector.h:81
LinearAlgebraDealII::BlockVector
BlockVector< double > BlockVector
Definition: generic_linear_algebra.h:48
PETScWrappers::MPI::BlockVector::has_ghost_elements
bool has_ghost_elements() const
Definition: petsc_block_vector.h:453
is_serial_vector
Definition: vector_type_traits.h:49
BlockVectorBase< Vector >::components
std::vector< Vector > components
Definition: block_vector_base.h:956
vector_type_traits.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
BlockVectorBase< Vector >::n_blocks
unsigned int n_blocks() const
BlockVectorBase
Definition: block_vector_base.h:450
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
PETScWrappers::MPI::BlockVector::pointer
BaseClass::pointer pointer
Definition: petsc_block_vector.h:78
BlockVectorBase< Vector >::BlockType
Vector BlockType
Definition: block_vector_base.h:456
internal::BlockVectorIterators::Iterator
Definition: block_vector_base.h:139
PETScWrappers::MPI::BlockVector::value_type
BaseClass::value_type value_type
Definition: petsc_block_vector.h:77
PETScWrappers::MPI::BlockVector::reference
BaseClass::reference reference
Definition: petsc_block_vector.h:80
BlockVectorBase< Vector >::const_pointer
const value_type * const_pointer
Definition: block_vector_base.h:470
BlockVectorBase< Vector >::size_type
types::global_dof_index size_type
Definition: block_vector_base.h:477
BlockVectorBase< Vector >::pointer
value_type * pointer
Definition: block_vector_base.h:469
PETScWrappers::MPI::BlockVector
Definition: petsc_block_vector.h:61
BlockVectorBase< Vector >::operator=
BlockVectorBase & operator=(const value_type s)
PETScWrappers::MPI::Vector
Definition: petsc_vector.h:158
BlockVectorBase< Vector >::const_iterator
::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > const_iterator
Definition: block_vector_base.h:474
BlockVectorBase< Vector >::block_indices
BlockIndices block_indices
Definition: block_vector_base.h:962
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
BlockVectorBase::get_block_indices
const BlockIndices & get_block_indices() const
BlockVectorBase< Vector >::iterator
::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > iterator
Definition: block_vector_base.h:472
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
PETScWrappers::MPI::BlockVector::BlockType
BaseClass::BlockType BlockType
Definition: petsc_block_vector.h:72
PETScWrappers::MPI::BlockVector::const_pointer
BaseClass::const_pointer const_pointer
Definition: petsc_block_vector.h:79
BlockIndices::reinit
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
Definition: block_indices.h:260
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
petsc_vector.h
config.h
internal
Definition: aligned_vector.h:369
PETScWrappers::MPI::BlockVector::print
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
Definition: petsc_block_vector.h:475
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
BlockVectorBase::block
BlockType & block(const unsigned int i)
internal::LinearOperatorImplementation::ReinitHelper
Definition: block_vector.h:502