Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
la_parallel_block_vector.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1999 - 2021 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.md at
12// the top level directory of deal.II.
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
23
29
30#include <cstdio>
31#include <vector>
32
34
35
36// Forward declarations
37#ifndef DOXYGEN
38# ifdef DEAL_II_WITH_PETSC
39namespace PETScWrappers
40{
41 namespace MPI
42 {
43 class BlockVector;
44 }
45} // namespace PETScWrappers
46# endif
47
48# ifdef DEAL_II_WITH_TRILINOS
49namespace TrilinosWrappers
50{
51 namespace MPI
52 {
53 class BlockVector;
54 }
55} // namespace TrilinosWrappers
56# endif
57#endif
58
59namespace LinearAlgebra
60{
61 namespace distributed
62 {
83 template <typename Number>
84 class BlockVector : public BlockVectorBase<Vector<Number>>,
85 public VectorSpaceVector<Number>
86 {
87 public:
100 static constexpr unsigned int communication_block_size = 20;
101
106
111
117 using pointer = typename BaseClass::pointer;
124
129
140 explicit BlockVector(const size_type num_blocks = 0,
141 const size_type block_size = 0);
142
148
161 template <typename OtherNumber>
163
168 BlockVector(const std::vector<size_type> &block_sizes);
169
174 BlockVector(const std::vector<IndexSet> &local_ranges,
175 const std::vector<IndexSet> &ghost_indices,
176 const MPI_Comm & communicator);
177
181 BlockVector(const std::vector<IndexSet> &local_ranges,
182 const MPI_Comm & communicator);
183
192 virtual ~BlockVector() override = default;
193
198 virtual BlockVector &
199 operator=(const value_type s) override;
200
207
212 template <class Number2>
215
221
222#ifdef DEAL_II_WITH_PETSC
232#endif
233
234#ifdef DEAL_II_WITH_TRILINOS
245#endif
246
260 void
261 reinit(const size_type num_blocks,
262 const size_type block_size = 0,
263 const bool omit_zeroing_entries = false);
264
284 void
285 reinit(const std::vector<size_type> &N,
286 const bool omit_zeroing_entries = false);
287
302 template <typename Number2>
303 void
305 const bool omit_zeroing_entries = false);
306
330 virtual void
332
341 void
343
356
357
366 void
368
372 bool
374
379 template <typename OtherNumber>
380 void
381 add(const std::vector<size_type> & indices,
382 const ::Vector<OtherNumber> &values);
383
388 void
389 sadd(const Number s, const BlockVector<Number> &V);
390
396 virtual bool
397 all_zero() const override;
398
402 virtual Number
403 mean_value() const override;
404
410 lp_norm(const real_type p) const;
411
429 void
432
437
442 virtual void
444 const bool omit_zeroing_entries = false) override;
445
449 virtual BlockVector<Number> &
450 operator*=(const Number factor) override;
451
455 virtual BlockVector<Number> &
456 operator/=(const Number factor) override;
457
461 virtual BlockVector<Number> &
463
467 virtual BlockVector<Number> &
469
478 virtual void
480 VectorOperation::values operation,
481 std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
482 communication_pattern = {}) override;
483
487 virtual Number
488 operator*(const VectorSpaceVector<Number> &V) const override;
489
506 template <typename FullMatrixType>
507 void
509 const BlockVector<Number> &V,
510 const bool symmetric = false) const;
511
528 template <typename FullMatrixType>
529 Number
531 const BlockVector<Number> &V,
532 const bool symmetric = false) const;
533
543 template <typename FullMatrixType>
544 void
546 const FullMatrixType &matrix,
547 const Number s = Number(0.),
548 const Number b = Number(1.)) const;
549
553 virtual void
554 add(const Number a) override;
555
559 virtual void
560 add(const Number a, const VectorSpaceVector<Number> &V) override;
561
565 virtual void
566 add(const Number a,
568 const Number b,
569 const VectorSpaceVector<Number> &W) override;
570
575 virtual void
576 add(const std::vector<size_type> &indices,
577 const std::vector<Number> & values);
578
583 virtual void
584 sadd(const Number s,
585 const Number a,
586 const VectorSpaceVector<Number> &V) override;
587
593 virtual void
594 scale(const VectorSpaceVector<Number> &scaling_factors) override;
595
599 virtual void
600 equ(const Number a, const VectorSpaceVector<Number> &V) override;
601
606 virtual real_type
607 l1_norm() const override;
608
613 virtual real_type
614 l2_norm() const override;
615
620 norm_sqr() const;
621
626 virtual real_type
627 linfty_norm() const override;
628
649 virtual Number
650 add_and_dot(const Number a,
652 const VectorSpaceVector<Number> &W) override;
653
658 virtual size_type
659 size() const override;
660
672 virtual ::IndexSet
673 locally_owned_elements() const override;
674
678 virtual void
679 print(std::ostream & out,
680 const unsigned int precision = 3,
681 const bool scientific = true,
682 const bool across = true) const override;
683
687 virtual std::size_t
688 memory_consumption() const override;
690
702
708 };
709
712 } // end of namespace distributed
713
714} // end of namespace LinearAlgebra
715
716
724template <typename Number>
725inline void
728{
729 u.swap(v);
730}
731
732
737template <typename Number>
739 : std::false_type
740{};
741
742
744
745#ifdef DEAL_II_MSVC
746# include <deal.II/lac/la_parallel_block_vector.templates.h>
747#endif
748
749#endif
void swap(LinearAlgebra::distributed::BlockVector< Number > &u, LinearAlgebra::distributed::BlockVector< Number > &v)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_DEPRECATED_EARLY
Definition: config.h:170
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define DeclException0(Exception0)
Definition: exceptions.h:470
::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > iterator
typename BlockType::const_reference const_reference
typename BlockType::value_type value_type
::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > const_iterator
typename BlockType::reference reference
typename BlockType::real_type real_type
virtual real_type linfty_norm() const override
void swap(BlockVector< Number > &v)
BlockVector(const size_type num_blocks=0, const size_type block_size=0)
virtual real_type l1_norm() const override
static ::ExceptionBase & ExcVectorTypeNotCompatible()
virtual void add(const Number a, const VectorSpaceVector< Number > &V, const Number b, const VectorSpaceVector< Number > &W) override
virtual void add(const std::vector< size_type > &indices, const std::vector< Number > &values)
virtual void reinit(const VectorSpaceVector< Number > &V, const bool omit_zeroing_entries=false) override
virtual ::IndexSet locally_owned_elements() const override
static constexpr unsigned int communication_block_size
virtual ~BlockVector() override=default
virtual void sadd(const Number s, const Number a, const VectorSpaceVector< Number > &V) override
virtual BlockVector< Number > & operator*=(const Number factor) override
virtual BlockVector< Number > & operator/=(const Number factor) override
BlockVector(const std::vector< size_type > &block_sizes)
virtual BlockVector< Number > & operator-=(const VectorSpaceVector< Number > &V) override
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()
virtual void scale(const VectorSpaceVector< Number > &scaling_factors) override
BlockVector(const BlockVector< OtherNumber > &v)
BlockVector(const BlockVector< Number > &V)
void reinit(const size_type num_blocks, const size_type block_size=0, const bool omit_zeroing_entries=false)
virtual std::size_t memory_consumption() const override
virtual Number mean_value() const override
BlockVector< Number > & operator=(const TrilinosWrappers::MPI::BlockVector &trilinos_vec)
BlockVector< Number > & operator=(const PETScWrappers::MPI::BlockVector &petsc_vec)
virtual bool all_zero() const override
virtual void compress(::VectorOperation::values operation) override
void add(const std::vector< size_type > &indices, const ::Vector< OtherNumber > &values)
BlockVector & operator=(const Vector< Number > &V)
virtual real_type l2_norm() const override
virtual Number operator*(const VectorSpaceVector< Number > &V) const override
Number multivector_inner_product_with_metric(const FullMatrixType &matrix, const BlockVector< Number > &V, const bool symmetric=false) const
virtual BlockVector< Number > & operator+=(const VectorSpaceVector< Number > &V) override
virtual void add(const Number a, const VectorSpaceVector< Number > &V) override
typename BaseClass::const_reference const_reference
real_type lp_norm(const real_type p) const
virtual BlockVector & operator=(const value_type s) override
BlockVector & operator=(const BlockVector< Number2 > &V)
void reinit(const BlockVector< Number2 > &V, const bool omit_zeroing_entries=false)
BlockVector & operator=(const BlockVector &V)
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const override
void sadd(const Number s, const BlockVector< Number > &V)
virtual void add(const Number a) override
BlockVector(const std::vector< IndexSet > &local_ranges, const std::vector< IndexSet > &ghost_indices, const MPI_Comm &communicator)
void reinit(const std::vector< size_type > &N, const bool omit_zeroing_entries=false)
virtual void equ(const Number a, const VectorSpaceVector< Number > &V) override
typename BaseClass::const_pointer const_pointer
typename BaseClass::const_iterator const_iterator
BlockVector(const std::vector< IndexSet > &local_ranges, const MPI_Comm &communicator)
void multivector_inner_product(FullMatrixType &matrix, const BlockVector< Number > &V, const bool symmetric=false) const
void mmult(BlockVector< Number > &V, const FullMatrixType &matrix, const Number s=Number(0.), const Number b=Number(1.)) const
virtual Number add_and_dot(const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W) override
virtual size_type size() const override
constexpr int block_size
Definition: cuda_size.h:29
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
static const char N
static const char V
BlockVector< double > BlockVector
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)