Reference documentation for deal.II version 8.5.1
trilinos_parallel_block_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 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__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 BlockVector;
35 
40 namespace TrilinosWrappers
41 {
42  // forward declaration
43  namespace MPI
44  {
45  class BlockVector;
46  }
47  class BlockVector;
48  class BlockSparseMatrix;
49 
50 
51  namespace MPI
52  {
71  class BlockVector : public BlockVectorBase<Vector>
72  {
73  public:
78 
83 
87  typedef BaseClass::value_type value_type;
88  typedef BaseClass::pointer pointer;
89  typedef BaseClass::const_pointer const_pointer;
90  typedef BaseClass::reference reference;
91  typedef BaseClass::const_reference const_reference;
92  typedef BaseClass::size_type size_type;
95 
99  BlockVector ();
100 
108  explicit BlockVector (const std::vector<Epetra_Map> &parallel_partitioning) DEAL_II_DEPRECATED;
109 
116  explicit BlockVector (const std::vector<IndexSet> &parallel_partitioning,
117  const MPI_Comm &communicator = MPI_COMM_WORLD);
118 
124  BlockVector (const std::vector<IndexSet> &parallel_partitioning,
125  const std::vector<IndexSet> &ghost_values,
126  const MPI_Comm &communicator,
127  const bool vector_writable = false);
128 
133  BlockVector (const BlockVector &v);
134 
135 #ifdef DEAL_II_WITH_CXX11
136 
143  BlockVector (BlockVector &&v);
144 #endif
145 
151  explicit BlockVector (const size_type num_blocks);
152 
156  ~BlockVector ();
157 
163 
168 
169 #ifdef DEAL_II_WITH_CXX11
170 
178 #endif
179 
183  BlockVector &
184  operator= (const ::TrilinosWrappers::BlockVector &v);
185 
196  template <typename Number>
197  BlockVector &operator= (const ::BlockVector<Number> &v);
198 
209  void reinit (const std::vector<Epetra_Map> &parallel_partitioning,
210  const bool omit_zeroing_entries = false) DEAL_II_DEPRECATED;
211 
220  void reinit (const std::vector<IndexSet> &parallel_partitioning,
221  const MPI_Comm &communicator = MPI_COMM_WORLD,
222  const bool omit_zeroing_entries = false);
223 
241  void reinit (const std::vector<IndexSet> &partitioning,
242  const std::vector<IndexSet> &ghost_values,
243  const MPI_Comm &communicator = MPI_COMM_WORLD,
244  const bool vector_writable = false);
245 
246 
261  void reinit (const BlockVector &V,
262  const bool omit_zeroing_entries = false);
263 
270  void reinit (const size_type num_blocks);
271 
290  const BlockVector &v);
291 
301  bool is_compressed () const DEAL_II_DEPRECATED;
302 
309  bool has_ghost_elements() const;
310 
328  void swap (BlockVector &v);
329 
333  void print (std::ostream &out,
334  const unsigned int precision = 3,
335  const bool scientific = true,
336  const bool across = true) const;
337 
342 
347  };
348 
349 
350 
351  /*----------------------- Inline functions ----------------------------------*/
352 
353 
354  inline
356  {}
357 
358 
359 
360  inline
361  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
362  const MPI_Comm &communicator)
363  {
364  reinit (parallel_partitioning, communicator, false);
365  }
366 
367 
368 
369  inline
370  BlockVector::BlockVector (const std::vector<IndexSet> &parallel_partitioning,
371  const std::vector<IndexSet> &ghost_values,
372  const MPI_Comm &communicator,
373  const bool vector_writable)
374  {
375  reinit(parallel_partitioning, ghost_values, communicator,
376  vector_writable);
377  }
378 
379 
380 
381  inline
382  BlockVector::BlockVector (const size_type num_blocks)
383  {
384  reinit (num_blocks);
385  }
386 
387 
388 
389  inline
391  :
393  {
394  this->components.resize (v.n_blocks());
395  this->block_indices = v.block_indices;
396 
397  for (size_type i=0; i<this->n_blocks(); ++i)
398  this->components[i] = v.components[i];
399  }
400 
401 
402 
403 #ifdef DEAL_II_WITH_CXX11
404  inline
406  {
407  // initialize a minimal, valid object and swap
408  reinit (0);
409  swap(v);
410  }
411 #endif
412 
413 
414 
415  template <typename Number>
416  BlockVector &
417  BlockVector::operator= (const ::BlockVector<Number> &v)
418  {
419  if (n_blocks() != v.n_blocks())
420  {
421  std::vector<size_type> block_sizes (v.n_blocks(), 0);
422  block_indices.reinit (block_sizes);
423  if (components.size() != n_blocks())
424  components.resize(n_blocks());
425  }
426 
427  for (size_type i=0; i<this->n_blocks(); ++i)
428  this->components[i] = v.block(i);
429 
430  collect_sizes();
431 
432  return *this;
433  }
434 
435 
436 
437  inline
438  bool
440  {
441  bool ghosted=block(0).has_ghost_elements();
442 #ifdef DEBUG
443  for (unsigned int i=0; i<this->n_blocks(); ++i)
444  Assert(block(i).has_ghost_elements()==ghosted, ExcInternalError());
445 #endif
446  return ghosted;
447  }
448 
449 
450 
451  inline
452  void
454  {
455  std::swap(this->components, v.components);
456 
458  }
459 
460 
461 
470  inline
471  void swap (BlockVector &u,
472  BlockVector &v)
473  {
474  u.swap (v);
475  }
476 
477  } /* namespace MPI */
478 
479 } /* namespace TrilinosWrappers */
480 
484 namespace internal
485 {
486  namespace LinearOperator
487  {
488  template <typename> class ReinitHelper;
489 
494  template<>
495  class ReinitHelper<TrilinosWrappers::MPI::BlockVector>
496  {
497  public:
498  template <typename Matrix>
499  static
500  void reinit_range_vector (const Matrix &matrix,
502  bool omit_zeroing_entries)
503  {
504  v.reinit(matrix.range_partitioner(), omit_zeroing_entries);
505  }
506 
507  template <typename Matrix>
508  static
509  void reinit_domain_vector(const Matrix &matrix,
511  bool omit_zeroing_entries)
512  {
513  v.reinit(matrix.domain_partitioner(), omit_zeroing_entries);
514  }
515  };
516 
517  } /* namespace LinearOperator */
518 } /* namespace internal */
519 
520 
521 DEAL_II_NAMESPACE_CLOSE
522 
523 #endif // DEAL_II_WITH_TRILINOS
524 
525 #endif
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
STL namespace.
static ::ExceptionBase & ExcNonMatchingBlockVectors()
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
#define Assert(cond, exc)
Definition: exceptions.h:313
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
#define DeclException0(Exception0)
Definition: exceptions.h:541
void reinit(const std::vector< Epetra_Map > &parallel_partitioning, const bool omit_zeroing_entries=false) 1
std::vector< Vector > components
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
BlockVector & operator=(const value_type s)
unsigned int n_blocks() const
void import_nonlocal_data_for_fe(const TrilinosWrappers::BlockSparseMatrix &m, const BlockVector &v)
void swap(BlockVector &u, 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
static ::ExceptionBase & ExcInternalError()