ginkgo/core/solver/gmres.hpp Source File

ginkgo/core/solver/gmres.hpp Source File#

Reference API: ginkgo/core/solver/gmres.hpp Source File
Reference API
gmres.hpp
1// SPDX-FileCopyrightText: 2017 - 2025 The Ginkgo authors
2//
3// SPDX-License-Identifier: BSD-3-Clause
4
5#ifndef GKO_PUBLIC_CORE_SOLVER_GMRES_HPP_
6#define GKO_PUBLIC_CORE_SOLVER_GMRES_HPP_
7
8
9#include <vector>
10
11#include <ginkgo/core/base/array.hpp>
12#include <ginkgo/core/base/exception_helpers.hpp>
13#include <ginkgo/core/base/lin_op.hpp>
14#include <ginkgo/core/base/math.hpp>
15#include <ginkgo/core/base/types.hpp>
16#include <ginkgo/core/config/config.hpp>
17#include <ginkgo/core/config/registry.hpp>
18#include <ginkgo/core/log/logger.hpp>
19#include <ginkgo/core/matrix/dense.hpp>
20#include <ginkgo/core/matrix/identity.hpp>
21#include <ginkgo/core/solver/solver_base.hpp>
22#include <ginkgo/core/stop/combined.hpp>
23#include <ginkgo/core/stop/criterion.hpp>
24
25
26namespace gko {
27namespace solver {
28
29
30[[deprecated]] constexpr size_type default_krylov_dim = 100u;
31
32constexpr size_type gmres_default_krylov_dim = 100u;
33
34namespace gmres {
38enum class ortho_method {
42 mgs,
46 cgs,
50 cgs2
51};
52
54std::ostream& operator<<(std::ostream& stream, ortho_method ortho);
55
56} // namespace gmres
57
70template <typename ValueType = default_precision>
71class Gmres
72 : public EnableLinOp<Gmres<ValueType>>,
73 public EnablePreconditionedIterativeSolver<ValueType, Gmres<ValueType>>,
74 public Transposable {
75 friend class EnableLinOp<Gmres>;
76 friend class EnablePolymorphicObject<Gmres, LinOp>;
77
78public:
79 using value_type = ValueType;
81
82 std::unique_ptr<LinOp> transpose() const override;
83
84 std::unique_ptr<LinOp> conj_transpose() const override;
85
91 bool apply_uses_initial_guess() const override { return true; }
92
98 size_type get_krylov_dim() const { return parameters_.krylov_dim; }
99
105 void set_krylov_dim(size_type other) { parameters_.krylov_dim = other; }
106
107
108 class Factory;
109
112 parameters_type, Factory> {
115
118
120 gmres::ortho_method GKO_FACTORY_PARAMETER_SCALAR(
121 ortho_method, gmres::ortho_method::mgs);
122 };
125
139 static parameters_type parse(const config::pnode& config,
140 const config::registry& context,
141 const config::type_descriptor& td_for_child =
142 config::make_type_descriptor<ValueType>());
143
144protected:
145 void apply_impl(const LinOp* b, LinOp* x) const override;
146
147 template <typename VectorType>
148 void apply_dense_impl(const VectorType* b, VectorType* x) const;
149
150 void apply_impl(const LinOp* alpha, const LinOp* b, const LinOp* beta,
151 LinOp* x) const override;
152
153 explicit Gmres(std::shared_ptr<const Executor> exec)
154 : EnableLinOp<Gmres>(std::move(exec))
155 {}
156
157 explicit Gmres(const Factory* factory,
158 std::shared_ptr<const LinOp> system_matrix)
159 : EnableLinOp<Gmres>(factory->get_executor(),
160 gko::transpose(system_matrix->get_size())),
161 EnablePreconditionedIterativeSolver<ValueType, Gmres<ValueType>>{
162 std::move(system_matrix), factory->get_parameters()},
163 parameters_{factory->get_parameters()}
164 {
165 if (!parameters_.krylov_dim) {
166 parameters_.krylov_dim = gmres_default_krylov_dim;
167 }
168 }
169};
170
171
172template <typename ValueType>
173struct workspace_traits<Gmres<ValueType>> {
174 using Solver = Gmres<ValueType>;
175 // number of vectors used by this workspace
176 static int num_vectors(const Solver&);
177 // number of arrays used by this workspace
178 static int num_arrays(const Solver&);
179 // array containing the num_vectors names for the workspace vectors
180 static std::vector<std::string> op_names(const Solver&);
181 // array containing the num_arrays names for the workspace vectors
182 static std::vector<std::string> array_names(const Solver&);
183 // array containing all varying scalar vectors (independent of problem size)
184 static std::vector<int> scalars(const Solver&);
185 // array containing all varying vectors (dependent on problem size)
186 static std::vector<int> vectors(const Solver&);
187
188 // residual vector
189 constexpr static int residual = 0;
190 // preconditioned vector
191 constexpr static int preconditioned_vector = 1;
192 // krylov basis multivector
193 constexpr static int krylov_bases = 2;
194 // hessenberg matrix
195 constexpr static int hessenberg = 3;
196 // auxiliary space for CGS2
197 constexpr static int hessenberg_aux = 4;
198 // givens sin parameters
199 constexpr static int givens_sin = 5;
200 // givens cos parameters
201 constexpr static int givens_cos = 6;
202 // coefficients of the residual in Krylov space
203 constexpr static int residual_norm_collection = 7;
204 // residual norm scalar
205 constexpr static int residual_norm = 8;
206 // solution of the least-squares problem in Krylov space
207 constexpr static int y = 9;
208 // solution of the least-squares problem mapped to the full space
209 constexpr static int before_preconditioner = 10;
210 // preconditioned solution of the least-squares problem
211 constexpr static int after_preconditioner = 11;
212 // constant 1.0 scalar
213 constexpr static int one = 12;
214 // constant -1.0 scalar
215 constexpr static int minus_one = 13;
216 // temporary norm vector of next_krylov to copy into hessenberg matrix
217 constexpr static int next_krylov_norm_tmp = 14;
218 // preconditioned krylov basis multivector
219 constexpr static int preconditioned_krylov_bases = 15;
220
221 // stopping status array
222 constexpr static int stop = 0;
223 // reduction tmp array
224 constexpr static int tmp = 1;
225 // final iteration number array
226 constexpr static int final_iter_nums = 2;
227};
228
229
230} // namespace solver
231} // namespace gko
232
233
234#endif // GKO_PUBLIC_CORE_SOLVER_GMRES_HPP_
Definition lin_op.hpp:878
Definition polymorphic_object.hpp:668
Definition lin_op.hpp:117
std::shared_ptr< const Executor > get_executor() const noexcept
Definition polymorphic_object.hpp:243
Definition lin_op.hpp:433
Definition property_tree.hpp:28
Definition registry.hpp:167
Definition type_descriptor.hpp:39
Definition gmres.hpp:123
Definition gmres.hpp:74
std::unique_ptr< LinOp > conj_transpose() const override
std::unique_ptr< LinOp > transpose() const override
size_type get_krylov_dim() const
Definition gmres.hpp:98
bool apply_uses_initial_guess() const override
Definition gmres.hpp:91
void set_krylov_dim(size_type other)
Definition gmres.hpp:105
static parameters_type parse(const config::pnode &config, const config::registry &context, const config::type_descriptor &td_for_child=config::make_type_descriptor< ValueType >())
#define GKO_FACTORY_PARAMETER_SCALAR(_name, _default)
Definition abstract_factory.hpp:445
#define GKO_ENABLE_BUILD_METHOD(_factory_name)
Definition abstract_factory.hpp:394
#define GKO_ENABLE_LIN_OP_FACTORY(_lin_op, _parameters_name, _factory_name)
Definition lin_op.hpp:1016
The Ginkgo namespace.
Definition abstract_factory.hpp:20
constexpr T one()
Definition math.hpp:630
std::size_t size_type
Definition types.hpp:89
STL namespace.
Definition gmres.hpp:112
size_type krylov_dim
Definition gmres.hpp:114
gmres::ortho_method ortho_method
Definition gmres.hpp:121
bool flexible
Definition gmres.hpp:117
Definition solver_base.hpp:238