Reference documentation for deal.II version 8.5.1
petsc_parallel_block_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__petsc_parallel_block_vector_h
17 #define dealii__petsc_parallel_block_vector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_PETSC
23 
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>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
33 namespace PETScWrappers
34 {
35  // forward declaration
36  class BlockVector;
37 
38  namespace MPI
39  {
40 
62  class BlockVector : public BlockVectorBase<Vector>
63  {
64  public:
69 
74 
78  typedef BaseClass::value_type value_type;
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;
86 
90  BlockVector ();
91 
98  explicit BlockVector (const unsigned int n_blocks,
99  const MPI_Comm &communicator,
100  const size_type block_size,
101  const size_type local_size);
102 
107  BlockVector (const BlockVector &V);
108 
116  BlockVector (const std::vector<size_type> &block_sizes,
117  const MPI_Comm &communicator,
118  const std::vector<size_type> &local_elements);
119 
124  explicit BlockVector (const std::vector<IndexSet> &parallel_partitioning,
125  const MPI_Comm &communicator = MPI_COMM_WORLD);
126 
130  BlockVector (const std::vector<IndexSet> &parallel_partitioning,
131  const std::vector<IndexSet> &ghost_indices,
132  const MPI_Comm &communicator);
133 
134 
135 
139  ~BlockVector ();
140 
146 
150  BlockVector &
151  operator= (const BlockVector &V);
152 
169  BlockVector &
171 
181  void reinit (const unsigned int n_blocks,
182  const MPI_Comm &communicator,
183  const size_type block_size,
184  const size_type local_size,
185  const bool omit_zeroing_entries = false);
186 
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);
211 
226  void reinit (const BlockVector &V,
227  const bool omit_zeroing_entries=false);
228 
233  void reinit (const std::vector<IndexSet> &parallel_partitioning,
234  const MPI_Comm &communicator);
235 
239  void reinit (const std::vector<IndexSet> &parallel_partitioning,
240  const std::vector<IndexSet> &ghost_entries,
241  const MPI_Comm &communicator);
242 
249  void reinit (const unsigned int num_blocks);
250 
254  bool has_ghost_elements() const;
255 
260  const MPI_Comm &get_mpi_communicator () const;
261 
279  void swap (BlockVector &v);
280 
284  void print (std::ostream &out,
285  const unsigned int precision = 3,
286  const bool scientific = true,
287  const bool across = true) const;
288 
297  };
298 
301  /*----------------------- Inline functions ----------------------------------*/
302 
303 
304  inline
306  {}
307 
308 
309 
310  inline
311  BlockVector::BlockVector (const unsigned int n_blocks,
312  const MPI_Comm &communicator,
313  const size_type block_size,
314  const size_type local_size)
315  {
316  reinit (n_blocks, communicator, block_size, local_size);
317  }
318 
319 
320 
321  inline
322  BlockVector::BlockVector (const std::vector<size_type> &block_sizes,
323  const MPI_Comm &communicator,
324  const std::vector<size_type> &local_elements)
325  {
326  reinit (block_sizes, communicator, local_elements, false);
327  }
328 
329 
330  inline
332  :
334  {
335  this->components.resize (v.n_blocks());
336  this->block_indices = v.block_indices;
337 
338  for (unsigned int i=0; i<this->n_blocks(); ++i)
339  this->components[i] = v.components[i];
340  }
341 
342  inline
343  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
344  const MPI_Comm &communicator)
345  {
346  reinit(parallel_partitioning, communicator);
347  }
348 
349  inline
350  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
351  const std::vector<IndexSet> &ghost_indices,
352  const MPI_Comm &communicator)
353  {
354  reinit(parallel_partitioning, ghost_indices, communicator);
355  }
356 
357  inline
358  BlockVector &
360  {
362  return *this;
363  }
364 
365  inline
366  BlockVector &
368  {
369  // we only allow assignment to vectors with the same number of blocks
370  // or to an empty BlockVector
371  Assert (n_blocks() == 0 || n_blocks() == v.n_blocks(),
373 
374  if (this->n_blocks() != v.n_blocks())
375  reinit(v.n_blocks());
376 
377  for (size_type i=0; i<this->n_blocks(); ++i)
378  this->components[i] = v.block(i);
379 
380  collect_sizes();
381 
382  return *this;
383  }
384 
385  inline
387  {}
388 
389 
390  inline
391  void
392  BlockVector::reinit (const unsigned int n_blocks,
393  const MPI_Comm &communicator,
394  const size_type block_size,
395  const size_type local_size,
396  const bool omit_zeroing_entries)
397  {
398  reinit(std::vector<size_type>(n_blocks, block_size),
399  communicator,
400  std::vector<size_type>(n_blocks, local_size),
401  omit_zeroing_entries);
402  }
403 
404 
405 
406  inline
407  void
408  BlockVector::reinit (const std::vector<size_type> &block_sizes,
409  const MPI_Comm &communicator,
410  const std::vector<size_type> &local_sizes,
411  const bool omit_zeroing_entries)
412  {
413  this->block_indices.reinit (block_sizes);
414  if (this->components.size() != this->n_blocks())
415  this->components.resize(this->n_blocks());
416 
417  for (unsigned int i=0; i<this->n_blocks(); ++i)
418  this->components[i].reinit(communicator, block_sizes[i],
419  local_sizes[i], omit_zeroing_entries);
420  }
421 
422 
423  inline
424  void
426  const bool omit_zeroing_entries)
427  {
428  this->block_indices = v.get_block_indices();
429  if (this->components.size() != this->n_blocks())
430  this->components.resize(this->n_blocks());
431 
432  for (unsigned int i=0; i<this->n_blocks(); ++i)
433  block(i).reinit(v.block(i), omit_zeroing_entries);
434  }
435 
436  inline
437  void
438  BlockVector::reinit (const std::vector<IndexSet> &parallel_partitioning,
439  const MPI_Comm &communicator)
440  {
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();
444 
445  this->block_indices.reinit(sizes);
446  if (this->components.size() != this->n_blocks())
447  this->components.resize(this->n_blocks());
448 
449  for (unsigned int i=0; i<this->n_blocks(); ++i)
450  block(i).reinit(parallel_partitioning[i], communicator);
451  }
452 
453  inline
454  void
455  BlockVector::reinit (const std::vector<IndexSet> &parallel_partitioning,
456  const std::vector<IndexSet> &ghost_entries,
457  const MPI_Comm &communicator)
458  {
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();
462 
463  this->block_indices.reinit(sizes);
464  if (this->components.size() != this->n_blocks())
465  this->components.resize(this->n_blocks());
466 
467  for (unsigned int i=0; i<this->n_blocks(); ++i)
468  block(i).reinit(parallel_partitioning[i], ghost_entries[i], communicator);
469  }
470 
471 
472 
473  inline
474  const MPI_Comm &
476  {
477  return block(0).get_mpi_communicator();
478  }
479 
480  inline
481  bool
483  {
484  bool ghosted=block(0).has_ghost_elements();
485 #ifdef DEBUG
486  for (unsigned int i=0; i<this->n_blocks(); ++i)
487  Assert(block(i).has_ghost_elements()==ghosted, ExcInternalError());
488 #endif
489  return ghosted;
490  }
491 
492 
493  inline
494  void
496  {
497  std::swap(this->components, v.components);
498 
500  }
501 
502 
503 
504  inline
505  void
506  BlockVector::print (std::ostream &out,
507  const unsigned int precision,
508  const bool scientific,
509  const bool across) const
510  {
511  for (unsigned int i=0; i<this->n_blocks(); ++i)
512  {
513  if (across)
514  out << 'C' << i << ':';
515  else
516  out << "Component " << i << std::endl;
517  this->components[i].print(out, precision, scientific, across);
518  }
519  }
520 
521 
522 
531  inline
532  void swap (BlockVector &u,
533  BlockVector &v)
534  {
535  u.swap (v);
536  }
537 
538  }
539 
540 }
541 
542 namespace internal
543 {
544  namespace LinearOperator
545  {
546  template <typename> class ReinitHelper;
547 
552  template<>
553  class ReinitHelper<PETScWrappers::MPI::BlockVector>
554  {
555  public:
556  template <typename Matrix>
557  static
558  void reinit_range_vector (const Matrix &matrix,
560  bool /*omit_zeroing_entries*/)
561  {
562  v.reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator());
563  }
564 
565  template <typename Matrix>
566  static
567  void reinit_domain_vector(const Matrix &matrix,
569  bool /*omit_zeroing_entries*/)
570  {
571  v.reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator());
572  }
573  };
574 
575  } /* namespace LinearOperator */
576 } /* namespace internal */
577 
578 
584 template <>
585 struct is_serial_vector< PETScWrappers::MPI::BlockVector > : std_cxx11::false_type
586 {
587 };
588 
589 
590 DEAL_II_NAMESPACE_CLOSE
591 
592 #endif // DEAL_II_WITH_PETSC
593 
594 #endif
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 &u, BlockVector &v)
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
static ::ExceptionBase & ExcNonMatchingBlockVectors()
std::size_t size() const
const BlockIndices & get_block_indices() const
#define Assert(cond, exc)
Definition: exceptions.h:313
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)
Definition: exceptions.h:541
std::vector< Vector > components
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
unsigned int n_blocks() const
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
BlockVector & operator=(const value_type s)
BlockType & block(const unsigned int i)
BlockVectorBase & operator=(const value_type s)
static ::ExceptionBase & ExcInternalError()