Reference documentation for deal.II version 8.5.1
la_parallel_block_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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__la_parallel_block_vector_h
17 #define dealii__la_parallel_block_vector_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/lac/block_indices.h>
23 #include <deal.II/lac/block_vector_base.h>
24 #include <deal.II/lac/la_parallel_vector.h>
25 #include <deal.II/lac/vector_type_traits.h>
26 
27 
28 #include <cstdio>
29 #include <vector>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
34 #ifdef DEAL_II_WITH_PETSC
35 namespace PETScWrappers
36 {
37  namespace MPI
38  {
39  class BlockVector;
40  }
41 }
42 #endif
43 
44 #ifdef DEAL_II_WITH_TRILINOS
45 namespace TrilinosWrappers
46 {
47  namespace MPI
48  {
49  class BlockVector;
50  }
51 }
52 #endif
53 
54 namespace LinearAlgebra
55 {
56  namespace distributed
57  {
58 
80  template <typename Number>
81  class BlockVector : public BlockVectorBase<Vector<Number> >,
82  public VectorSpaceVector<Number>
83  {
84  public:
89 
93  typedef typename BaseClass::BlockType BlockType;
94 
98  typedef typename BaseClass::value_type value_type;
99  typedef typename BaseClass::real_type real_type;
100  typedef typename BaseClass::pointer pointer;
101  typedef typename BaseClass::const_pointer const_pointer;
102  typedef typename BaseClass::reference reference;
103  typedef typename BaseClass::const_reference const_reference;
104  typedef typename BaseClass::size_type size_type;
105  typedef typename BaseClass::iterator iterator;
106  typedef typename BaseClass::const_iterator const_iterator;
107 
112 
123  explicit BlockVector (const size_type num_blocks = 0,
124  const size_type block_size = 0);
125 
130  BlockVector (const BlockVector<Number> &V);
131 
132 
133 #ifndef DEAL_II_EXPLICIT_CONSTRUCTOR_BUG
134 
146  template <typename OtherNumber>
147  explicit
149 #endif
150 
155  BlockVector (const std::vector<size_type> &block_sizes);
156 
161  BlockVector (const std::vector<IndexSet> &local_ranges,
162  const std::vector<IndexSet> &ghost_indices,
163  const MPI_Comm communicator);
164 
168  BlockVector (const std::vector<IndexSet> &local_ranges,
169  const MPI_Comm communicator);
170 
176 
181  BlockVector &
182  operator= (const BlockVector &V);
183 
188  template <class Number2>
189  BlockVector &
191 
195  BlockVector &
196  operator= (const Vector<Number> &V);
197 
198 #ifdef DEAL_II_WITH_PETSC
199 
208 #endif
209 
210 #ifdef DEAL_II_WITH_TRILINOS
211 
220  operator = (const TrilinosWrappers::MPI::BlockVector &trilinos_vec);
221 #endif
222 
236  void reinit (const size_type num_blocks,
237  const size_type block_size = 0,
238  const bool omit_zeroing_entries = false);
239 
259  void reinit (const std::vector<size_type> &N,
260  const bool omit_zeroing_entries=false);
261 
276  template <typename Number2>
277  void reinit (const BlockVector<Number2> &V,
278  const bool omit_zeroing_entries=false);
279 
304  void compress (::VectorOperation::values operation);
305 
314  void update_ghost_values () const;
315 
324  void zero_out_ghosts ();
325 
329  bool has_ghost_elements() const;
330 
335  template <typename OtherNumber>
336  void add (const std::vector<size_type> &indices,
337  const ::Vector<OtherNumber> &values);
338 
343  void sadd (const Number s,
344  const BlockVector<Number> &V);
345 
351  void equ (const Number a, const BlockVector<Number> &u,
352  const Number b, const BlockVector<Number> &v) DEAL_II_DEPRECATED;
353 
359  void sadd (const Number s,
360  const Number a,
361  const BlockVector<Number> &V,
362  const Number b,
363  const BlockVector<Number> &W) DEAL_II_DEPRECATED;
364 
370  bool all_zero () const;
371 
375  Number mean_value () const;
376 
381  real_type lp_norm (const real_type p) const;
382 
400  void swap (BlockVector<Number> &v);
402 
407 
411  virtual BlockVector<Number> &operator*= (const Number factor);
412 
416  virtual BlockVector<Number> &operator/= (const Number factor);
417 
422 
427 
436  virtual void import(const LinearAlgebra::ReadWriteVector<Number> &V,
437  VectorOperation::values operation,
438  std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
439  std_cxx11::shared_ptr<const CommunicationPatternBase> ());
440 
444  virtual Number operator* (const VectorSpaceVector<Number> &V) const;
445 
449  virtual void add(const Number a);
450 
454  virtual void add(const Number a, const VectorSpaceVector<Number> &V);
455 
459  virtual void add(const Number a, const VectorSpaceVector<Number> &V,
460  const Number b, const VectorSpaceVector<Number> &W);
461 
466  virtual void add (const std::vector<size_type> &indices,
467  const std::vector<Number> &values);
468 
473  virtual void sadd(const Number s, const Number a,
474  const VectorSpaceVector<Number> &V);
475 
481  virtual void scale(const VectorSpaceVector<Number> &scaling_factors);
482 
486  virtual void equ(const Number a, const VectorSpaceVector<Number> &V);
487 
492  virtual real_type l1_norm() const;
493 
498  virtual real_type l2_norm() const;
499 
504  virtual real_type linfty_norm() const;
505 
522  virtual Number add_and_dot(const Number a,
523  const VectorSpaceVector<Number> &V,
524  const VectorSpaceVector<Number> &W);
525 
530  virtual size_type size() const;
531 
543  virtual ::IndexSet locally_owned_elements() const;
544 
548  virtual void print(std::ostream &out,
549  const unsigned int precision=3,
550  const bool scientific=true,
551  const bool across=true) const;
552 
556  virtual std::size_t memory_consumption() const;
558 
570 
576  };
577 
580  } // end of namespace distributed
581 
582 } // end of namespace LinearAlgebra
583 
584 
593 template <typename Number>
594 inline
597 {
598  u.swap (v);
599 }
600 
601 
607 template <typename Number>
608 struct is_serial_vector< LinearAlgebra::distributed::BlockVector< Number > > : std_cxx11::false_type
609 {
610 };
611 
612 
613 DEAL_II_NAMESPACE_CLOSE
614 
615 #ifdef DEAL_II_MSVC
616 #include <deal.II/lac/la_parallel_block_vector.templates.h>
617 #endif
618 
619 #endif
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
real_type lp_norm(const real_type p) const
virtual ::IndexSet locally_owned_elements() const
virtual BlockVector< Number > & operator-=(const VectorSpaceVector< Number > &V)
virtual size_type size() const
virtual std::size_t memory_consumption() const
virtual BlockVector< Number > & operator*=(const Number factor)
virtual real_type l2_norm() const
void reinit(const size_type num_blocks, const size_type block_size=0, const bool omit_zeroing_entries=false)
BlockVectorBase< Vector< Number > > BaseClass
static ::ExceptionBase & ExcVectorTypeNotCompatible()
virtual void scale(const VectorSpaceVector< Number > &scaling_factors)
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
virtual Number add_and_dot(const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W)
virtual BlockVector< Number > & operator/=(const Number factor)
virtual Number operator*(const VectorSpaceVector< Number > &V) const
virtual real_type l1_norm() const
void swap(LinearAlgebra::distributed::BlockVector< Number > &u, LinearAlgebra::distributed::BlockVector< Number > &v)
#define DeclException0(Exception0)
Definition: exceptions.h:541
void swap(BlockVector< Number > &v)
virtual real_type linfty_norm() const
BlockVector(const size_type num_blocks=0, const size_type block_size=0)
void equ(const Number a, const BlockVector< Number > &u, const Number b, const BlockVector< Number > &v) 1
void add(const std::vector< size_type > &indices, const ::Vector< OtherNumber > &values)
BlockVector & operator=(const value_type s)
void sadd(const Number s, const BlockVector< Number > &V)
void compress(::VectorOperation::values operation)
virtual BlockVector< Number > & operator+=(const VectorSpaceVector< Number > &V)