Reference documentation for deal.II version 8.4.2
petsc_parallel_block_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2016 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 
29 DEAL_II_NAMESPACE_OPEN
30 
31 
32 namespace PETScWrappers
33 {
34  // forward declaration
35  class BlockVector;
36 
37  namespace MPI
38  {
39 
61  class BlockVector : public BlockVectorBase<Vector>
62  {
63  public:
68 
73 
77  typedef BaseClass::value_type value_type;
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;
85 
89  BlockVector ();
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 ();
139 
144  BlockVector &operator= (const value_type s);
145 
149  BlockVector &
150  operator= (const BlockVector &V);
151 
168  BlockVector &
170 
180  void reinit (const unsigned int n_blocks,
181  const MPI_Comm &communicator,
182  const size_type block_size,
183  const size_type local_size,
184  const bool omit_zeroing_entries = false);
185 
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);
210 
225  void reinit (const BlockVector &V,
226  const bool omit_zeroing_entries=false);
227 
232  void reinit (const std::vector<IndexSet> &parallel_partitioning,
233  const MPI_Comm &communicator);
234 
238  void reinit (const std::vector<IndexSet> &parallel_partitioning,
239  const std::vector<IndexSet> &ghost_entries,
240  const MPI_Comm &communicator);
241 
248  void reinit (const unsigned int num_blocks);
249 
253  bool has_ghost_elements() const;
254 
259  const MPI_Comm &get_mpi_communicator () const;
260 
278  void swap (BlockVector &v);
279 
283  void print (std::ostream &out,
284  const unsigned int precision = 3,
285  const bool scientific = true,
286  const bool across = true) const;
287 
291  DeclException0 (ExcIteratorRangeDoesNotMatchVectorSize);
295  DeclException0 (ExcNonMatchingBlockVectors);
296  };
297 
300  /*----------------------- Inline functions ----------------------------------*/
301 
302 
303  inline
305  {}
306 
307 
308 
309  inline
310  BlockVector::BlockVector (const unsigned int n_blocks,
311  const MPI_Comm &communicator,
312  const size_type block_size,
313  const size_type local_size)
314  {
315  reinit (n_blocks, communicator, block_size, local_size);
316  }
317 
318 
319 
320  inline
321  BlockVector::BlockVector (const std::vector<size_type> &block_sizes,
322  const MPI_Comm &communicator,
323  const std::vector<size_type> &local_elements)
324  {
325  reinit (block_sizes, communicator, local_elements, false);
326  }
327 
328 
329  inline
331  :
333  {
334  this->components.resize (v.n_blocks());
335  this->block_indices = v.block_indices;
336 
337  for (unsigned int i=0; i<this->n_blocks(); ++i)
338  this->components[i] = v.components[i];
339  }
340 
341  inline
342  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
343  const MPI_Comm &communicator)
344  {
345  reinit(parallel_partitioning, communicator);
346  }
347 
348  inline
349  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
350  const std::vector<IndexSet> &ghost_indices,
351  const MPI_Comm &communicator)
352  {
353  reinit(parallel_partitioning, ghost_indices, communicator);
354  }
355 
356  inline
357  BlockVector &
359  {
361  return *this;
362  }
363 
364  inline
365  BlockVector &
367  {
368  // we only allow assignment to vectors with the same number of blocks
369  // or to an empty BlockVector
370  Assert (n_blocks() == 0 || n_blocks() == v.n_blocks(),
371  ExcDimensionMismatch(n_blocks(), v.n_blocks()));
372 
373  if (this->n_blocks() != v.n_blocks())
374  reinit(v.n_blocks());
375 
376  for (size_type i=0; i<this->n_blocks(); ++i)
377  this->components[i] = v.block(i);
378 
379  collect_sizes();
380 
381  return *this;
382  }
383 
384  inline
386  {}
387 
388 
389  inline
390  void
391  BlockVector::reinit (const unsigned int n_blocks,
392  const MPI_Comm &communicator,
393  const size_type block_size,
394  const size_type local_size,
395  const bool omit_zeroing_entries)
396  {
397  reinit(std::vector<size_type>(n_blocks, block_size),
398  communicator,
399  std::vector<size_type>(n_blocks, local_size),
400  omit_zeroing_entries);
401  }
402 
403 
404 
405  inline
406  void
407  BlockVector::reinit (const std::vector<size_type> &block_sizes,
408  const MPI_Comm &communicator,
409  const std::vector<size_type> &local_sizes,
410  const bool omit_zeroing_entries)
411  {
412  this->block_indices.reinit (block_sizes);
413  if (this->components.size() != this->n_blocks())
414  this->components.resize(this->n_blocks());
415 
416  for (unsigned int i=0; i<this->n_blocks(); ++i)
417  this->components[i].reinit(communicator, block_sizes[i],
418  local_sizes[i], omit_zeroing_entries);
419  }
420 
421 
422  inline
423  void
425  const bool omit_zeroing_entries)
426  {
427  this->block_indices = v.get_block_indices();
428  if (this->components.size() != this->n_blocks())
429  this->components.resize(this->n_blocks());
430 
431  for (unsigned int i=0; i<this->n_blocks(); ++i)
432  block(i).reinit(v.block(i), omit_zeroing_entries);
433  }
434 
435  inline
436  void
437  BlockVector::reinit (const std::vector<IndexSet> &parallel_partitioning,
438  const MPI_Comm &communicator)
439  {
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();
443 
444  this->block_indices.reinit(sizes);
445  if (this->components.size() != this->n_blocks())
446  this->components.resize(this->n_blocks());
447 
448  for (unsigned int i=0; i<this->n_blocks(); ++i)
449  block(i).reinit(parallel_partitioning[i], communicator);
450  }
451 
452  inline
453  void
454  BlockVector::reinit (const std::vector<IndexSet> &parallel_partitioning,
455  const std::vector<IndexSet> &ghost_entries,
456  const MPI_Comm &communicator)
457  {
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();
461 
462  this->block_indices.reinit(sizes);
463  if (this->components.size() != this->n_blocks())
464  this->components.resize(this->n_blocks());
465 
466  for (unsigned int i=0; i<this->n_blocks(); ++i)
467  block(i).reinit(parallel_partitioning[i], ghost_entries[i], communicator);
468  }
469 
470 
471 
472  inline
473  const MPI_Comm &
475  {
476  return block(0).get_mpi_communicator();
477  }
478 
479  inline
480  bool
482  {
483  bool ghosted=block(0).has_ghost_elements();
484 #ifdef DEBUG
485  for (unsigned int i=0; i<this->n_blocks(); ++i)
486  Assert(block(i).has_ghost_elements()==ghosted, ExcInternalError());
487 #endif
488  return ghosted;
489  }
490 
491 
492  inline
493  void
495  {
496  std::swap(this->components, v.components);
497 
499  }
500 
501 
502 
503  inline
504  void
505  BlockVector::print (std::ostream &out,
506  const unsigned int precision,
507  const bool scientific,
508  const bool across) const
509  {
510  for (unsigned int i=0; i<this->n_blocks(); ++i)
511  {
512  if (across)
513  out << 'C' << i << ':';
514  else
515  out << "Component " << i << std::endl;
516  this->components[i].print(out, precision, scientific, across);
517  }
518  }
519 
520 
521 
530  inline
531  void swap (BlockVector &u,
532  BlockVector &v)
533  {
534  u.swap (v);
535  }
536 
537  }
538 
539 }
540 
541 DEAL_II_NAMESPACE_CLOSE
542 
543 #endif // DEAL_II_WITH_PETSC
544 
545 #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)
DeclException0(ExcIteratorRangeDoesNotMatchVectorSize)
std::size_t size() const
const BlockIndices & get_block_indices() const
#define Assert(cond, exc)
Definition: exceptions.h:294
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
std::vector< Vector > components
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)