Reference documentation for deal.II version GIT relicensing-422-gb369f187d8 2024-04-17 18:10:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
la_parallel_block_vector.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_la_parallel_block_vector_h
16#define dealii_la_parallel_block_vector_h
17
18
19#include <deal.II/base/config.h>
20
22
28
29#include <cstdio>
30#include <vector>
31
33
34
35// Forward declarations
36#ifndef DOXYGEN
37# ifdef DEAL_II_WITH_PETSC
38namespace PETScWrappers
39{
40 namespace MPI
41 {
42 class BlockVector;
43 }
44} // namespace PETScWrappers
45# endif
46
47# ifdef DEAL_II_WITH_TRILINOS
48namespace TrilinosWrappers
49{
50 namespace MPI
51 {
52 class BlockVector;
53 }
54} // namespace TrilinosWrappers
55# endif
56#endif
57
58namespace LinearAlgebra
59{
60 namespace distributed
61 {
83 template <typename Number>
84 class BlockVector : public BlockVectorBase<Vector<Number>>
85 {
86 public:
99 static constexpr unsigned int communication_block_size = 20;
100
105
110
116 using pointer = typename BaseClass::pointer;
123
139 explicit BlockVector(const size_type num_blocks = 0,
140 const size_type block_size = 0);
141
147
154 template <typename OtherNumber>
156
161 BlockVector(const std::vector<size_type> &block_sizes);
162
167 BlockVector(const std::vector<IndexSet> &local_ranges,
168 const std::vector<IndexSet> &ghost_indices,
169 const MPI_Comm communicator);
170
174 BlockVector(const std::vector<IndexSet> &local_ranges,
175 const MPI_Comm communicator);
176
189 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
190 &partitioners,
191 const MPI_Comm &comm_sm = MPI_COMM_SELF);
192
201 ~BlockVector() = default;
202
209
222
233 template <class Number2>
236
242
243#ifdef DEAL_II_WITH_PETSC
253#endif
254
255#ifdef DEAL_II_WITH_TRILINOS
266#endif
267
281 void
282 reinit(const size_type num_blocks,
283 const size_type block_size = 0,
284 const bool omit_zeroing_entries = false);
285
305 void
306 reinit(const std::vector<size_type> &block_sizes,
307 const bool omit_zeroing_entries = false);
308
323 template <typename Number2>
324 void
326 const bool omit_zeroing_entries = false);
327
345 void
346 reinit(const std::vector<IndexSet> &local_ranges,
347 const std::vector<IndexSet> &ghost_indices,
348 const MPI_Comm communicator);
349
353 void
354 reinit(const std::vector<IndexSet> &local_ranges,
355 const MPI_Comm communicator);
356
370 void
372 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
373 &partitioners,
374 const MPI_Comm &comm_sm = MPI_COMM_SELF);
375
382 void
384 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
385 &partitioners,
386 const bool make_ghosted,
387 const MPI_Comm &comm_sm = MPI_COMM_SELF);
388
412 void
414
423 void
425
434 void
436
440 bool
442
446 void
447 set_ghost_state(const bool ghosted) const;
448
463 template <typename Number2>
464 void
466
471 template <typename OtherNumber>
472 void
473 add(const std::vector<size_type> &indices,
474 const ::Vector<OtherNumber> &values);
475
480 void
481 sadd(const Number s, const BlockVector<Number> &V);
482
488 bool
489 all_zero() const;
490
494 Number
495 mean_value() const;
496
502 lp_norm(const real_type p) const;
503
521 void
534 operator*=(const Number factor);
535
540 operator/=(const Number factor);
541
547
553
562 void
565 const VectorOperation::values operation,
566 const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
567 &communication_pattern = {});
568
574 VectorOperation::values operation,
575 std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
576 communication_pattern = {})
577 {
578 import_elements(V, operation, communication_pattern);
579 }
580
584 Number
586
603 template <typename FullMatrixType>
604 void
605 multivector_inner_product(FullMatrixType &matrix,
606 const BlockVector<Number> &V,
607 const bool symmetric = false) const;
608
625 template <typename FullMatrixType>
626 Number
627 multivector_inner_product_with_metric(const FullMatrixType &matrix,
628 const BlockVector<Number> &V,
629 const bool symmetric = false) const;
630
640 template <typename FullMatrixType>
641 void
643 const FullMatrixType &matrix,
644 const Number s = Number(0.),
645 const Number b = Number(1.)) const;
646
650 void
651 add(const Number a);
652
656 void
657 add(const Number a, const BlockVector<Number> &V);
658
662 void
663 add(const Number a,
664 const BlockVector<Number> &V,
665 const Number b,
666 const BlockVector<Number> &W);
667
672 void
673 add(const std::vector<size_type> &indices,
674 const std::vector<Number> &values);
675
680 void
681 sadd(const Number s, const Number a, const BlockVector<Number> &V);
682
688 void
689 scale(const BlockVector<Number> &scaling_factors);
690
694 void
695 equ(const Number a, const BlockVector<Number> &V);
696
702 l1_norm() const;
703
709 l2_norm() const;
710
715 norm_sqr() const;
716
722 linfty_norm() const;
723
744 Number
745 add_and_dot(const Number a,
746 const BlockVector<Number> &V,
747 const BlockVector<Number> &W);
748
753 virtual size_type
754 size() const override;
755
769
773 void
774 print(std::ostream &out,
775 const unsigned int precision = 3,
776 const bool scientific = true,
777 const bool across = true) const;
778
782 std::size_t
797
803 };
804
807 } // end of namespace distributed
808
809} // end of namespace LinearAlgebra
810
811
819template <typename Number>
820inline void
826
827
832template <typename Number>
833struct is_serial_vector<LinearAlgebra::distributed::BlockVector<Number>>
834 : std::false_type
835{};
836
837
839
840#ifdef DEAL_II_MSVC
841# include <deal.II/lac/la_parallel_block_vector.templates.h>
842#endif
843
844#endif
::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
void swap(LinearAlgebra::distributed::BlockVector< Number > &u, LinearAlgebra::distributed::BlockVector< Number > &v)
BlockVector< Number > & operator-=(const BlockVector< Number > &V)
void swap(BlockVector< Number > &v)
Number operator*(const BlockVector< Number > &V) const
BlockVector(const size_type num_blocks=0, const size_type block_size=0)
void add(const std::vector< size_type > &indices, const std::vector< Number > &values)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
void reinit(const std::vector< IndexSet > &local_ranges, const MPI_Comm communicator)
void sadd(const Number s, const Number a, const BlockVector< Number > &V)
static constexpr unsigned int communication_block_size
void add(const Number a, const BlockVector< Number > &V, const Number b, const BlockVector< Number > &W)
void reinit(const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &partitioners, const bool make_ghosted, const MPI_Comm &comm_sm=MPI_COMM_SELF)
BlockVector(const std::vector< size_type > &block_sizes)
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)
BlockVector< Number > & operator=(const TrilinosWrappers::MPI::BlockVector &trilinos_vec)
BlockVector< Number > & operator=(const PETScWrappers::MPI::BlockVector &petsc_vec)
BlockVector< Number > & operator/=(const Number factor)
void set_ghost_state(const bool ghosted) const
void add(const Number a, const BlockVector< Number > &V)
void add(const std::vector< size_type > &indices, const ::Vector< OtherNumber > &values)
BlockVector & operator=(const Vector< Number > &V)
BlockVector< Number > & operator*=(const Number factor)
BlockVector(const std::vector< IndexSet > &local_ranges, const MPI_Comm communicator)
void scale(const BlockVector< Number > &scaling_factors)
void reinit(const std::vector< IndexSet > &local_ranges, const std::vector< IndexSet > &ghost_indices, const MPI_Comm communicator)
Number multivector_inner_product_with_metric(const FullMatrixType &matrix, const BlockVector< Number > &V, const bool symmetric=false) const
typename BaseClass::const_reference const_reference
real_type lp_norm(const real_type p) const
BlockVector & operator=(const BlockVector< Number2 > &V)
void reinit(const BlockVector< Number2 > &V, const bool omit_zeroing_entries=false)
BlockVector & operator=(const BlockVector &V)
void reinit(const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &partitioners, const MPI_Comm &comm_sm=MPI_COMM_SELF)
Number add_and_dot(const Number a, const BlockVector< Number > &V, const BlockVector< Number > &W)
void equ(const Number a, const BlockVector< Number > &V)
void sadd(const Number s, const BlockVector< Number > &V)
void copy_locally_owned_data_from(const BlockVector< Number2 > &src)
BlockVector(const std::vector< IndexSet > &local_ranges, const std::vector< IndexSet > &ghost_indices, const MPI_Comm communicator)
BlockVector(const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &partitioners, const MPI_Comm &comm_sm=MPI_COMM_SELF)
typename BaseClass::const_pointer const_pointer
typename BaseClass::const_iterator const_iterator
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
void reinit(const std::vector< size_type > &block_sizes, const bool omit_zeroing_entries=false)
void compress(VectorOperation::values operation)
void import_elements(const LinearAlgebra::ReadWriteVector< Number > &V, const VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
BlockVector< Number > & operator+=(const BlockVector< Number > &V)
BlockVector & operator=(const value_type s)
virtual size_type size() const override
#define DEAL_II_DEPRECATED
Definition config.h:207
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define DeclException0(Exception0)
Definition exceptions.h:471
static ::ExceptionBase & ExcVectorTypeNotCompatible()
static ::ExceptionBase & ExcIteratorRangeDoesNotMatchVectorSize()