23#include <ginkgo/ginkgo.hpp>
32#define MFEM_GINKGO_VERSION \
33 ((GKO_VERSION_MAJOR*100 + GKO_VERSION_MINOR)*100 + GKO_VERSION_PATCH)
42#if MFEM_GINKGO_VERSION < 10500
43template <
typename T>
using gko_array = gko::Array<T>;
45template <
typename T>
using gko_array = gko::array<T>;
76 gko::size_type size,
Vector *mfem_vec,
77 bool ownership =
false)
78 : gko::matrix::Dense<
real_t>(
80 gko::
dim<2> {size, 1},
84 exec != exec->get_master() ?
true : false)),
95 wrapped_vec = std::unique_ptr<
Vector,
96 std::function<void(
Vector *)>>(
101 using deleter = gko::null_deleter<Vector>;
102 wrapped_vec = std::unique_ptr<Vector,
103 std::function<void(Vector *)>>(
104 mfem_vec, deleter{});
108 static std::unique_ptr<VectorWrapper>
create(
109 std::shared_ptr<const gko::Executor> exec,
112 bool ownership =
false)
114 return std::unique_ptr<VectorWrapper>(
126 virtual std::unique_ptr<gko::matrix::Dense<real_t>>
131 this->wrapped_vec.get()->GetMemory().GetMemoryType());
133 mfem_vec->
UseDevice(this->wrapped_vec.get()->UseDevice());
149 std::shared_ptr<const gko::Executor> exec,
150 const gko::dim<2> &size,
151 gko::size_type stride)
const override
158 "VectorWrapper cannot be created with stride > 1");
161 gko::size_type total_size = size[0]*size[1];
164 this->wrapped_vec.get()->GetMemory().GetMemoryType());
166 mfem_vec->
UseDevice(this->wrapped_vec.get()->UseDevice());
172 this->get_executor(), total_size, mfem_vec,
179 const gko::span &rows,
180 const gko::span &columns,
181 const gko::size_type stride)
override
184 gko::size_type num_rows = rows.end - rows.begin;
185 gko::size_type num_cols = columns.end - columns.begin;
189 if (num_cols > 1 || stride > 1)
191 throw gko::BadDimension(
192 __FILE__, __LINE__, __func__,
"new_submatrix", num_rows,
194 "VectorWrapper submatrix must have one column and stride = 1");
196 int data_size =
static_cast<int>(num_rows * num_cols);
197 int start =
static_cast<int>(rows.begin);
200 mfem_vec->
MakeRef(*(this->wrapped_vec.get()), start, data_size);
201 mfem_vec->
UseDevice(this->wrapped_vec.get()->UseDevice());
208 this->get_executor(), data_size, mfem_vec,
213 std::unique_ptr<
Vector, std::function<void(
Vector *)>> wrapped_vec;
225 :
public gko::EnableLinOp<OperatorWrapper>,
226 public gko::EnableCreateMethod<OperatorWrapper>
230 gko::size_type size = 0,
233 gko::EnableCreateMethod<OperatorWrapper>()
235 this->wrapped_oper = oper;
239 void apply_impl(
const gko::LinOp *
b, gko::LinOp *x)
const override;
241 const gko::LinOp *
beta, gko::LinOp *x)
const override;
249template <
typename ValueType=real_t>
253 auto cpu_norm = clone(norm->get_executor()->get_master(), norm);
255 return cpu_norm->at(0, 0);
260template <
typename ValueType=real_t>
264 auto exec =
b->get_executor();
266 auto b_norm = gko::initialize<gko::matrix::Dense<ValueType>>({0.0}, exec);
268#if MFEM_GINKGO_VERSION < 10600
269 b->compute_norm2(gko::lend(b_norm));
271 b->compute_norm2(b_norm);
274 return std::pow(
get_norm(b_norm.get()),2);
288template <
typename ValueType=real_t>
295 if (compute_real_residual)
297 mfem::out <<
"Iteration log with real residual norms:" << std::endl;
301 mfem::out <<
"Iteration log with residual norms:" << std::endl;
303 mfem::out <<
'|' << std::setw(10) <<
"Iteration" <<
'|' << std::setw(25)
304 <<
"Residual Norm" <<
'|' << std::endl;
307 mfem::out <<
'|' << std::setfill(
'-') << std::setw(11) <<
'|' <<
308 std::setw(26) <<
'|' << std::setfill(
' ') << std::endl;
311 for (std::size_t i = 0; i < iterations.size(); i++)
313 mfem::out <<
'|' << std::setw(10) << iterations[i] <<
'|'
314 << std::setw(25) << residual_norms[i] <<
'|' << std::endl;
318 mfem::out.unsetf(std::ios_base::floatfield);
320 mfem::out <<
'|' << std::setfill(
'-') << std::setw(11) <<
'|' <<
321 std::setw(26) <<
'|' << std::setfill(
' ') << std::endl;
328 const gko::size_type &iteration,
329 const gko::LinOp *residual,
330 const gko::LinOp *solution,
331 const gko::LinOp *residual_norm)
const override
333 iteration_complete_core(iteration, residual, solution, residual_norm,
338 const gko::size_type &iteration,
339 const gko::LinOp *residual,
340 const gko::LinOp *solution,
341 const gko::LinOp *residual_norm,
342 const gko::LinOp *implicit_sq_residual_norm)
const override
344 iteration_complete_core(iteration, residual, solution, residual_norm,
345 implicit_sq_residual_norm);
347#if MFEM_GINKGO_VERSION > 10500
350 const gko::LinOp *rhs,
351 const gko::LinOp *solution,
352 const gko::size_type &iteration,
353 const gko::LinOp *residual,
354 const gko::LinOp *residual_norm,
355 const gko::LinOp *implicit_sq_residual_norm,
356 const gko::array<gko::stopping_status>* status,
357 bool stopped)
const override
359 iteration_complete_core(iteration, residual, solution, residual_norm,
360 implicit_sq_residual_norm);
366 const gko::LinOp *matrix,
const gko_dense *b,
367 bool compute_real_residual=
false)
369#if MFEM_GINKGO_VERSION < 10500
370 gko::log::Logger(exec,
371 gko::log::Logger::iteration_complete_mask),
373 gko::log::Logger(gko::log::Logger::iteration_complete_mask),
377 compute_real_residual{compute_real_residual}
379 if (compute_real_residual ==
true)
384 res = std::move(gko_dense::create_with_config_of(b_cast).release());
388 res = std::move(gko::clone(
b).release());
396 void iteration_complete_core(
const gko::size_type &iteration,
397 const gko::LinOp *residual,
398 const gko::LinOp *solution,
399 const gko::LinOp *residual_norm,
400 const gko::LinOp *implicit_sq_residual_norm)
const
404 if (solution && compute_real_residual)
407 auto exec = matrix->get_executor();
410#if MFEM_GINKGO_VERSION < 10600
411 matrix->apply(gko::lend(solution), gko::lend(res));
413 matrix->apply(solution, res);
419 const VectorWrapper *b_cast = gko::as<const VectorWrapper>(b);
427 auto neg_one = gko::initialize<gko_dense>({-1.0}, exec);
428#if MFEM_GINKGO_VERSION < 10600
429 res->add_scaled(gko::lend(neg_one), gko::lend(b));
431 res->add_scaled(neg_one, b);
437#if MFEM_GINKGO_VERSION < 10600
446 if (implicit_sq_residual_norm)
448 auto dense_norm = gko::as<gko_dense>(implicit_sq_residual_norm);
450 residual_norms.push_back(
get_norm(dense_norm));
453 else if (residual_norm)
455 auto dense_norm = gko::as<gko_dense>(residual_norm);
457 residual_norms.push_back(
get_norm(dense_norm));
462 auto dense_residual = gko::as<gko_dense>(residual);
464#if MFEM_GINKGO_VERSION < 10600
470 residual_norms.push_back(norm);
474 iterations.push_back(iteration);
478 const gko::LinOp *matrix;
484 mutable std::vector<ValueType> residual_norms{};
486 mutable std::vector<std::size_t> iterations{};
489 const bool compute_real_residual;
559 return this->executor;
563 std::shared_ptr<gko::Executor> executor;
766 std::shared_ptr<gko::stop::ResidualNorm<>::Factory>
773 std::shared_ptr<gko::stop::ResidualNorm<>::Factory>
780 std::shared_ptr<gko::stop::ImplicitResidualNorm<>::Factory>
787 std::shared_ptr<gko::stop::ImplicitResidualNorm<>::Factory>
838 initialize_ginkgo_log(gko::matrix::Dense<real_t>*
b)
const;
844 std::shared_ptr<gko::LinOp> system_oper;
854template<
typename SolverType>
865 auto current_params = gko::as<typename SolverType::Factory>
879 auto current_params = gko::as<typename SolverType::Factory>
893 auto current_params = gko::as<typename SolverType::Factory>
1059using gko::solver::cb_gmres::storage_precision;
1090 storage_precision prec = storage_precision::reduce1);
1109 storage_precision prec = storage_precision::reduce1);
1170 const std::string &storage_opt =
"none",
1171 const real_t accuracy = 1.e-1,
1172 const int max_block_size = 32
1204 const std::string &factorization_type =
"exact",
1205 const int sweeps = 0,
1206 const bool skip_sort =
false
1246 const std::string &factorization_type =
"exact",
1247 const int sweeps = 0,
1248 const int sparsity_power = 1,
1249 const bool skip_sort =
false
1281 const std::string &factorization_type =
"exact",
1282 const int sweeps = 0,
1283 const bool skip_sort =
false
1323 const std::string &factorization_type =
"exact",
1324 const int sweeps = 0,
1325 const int sparsity_power = 1,
1326 const bool skip_sort =
false
1346 const Solver &mfem_precond
1356 MFEM_ABORT(
"Ginkgo::MFEMPreconditioner must be constructed "
1357 "with the MFEM Operator that it will wrap as an argument;\n"
1358 "calling SetOperator() is not allowed.");
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
BICGSTABSolver(GinkgoExecutor &exec)
CBGMRESSolver(GinkgoExecutor &exec, int dim=0, storage_precision prec=storage_precision::reduce1)
CGSSolver(GinkgoExecutor &exec)
CGSolver(GinkgoExecutor &exec)
void SetRelTol(real_t rtol)
void SetAbsTol(real_t atol)
void SetMaxIter(int max_it)
EnableGinkgoSolver(GinkgoExecutor &exec, bool use_implicit_res_norm)
FCGSolver(GinkgoExecutor &exec)
GMRESSolver(GinkgoExecutor &exec, int dim=0)
std::shared_ptr< gko::Executor > GetExecutor() const
GinkgoExecutor(ExecType exec_type)
virtual ~GinkgoExecutor()=default
@ REFERENCE
Reference CPU Executor.
@ OMP
OpenMP CPU Executor.
bool sub_op_needs_wrapped_vecs
bool UsesVectorWrappers() const
std::shared_ptr< gko::stop::ImplicitResidualNorm<>::Factory > imp_abs_criterion
virtual void Mult(const Vector &x, Vector &y) const
std::shared_ptr< ResidualLogger<> > residual_logger
void update_stop_factory()
bool use_implicit_res_norm
std::shared_ptr< gko::stop::ResidualNorm<>::Factory > abs_criterion
std::shared_ptr< gko::stop::Combined::Factory > combined_factory
real_t GetFinalNorm() const
int GetNumIterations() const
const std::shared_ptr< gko::LinOpFactory > GetFactory() const
std::shared_ptr< gko::log::Convergence<> > convergence_logger
GinkgoIterativeSolver(GinkgoExecutor &exec, bool use_implicit_res_norm)
void SetPrintLevel(int print_lvl)
std::shared_ptr< gko::LinOpFactory > solver_gen
std::shared_ptr< gko::stop::ResidualNorm<>::Factory > rel_criterion
virtual ~GinkgoIterativeSolver()=default
std::shared_ptr< gko::stop::ImplicitResidualNorm<>::Factory > imp_rel_criterion
std::shared_ptr< gko::Executor > executor
virtual void SetOperator(const Operator &op)
std::shared_ptr< gko::LinOp > solver
std::shared_ptr< gko::Executor > executor
const std::shared_ptr< gko::LinOpFactory > GetFactory() const
virtual ~GinkgoPreconditioner()=default
std::shared_ptr< gko::LinOp > generated_precond
virtual void SetOperator(const Operator &op)
virtual void Mult(const Vector &x, Vector &y) const
bool HasGeneratedPreconditioner() const
bool has_generated_precond
GinkgoPreconditioner(GinkgoExecutor &exec)
std::shared_ptr< gko::LinOpFactory > precond_gen
const std::shared_ptr< gko::LinOp > GetGeneratedPreconditioner() const
IRSolver(GinkgoExecutor &exec)
IcIsaiPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const int sparsity_power=1, const bool skip_sort=false)
IcPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const bool skip_sort=false)
IluIsaiPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const int sparsity_power=1, const bool skip_sort=false)
IluPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const bool skip_sort=false)
JacobiPreconditioner(GinkgoExecutor &exec, const std::string &storage_opt="none", const real_t accuracy=1.e-1, const int max_block_size=32)
MFEMPreconditioner(GinkgoExecutor &exec, const Solver &mfem_precond)
virtual void SetOperator(const Operator &op)
OperatorWrapper(std::shared_ptr< const gko::Executor > exec, gko::size_type size=0, const Operator *oper=NULL)
void apply_impl(const gko::LinOp *b, gko::LinOp *x) const override
static std::unique_ptr< VectorWrapper > create(std::shared_ptr< const gko::Executor > exec, gko::size_type size, Vector *mfem_vec, bool ownership=false)
Vector & get_mfem_vec_ref()
virtual std::unique_ptr< gko::matrix::Dense< real_t > > create_with_type_of_impl(std::shared_ptr< const gko::Executor > exec, const gko::dim< 2 > &size, gko::size_type stride) const override
VectorWrapper(std::shared_ptr< const gko::Executor > exec, gko::size_type size, Vector *mfem_vec, bool ownership=false)
virtual std::unique_ptr< gko::matrix::Dense< real_t > > create_with_same_config() const override
const Vector & get_mfem_vec_const_ref() const
virtual std::unique_ptr< gko::matrix::Dense< real_t > > create_submatrix_impl(const gko::span &rows, const gko::span &columns, const gko::size_type stride) override
void operator()(pointer ptr) const noexcept
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
real_t get_norm(const gko::matrix::Dense< ValueType > *norm)
real_t compute_norm(const gko::matrix::Dense< ValueType > *b)
gko::Array< T > gko_array
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void on_iteration_complete(const gko::LinOp *op, const gko::size_type &iteration, const gko::LinOp *residual, const gko::LinOp *solution, const gko::LinOp *residual_norm) const override
void on_iteration_complete(const gko::LinOp *op, const gko::LinOp *rhs, const gko::LinOp *solution, const gko::size_type &iteration, const gko::LinOp *residual, const gko::LinOp *residual_norm, const gko::LinOp *implicit_sq_residual_norm, const gko::array< gko::stopping_status > *status, bool stopped) const override
ResidualLogger(std::shared_ptr< const gko::Executor > exec, const gko::LinOp *matrix, const gko_dense *b, bool compute_real_residual=false)
gko::matrix::Dense< ValueType > gko_dense
void on_iteration_complete(const gko::LinOp *op, const gko::size_type &iteration, const gko::LinOp *residual, const gko::LinOp *solution, const gko::LinOp *residual_norm, const gko::LinOp *implicit_sq_residual_norm) const override