Reference documentation for deal.II version 9.0.0
trilinos_parallel_block_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_trilinos_parallel_block_vector_h
17 #define dealii_trilinos_parallel_block_vector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
24 # include <deal.II/lac/trilinos_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 # include <functional>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 // forward declaration
34 template <typename Number> class BlockVectorBase;
35 
40 namespace TrilinosWrappers
41 {
42  // forward declaration
43  namespace MPI
44  {
45  class BlockVector;
46  }
47  class BlockSparseMatrix;
48 
49 
50  namespace MPI
51  {
70  class BlockVector : public ::BlockVectorBase<MPI::Vector>
71  {
72  public:
76  typedef ::BlockVectorBase<MPI::Vector> BaseClass;
77 
82 
86  typedef BaseClass::value_type value_type;
87  typedef BaseClass::pointer pointer;
88  typedef BaseClass::const_pointer const_pointer;
89  typedef BaseClass::reference reference;
90  typedef BaseClass::const_reference const_reference;
91  typedef BaseClass::size_type size_type;
94 
98  BlockVector () = default;
99 
106  explicit BlockVector (const std::vector<IndexSet> &parallel_partitioning,
107  const MPI_Comm &communicator = MPI_COMM_WORLD);
108 
114  BlockVector (const std::vector<IndexSet> &parallel_partitioning,
115  const std::vector<IndexSet> &ghost_values,
116  const MPI_Comm &communicator,
117  const bool vector_writable = false);
118 
123  BlockVector (const BlockVector &v);
124 
129  BlockVector (BlockVector &&v) noexcept;
130 
136  explicit BlockVector (const size_type num_blocks);
137 
141  ~BlockVector () = default;
142 
148 
153 
158  BlockVector &operator= (BlockVector &&v) noexcept;
159 
170  template <typename Number>
171  BlockVector &operator= (const ::BlockVector<Number> &v);
172 
181  void reinit (const std::vector<IndexSet> &parallel_partitioning,
182  const MPI_Comm &communicator = MPI_COMM_WORLD,
183  const bool omit_zeroing_entries = false);
184 
202  void reinit (const std::vector<IndexSet> &partitioning,
203  const std::vector<IndexSet> &ghost_values,
204  const MPI_Comm &communicator = MPI_COMM_WORLD,
205  const bool vector_writable = false);
206 
207 
222  void reinit (const BlockVector &V,
223  const bool omit_zeroing_entries = false);
224 
231  void reinit (const size_type num_blocks);
232 
251  const BlockVector &v);
252 
259  bool has_ghost_elements() 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 
292 
297  };
298 
299 
300 
301  /*----------------------- Inline functions ----------------------------------*/
302  inline
303  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
304  const MPI_Comm &communicator)
305  {
306  reinit (parallel_partitioning, communicator, false);
307  }
308 
309 
310 
311  inline
312  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
313  const std::vector<IndexSet> &ghost_values,
314  const MPI_Comm &communicator,
315  const bool vector_writable)
316  {
317  reinit(parallel_partitioning, ghost_values, communicator,
318  vector_writable);
319  }
320 
321 
322 
323  inline
324  BlockVector::BlockVector (const size_type num_blocks)
325  {
326  reinit (num_blocks);
327  }
328 
329 
330 
331  inline
333  :
334  ::BlockVectorBase<MPI::Vector> ()
335  {
336  this->components.resize (v.n_blocks());
337  this->block_indices = v.block_indices;
338 
339  for (size_type i=0; i<this->n_blocks(); ++i)
340  this->components[i] = v.components[i];
341  }
342 
343 
344 
345  inline
347  {
348  // initialize a minimal, valid object and swap
349  reinit (0);
350  swap(v);
351  }
352 
353 
354 
355  template <typename Number>
356  BlockVector &
357  BlockVector::operator= (const ::BlockVector<Number> &v)
358  {
359  if (n_blocks() != v.n_blocks())
360  {
361  std::vector<size_type> block_sizes (v.n_blocks(), 0);
362  block_indices.reinit (block_sizes);
363  if (components.size() != n_blocks())
364  components.resize(n_blocks());
365  }
366 
367  for (size_type i=0; i<this->n_blocks(); ++i)
368  this->components[i] = v.block(i);
369 
370  collect_sizes();
371 
372  return *this;
373  }
374 
375 
376 
377  inline
378  bool
380  {
381  bool ghosted=block(0).has_ghost_elements();
382 #ifdef DEBUG
383  for (unsigned int i=0; i<this->n_blocks(); ++i)
384  Assert(block(i).has_ghost_elements()==ghosted, ExcInternalError());
385 #endif
386  return ghosted;
387  }
388 
389 
390 
391  inline
392  void
394  {
395  std::swap(this->components, v.components);
396 
398  }
399 
400 
401 
410  inline
411  void swap (BlockVector &u,
412  BlockVector &v)
413  {
414  u.swap (v);
415  }
416 
417  } /* namespace MPI */
418 
419 } /* namespace TrilinosWrappers */
420 
424 namespace internal
425 {
426  namespace LinearOperatorImplementation
427  {
428  template <typename> class ReinitHelper;
429 
434  template <>
435  class ReinitHelper<TrilinosWrappers::MPI::BlockVector>
436  {
437  public:
438  template <typename Matrix>
439  static
440  void reinit_range_vector (const Matrix &matrix,
442  bool omit_zeroing_entries)
443  {
444  v.reinit(matrix.locally_owned_range_indices(),
445  matrix.get_mpi_communicator(), omit_zeroing_entries);
446  }
447 
448  template <typename Matrix>
449  static
450  void reinit_domain_vector(const Matrix &matrix,
452  bool omit_zeroing_entries)
453  {
454  v.reinit(matrix.locally_owned_domain_indices(),
455  matrix.get_mpi_communicator(), omit_zeroing_entries);
456  }
457  };
458 
459  } /* namespace LinearOperator */
460 } /* namespace internal */
461 
462 
468 template <>
469 struct is_serial_vector< TrilinosWrappers::MPI::BlockVector > : std::false_type
470 {
471 };
472 
473 DEAL_II_NAMESPACE_CLOSE
474 
475 #endif // DEAL_II_WITH_TRILINOS
476 
477 #endif
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
static ::ExceptionBase & ExcNonMatchingBlockVectors()
void reinit(const std::vector< IndexSet > &parallel_partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
void swap(BlockVector &u, BlockVector &v)
#define Assert(cond, exc)
Definition: exceptions.h:1142
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
#define DeclException0(Exception0)
Definition: exceptions.h:323
std::vector< MPI::Vector > components
void swap(Vector< Number > &u, Vector< Number > &v)
Definition: vector.h:1312
BlockVector & operator=(const value_type s)
unsigned int n_blocks() const
void import_nonlocal_data_for_fe(const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v)
BlockType & block(const unsigned int i)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
::BlockVectorBase< MPI::Vector > BaseClass
static ::ExceptionBase & ExcInternalError()