ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
diagonal_preconditioner.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2022 - 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
10#include "state_vector.h"
11
12#include <deal.II/lac/la_parallel_block_vector.h>
13
14namespace ryujin
15{
22 template <typename Number>
24 {
25 public:
30
35
40
44 void reinit(const std::shared_ptr<const dealii::Utilities::MPI::Partitioner>
45 &scalar_partitioner)
46 {
47 diagonal_.reinit(scalar_partitioner);
48 }
49
54 {
55 return diagonal_;
56 }
57
61 void vmult(ScalarVector &dst, const ScalarVector &src) const
62 {
63 const auto n_owned = diagonal_.get_partitioner()->locally_owned_size();
64 AssertDimension(n_owned, src.get_partitioner()->locally_owned_size());
65 AssertDimension(n_owned, dst.get_partitioner()->locally_owned_size());
66
67 DEAL_II_OPENMP_SIMD_PRAGMA
68 for (unsigned int i = 0; i < n_owned; ++i)
69 dst.local_element(i) =
70 diagonal_.local_element(i) * src.local_element(i);
71 }
72
76 void vmult(BlockVector &dst, const BlockVector &src) const
77 {
78 const auto n_blocks = src.n_blocks();
79 AssertDimension(n_blocks, dst.n_blocks());
80
81 const auto n_owned = diagonal_.get_partitioner()->locally_owned_size();
82
83 for (unsigned int d = 0; d < n_blocks; ++d) {
84 AssertDimension(n_owned,
85 src.block(d).get_partitioner()->locally_owned_size());
86 AssertDimension(n_owned,
87 dst.block(d).get_partitioner()->locally_owned_size());
88
89 DEAL_II_OPENMP_SIMD_PRAGMA
90 for (unsigned int i = 0; i < n_owned; ++i)
91 dst.block(d).local_element(i) =
92 diagonal_.local_element(i) * src.block(d).local_element(i);
93 }
94 }
95
96 private:
97 ScalarVector diagonal_;
98 };
99
100} /* namespace ryujin */
void reinit(const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > &scalar_partitioner)
typename Vectors::ScalarVector< Number > ScalarVector
void vmult(BlockVector &dst, const BlockVector &src) const
void vmult(ScalarVector &dst, const ScalarVector &src) const
typename Vectors::BlockVector< Number > BlockVector
dealii::LinearAlgebra::distributed::Vector< Number > ScalarVector
Definition: state_vector.h:33
dealii::LinearAlgebra::distributed::BlockVector< Number > BlockVector
Definition: state_vector.h:39