5#ifndef GKO_PUBLIC_CORE_SOLVER_GMRES_HPP_
6#define GKO_PUBLIC_CORE_SOLVER_GMRES_HPP_
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>
30[[deprecated]]
constexpr size_type default_krylov_dim = 100u;
32constexpr size_type gmres_default_krylov_dim = 100u;
38enum class ortho_method {
54std::ostream& operator<<(std::ostream& stream, ortho_method ortho);
71template <
typename ValueType = default_precision>
74 public EnablePreconditionedIterativeSolver<ValueType, Gmres<ValueType>>,
78 GKO_ASSERT_SUPPORTED_VALUE_TYPE;
81 using value_type = ValueType;
82 using transposed_type = Gmres<ValueType>;
114 parameters_type, Factory> {
144 config::make_type_descriptor<ValueType>());
147 void apply_impl(
const LinOp* b,
LinOp* x)
const override;
149 template <
typename VectorType>
150 void apply_dense_impl(
const VectorType* b, VectorType* x)
const;
153 LinOp* x)
const override;
155 explicit Gmres(std::shared_ptr<const Executor> exec)
160 std::shared_ptr<const LinOp> system_matrix)
164 std::move(system_matrix), factory->get_parameters()},
165 parameters_{factory->get_parameters()}
168 parameters_.
krylov_dim = gmres_default_krylov_dim;
174template <
typename ValueType>
178 static int num_vectors(
const Solver&);
180 static int num_arrays(
const Solver&);
182 static std::vector<std::string> op_names(
const Solver&);
184 static std::vector<std::string> array_names(
const Solver&);
186 static std::vector<int> scalars(
const Solver&);
188 static std::vector<int> vectors(
const Solver&);
191 constexpr static int residual = 0;
193 constexpr static int preconditioned_vector = 1;
195 constexpr static int krylov_bases = 2;
197 constexpr static int hessenberg = 3;
199 constexpr static int hessenberg_aux = 4;
201 constexpr static int givens_sin = 5;
203 constexpr static int givens_cos = 6;
205 constexpr static int residual_norm_collection = 7;
207 constexpr static int residual_norm = 8;
209 constexpr static int y = 9;
211 constexpr static int before_preconditioner = 10;
213 constexpr static int after_preconditioner = 11;
215 constexpr static int one = 12;
217 constexpr static int minus_one = 13;
219 constexpr static int next_krylov_norm_tmp = 14;
221 constexpr static int preconditioned_krylov_bases = 15;
224 constexpr static int stop = 0;
226 constexpr static int tmp = 1;
228 constexpr static int final_iter_nums = 2;
The EnableLinOp mixin can be used to provide sensible default implementations of the majority of the ...
Definition lin_op.hpp:879
This mixin inherits from (a subclass of) PolymorphicObject and provides a base implementation of a ne...
Definition polymorphic_object.hpp:668
Definition lin_op.hpp:117
Linear operators which support transposition should implement the Transposable interface.
Definition lin_op.hpp:433
pnode describes a tree of properties.
Definition property_tree.hpp:28
This class stores additional context for creating Ginkgo objects from configuration files.
Definition registry.hpp:167
This class describes the value and index types to be used when building a Ginkgo type from a configur...
Definition type_descriptor.hpp:39
A LinOp implementing this interface stores a system matrix and stopping criterion factory.
Definition solver_base.hpp:802
GMRES or the generalized minimal residual method is an iterative type Krylov subspace method which is...
Definition gmres.hpp:75
std::unique_ptr< LinOp > conj_transpose() const override
Returns a LinOp representing the conjugate transpose of the Transposable object.
std::unique_ptr< LinOp > transpose() const override
Returns a LinOp representing the transpose of the Transposable object.
size_type get_krylov_dim() const
Gets the Krylov dimension of the solver.
Definition gmres.hpp:100
bool apply_uses_initial_guess() const override
Return true as iterative solvers use the data in x as an initial guess.
Definition gmres.hpp:93
void set_krylov_dim(size_type other)
Sets the Krylov dimension.
Definition gmres.hpp:107
static parameters_type parse(const config::pnode &config, const config::registry &context, const config::type_descriptor &td_for_child=config::make_type_descriptor< ValueType >())
Create the parameters from the property_tree.
#define GKO_FACTORY_PARAMETER_SCALAR(_name, _default)
Creates a scalar factory parameter in the factory parameters structure.
Definition abstract_factory.hpp:445
#define GKO_ENABLE_BUILD_METHOD(_factory_name)
Defines a build method for the factory, simplifying its construction by removing the repetitive typin...
Definition abstract_factory.hpp:394
#define GKO_ENABLE_LIN_OP_FACTORY(_lin_op, _parameters_name, _factory_name)
This macro will generate a default implementation of a LinOpFactory for the LinOp subclass it is defi...
Definition lin_op.hpp:1017
The ginkgo Solve namespace.
Definition bicg.hpp:28
The Ginkgo namespace.
Definition abstract_factory.hpp:20
std::size_t size_type
Integral type used for allocation quantities.
Definition types.hpp:90
batch_dim< 2, DimensionType > transpose(const batch_dim< 2, DimensionType > &input)
Returns a batch_dim object with its dimensions swapped for batched operators.
Definition batch_dim.hpp:119
size_type krylov_dim
Krylov subspace dimension/restart value.
Definition gmres.hpp:116
gmres::ortho_method ortho_method
Orthogonalization method.
Definition gmres.hpp:123
bool flexible
Flexible GMRES.
Definition gmres.hpp:119
Definition solver_base.hpp:855
Traits class providing information on the type and location of workspace vectors inside a solver.
Definition solver_base.hpp:238