MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ginkgo.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_GINKGO
13#define MFEM_GINKGO
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_GINKGO
18
19#include "operator.hpp"
20#include "sparsemat.hpp"
21#include "solvers.hpp"
22
23#include <ginkgo/ginkgo.hpp>
24
25#include <iomanip>
26#include <ios>
27#include <string>
28#include <vector>
29#include <fstream>
30#include <iostream>
31
32#define MFEM_GINKGO_VERSION \
33 ((GKO_VERSION_MAJOR*100 + GKO_VERSION_MINOR)*100 + GKO_VERSION_PATCH)
34
35namespace mfem
36{
37namespace Ginkgo
38{
39
40/// The alias 'gko_array' refers to 'gko::Array' (Ginkgo < 1.5.0) or
41/// 'gko::array' (Ginkgo >= 1.5.0).
42#if MFEM_GINKGO_VERSION < 10500
43template <typename T> using gko_array = gko::Array<T>;
44#else
45template <typename T> using gko_array = gko::array<T>;
46#endif
47
48/**
49* Helper class for a case where a wrapped MFEM Vector
50* should be owned by Ginkgo, and deleted when the wrapper
51* object goes out of scope.
52*/
53template <typename T>
55{
56public:
57 using pointer = T *;
58
59 // Destroys an MFEM object. Requires object to have a Destroy() method.
60 void operator()(pointer ptr) const noexcept { ptr->Destroy(); }
61};
62
63/**
64* This class wraps an MFEM vector object for Ginkgo's use. It
65* allows Ginkgo and MFEM to operate directly on the same
66* data, and is necessary to use MFEM Operators with Ginkgo
67* solvers.
68*
69* @ingroup Ginkgo
70*/
71
72class VectorWrapper : public gko::matrix::Dense<real_t>
73{
74public:
75 VectorWrapper(std::shared_ptr<const gko::Executor> exec,
76 gko::size_type size, Vector *mfem_vec,
77 bool ownership = false)
78 : gko::matrix::Dense<real_t>(
79 exec,
80 gko::dim<2> {size, 1},
82 size,
83 mfem_vec->ReadWrite(
84 exec != exec->get_master() ? true : false)),
85 1)
86 {
87 // This controls whether or not we want Ginkgo to own its MFEM Vector.
88 // Normally, when we are wrapping an MFEM Vector created outside
89 // Ginkgo, we do not want ownership to be true. However, Ginkgo
90 // creates its own temporary vectors as part of its solvers, and
91 // these will be owned (and deleted) by Ginkgo.
92 if (ownership)
93 {
94 using deleter = gko_mfem_destroy<Vector>;
95 wrapped_vec = std::unique_ptr<Vector,
96 std::function<void(Vector *)>>(
97 mfem_vec, deleter{});
98 }
99 else
100 {
101 using deleter = gko::null_deleter<Vector>;
102 wrapped_vec = std::unique_ptr<Vector,
103 std::function<void(Vector *)>>(
104 mfem_vec, deleter{});
105 }
106 }
107
108 static std::unique_ptr<VectorWrapper> create(
109 std::shared_ptr<const gko::Executor> exec,
110 gko::size_type size,
111 Vector *mfem_vec,
112 bool ownership = false)
113 {
114 return std::unique_ptr<VectorWrapper>(
115 new VectorWrapper(exec, size, mfem_vec, ownership));
116 }
117
118 // Return reference to MFEM Vector object
119 Vector &get_mfem_vec_ref() { return *(this->wrapped_vec.get()); }
120
121 // Return const reference to MFEM Vector object
122 const Vector &get_mfem_vec_const_ref() const { return *(this->wrapped_vec.get()); }
123
124 // Override base Dense class implementation for creating new vectors
125 // with same executor and size as self
126 virtual std::unique_ptr<gko::matrix::Dense<real_t>>
127 create_with_same_config() const override
128 {
129 Vector *mfem_vec = new Vector(
130 this->get_size()[0],
131 this->wrapped_vec.get()->GetMemory().GetMemoryType());
132
133 mfem_vec->UseDevice(this->wrapped_vec.get()->UseDevice());
134
135 // If this function is called, Ginkgo is creating this
136 // object and should control the memory, so ownership is
137 // set to true
138 return VectorWrapper::create(this->get_executor(),
139 this->get_size()[0],
140 mfem_vec,
141 true);
142 }
143
144 // Override base Dense class implementation for creating new vectors
145 // with same executor and type as self, but with a different size.
146 // This function will create "one large VectorWrapper" of size
147 // size[0] * size[1], since MFEM Vectors only have one dimension.
148 virtual std::unique_ptr<gko::matrix::Dense<real_t>> create_with_type_of_impl(
149 std::shared_ptr<const gko::Executor> exec,
150 const gko::dim<2> &size,
151 gko::size_type stride) const override
152 {
153 // Only stride of 1 is allowed for VectorWrapper type
154 if (stride > 1)
155 {
156 throw gko::Error(
157 __FILE__, __LINE__,
158 "VectorWrapper cannot be created with stride > 1");
159 }
160 // Compute total size of new Vector
161 gko::size_type total_size = size[0]*size[1];
162 Vector *mfem_vec = new Vector(
163 total_size,
164 this->wrapped_vec.get()->GetMemory().GetMemoryType());
165
166 mfem_vec->UseDevice(this->wrapped_vec.get()->UseDevice());
167
168 // If this function is called, Ginkgo is creating this
169 // object and should control the memory, so ownership is
170 // set to true
172 this->get_executor(), total_size, mfem_vec,
173 true);
174 }
175
176 // Override base Dense class implementation for creating new sub-vectors
177 // from a larger vector.
178 virtual std::unique_ptr<gko::matrix::Dense<real_t>> create_submatrix_impl(
179 const gko::span &rows,
180 const gko::span &columns,
181 const gko::size_type stride) override
182 {
183
184 gko::size_type num_rows = rows.end - rows.begin;
185 gko::size_type num_cols = columns.end - columns.begin;
186 // Data in the Dense matrix will be stored in row-major format.
187 // Check that we only have one column, and that the stride = 1
188 // (only allowed value for VectorWrappers).
189 if (num_cols > 1 || stride > 1)
190 {
191 throw gko::BadDimension(
192 __FILE__, __LINE__, __func__, "new_submatrix", num_rows,
193 num_cols,
194 "VectorWrapper submatrix must have one column and stride = 1");
195 }
196 int data_size = static_cast<int>(num_rows * num_cols);
197 int start = static_cast<int>(rows.begin);
198 // Create a new MFEM Vector pointing to this starting point in the data
199 Vector *mfem_vec = new Vector();
200 mfem_vec->MakeRef(*(this->wrapped_vec.get()), start, data_size);
201 mfem_vec->UseDevice(this->wrapped_vec.get()->UseDevice());
202
203 // If this function is called, Ginkgo is creating this
204 // object and should control the memory, so ownership is
205 // set to true (but MFEM doesn't own and won't delete
206 // the data, at it's only a reference to the parent Vector)
208 this->get_executor(), data_size, mfem_vec,
209 true);
210 }
211
212private:
213 std::unique_ptr<Vector, std::function<void(Vector *)>> wrapped_vec;
214};
215
216/**
217* This class wraps an MFEM Operator for Ginkgo, to make its Mult()
218* function available to Ginkgo, provided the input and output vectors
219* are of the VectorWrapper type.
220* Note that this class does NOT take ownership of the MFEM Operator.
221*
222* @ingroup Ginkgo
223*/
225 : public gko::EnableLinOp<OperatorWrapper>,
226 public gko::EnableCreateMethod<OperatorWrapper>
227{
228public:
229 OperatorWrapper(std::shared_ptr<const gko::Executor> exec,
230 gko::size_type size = 0,
231 const Operator *oper = NULL)
232 : gko::EnableLinOp<OperatorWrapper>(exec, gko::dim<2> {size, size}),
233 gko::EnableCreateMethod<OperatorWrapper>()
234 {
235 this->wrapped_oper = oper;
236 }
237
238protected:
239 void apply_impl(const gko::LinOp *b, gko::LinOp *x) const override;
240 void apply_impl(const gko::LinOp *alpha, const gko::LinOp *b,
241 const gko::LinOp *beta, gko::LinOp *x) const override;
242
243private:
244 const Operator *wrapped_oper;
245};
246
247// Utility function which gets the scalar value of a Ginkgo gko::matrix::Dense
248// matrix representing the norm of a vector.
249template <typename ValueType=real_t>
250real_t get_norm(const gko::matrix::Dense<ValueType> *norm)
251{
252 // Put the value on CPU thanks to the master executor
253 auto cpu_norm = clone(norm->get_executor()->get_master(), norm);
254 // Return the scalar value contained at position (0, 0)
255 return cpu_norm->at(0, 0);
256}
257
258// Utility function which computes the norm of a Ginkgo gko::matrix::Dense
259// vector.
260template <typename ValueType=real_t>
261real_t compute_norm(const gko::matrix::Dense<ValueType> *b)
262{
263 // Get the executor of the vector
264 auto exec = b->get_executor();
265 // Initialize a result scalar containing the value 0.0.
266 auto b_norm = gko::initialize<gko::matrix::Dense<ValueType>>({0.0}, exec);
267 // Use the dense `compute_norm2` function to compute the norm.
268#if MFEM_GINKGO_VERSION < 10600
269 b->compute_norm2(gko::lend(b_norm));
270#else
271 b->compute_norm2(b_norm);
272#endif
273 // Use the other utility function to return the norm contained in `b_norm``
274 return std::pow(get_norm(b_norm.get()),2);
275}
276
277/**
278 * Custom logger class which intercepts the residual norm scalar and solution
279 * vector in order to print a table of real vs recurrent (internal to the
280 * solvers) residual norms.
281 *
282 * This has been taken from the custom-logger example of Ginkgo. See the
283 * custom-logger example to understand how to write and modify your own loggers
284 * with Ginkgo.
285 *
286 * @ingroup Ginkgo
287 */
288template <typename ValueType=real_t>
289struct ResidualLogger : gko::log::Logger
290{
291 // Output the logger's data in a table format
292 void write() const
293 {
294 // Print a header for the table
295 if (compute_real_residual)
296 {
297 mfem::out << "Iteration log with real residual norms:" << std::endl;
298 }
299 else
300 {
301 mfem::out << "Iteration log with residual norms:" << std::endl;
302 }
303 mfem::out << '|' << std::setw(10) << "Iteration" << '|' << std::setw(25)
304 << "Residual Norm" << '|' << std::endl;
305 // Print a separation line. Note that for creating `10` characters
306 // `std::setw()` should be set to `11`.
307 mfem::out << '|' << std::setfill('-') << std::setw(11) << '|' <<
308 std::setw(26) << '|' << std::setfill(' ') << std::endl;
309 // Print the data one by one in the form
310 mfem::out << std::scientific;
311 for (std::size_t i = 0; i < iterations.size(); i++)
312 {
313 mfem::out << '|' << std::setw(10) << iterations[i] << '|'
314 << std::setw(25) << residual_norms[i] << '|' << std::endl;
315 }
316 // std::defaultfloat could be used here but some compilers do not support
317 // it properly, e.g. the Intel compiler
318 mfem::out.unsetf(std::ios_base::floatfield);
319 // Print a separation line
320 mfem::out << '|' << std::setfill('-') << std::setw(11) << '|' <<
321 std::setw(26) << '|' << std::setfill(' ') << std::endl;
322 }
323
324 using gko_dense = gko::matrix::Dense<ValueType>;
325
326 // Ginkgo 1.5 and older: version for solver that doesn't log implicit res norm
327 void on_iteration_complete(const gko::LinOp *op,
328 const gko::size_type &iteration,
329 const gko::LinOp *residual,
330 const gko::LinOp *solution,
331 const gko::LinOp *residual_norm) const override
332 {
333 iteration_complete_core(iteration, residual, solution, residual_norm,
334 nullptr);
335 }
336 // Ginkgo 1.5 and older: version with implicit residual norm
337 void on_iteration_complete(const gko::LinOp *op,
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
343 {
344 iteration_complete_core(iteration, residual, solution, residual_norm,
345 implicit_sq_residual_norm);
346 }
347#if MFEM_GINKGO_VERSION > 10500
348 // Ginkgo 1.6 and newer
349 void on_iteration_complete(const gko::LinOp *op,
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
358 {
359 iteration_complete_core(iteration, residual, solution, residual_norm,
360 implicit_sq_residual_norm);
361 }
362#endif
363
364 // Construct the logger and store the system matrix and b vectors
365 ResidualLogger(std::shared_ptr<const gko::Executor> exec,
366 const gko::LinOp *matrix, const gko_dense *b,
367 bool compute_real_residual=false)
368 :
369#if MFEM_GINKGO_VERSION < 10500
370 gko::log::Logger(exec,
371 gko::log::Logger::iteration_complete_mask),
372#else
373 gko::log::Logger(gko::log::Logger::iteration_complete_mask),
374#endif
375 matrix {matrix},
376 b{b},
377 compute_real_residual{compute_real_residual}
378 {
379 if (compute_real_residual == true)
380 {
381 if (dynamic_cast<const VectorWrapper*>(b))
382 {
383 const VectorWrapper *b_cast = gko::as<const VectorWrapper>(b);
384 res = std::move(gko_dense::create_with_config_of(b_cast).release());
385 }
386 else
387 {
388 res = std::move(gko::clone(b).release());
389 }
390 }
391 }
392
393private:
394 // Customize the logging hook which is called every time an iteration is
395 // completed.
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
401 {
402 // If the solver shares the current solution vector and we want to
403 // compute the residual from that
404 if (solution && compute_real_residual)
405 {
406 // Store the matrix's executor
407 auto exec = matrix->get_executor();
408 // Compute the real residual vector by calling apply on the system
409 // First, compute res = A * x
410#if MFEM_GINKGO_VERSION < 10600
411 matrix->apply(gko::lend(solution), gko::lend(res));
412#else
413 matrix->apply(solution, res);
414#endif
415 // Now do res = res - b, depending on which vector/oper type
416 // Check if b is a Ginkgo vector or wrapped MFEM Vector
417 if (dynamic_cast<const VectorWrapper*>(b))
418 {
419 const VectorWrapper *b_cast = gko::as<const VectorWrapper>(b);
420 // Copy the MFEM Vector stored in b
421 VectorWrapper *res_cast = gko::as<VectorWrapper>(res);
422 res_cast->get_mfem_vec_ref() -= b_cast->get_mfem_vec_const_ref();
423 }
424 else
425 {
426 // Create a scalar containing the value -1.0
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));
430#else
431 res->add_scaled(neg_one, b);
432#endif
433 }
434
435 // Compute the norm of the residual vector and add it to the
436 // `residual_norms` vector
437#if MFEM_GINKGO_VERSION < 10600
438 residual_norms.push_back(compute_norm(gko::lend(res)));
439#else
440 residual_norms.push_back(compute_norm(res));
441#endif
442 }
443 else
444 {
445 // If the solver shares an implicit or recurrent residual norm, log its value
446 if (implicit_sq_residual_norm)
447 {
448 auto dense_norm = gko::as<gko_dense>(implicit_sq_residual_norm);
449 // Add the norm to the `residual_norms` vector
450 residual_norms.push_back(get_norm(dense_norm));
451 // Otherwise, use the recurrent residual vector
452 }
453 else if (residual_norm)
454 {
455 auto dense_norm = gko::as<gko_dense>(residual_norm);
456 // Add the norm to the `residual_norms` vector
457 residual_norms.push_back(get_norm(dense_norm));
458 // Otherwise, use the recurrent residual vector
459 }
460 else
461 {
462 auto dense_residual = gko::as<gko_dense>(residual);
463 // Compute the residual vector's norm
464#if MFEM_GINKGO_VERSION < 10600
465 auto norm = compute_norm(gko::lend(dense_residual));
466#else
467 auto norm = compute_norm(dense_residual);
468#endif
469 // Add the computed norm to the `residual_norms` vector
470 residual_norms.push_back(norm);
471 }
472 }
473 // Add the current iteration number to the `iterations` vector
474 iterations.push_back(iteration);
475 }
476
477 // Pointer to the system matrix
478 const gko::LinOp *matrix;
479 // Pointer to the right hand sides
480 const gko_dense *b;
481 // Pointer to the residual workspace vector
482 gko_dense *res;
483 // Vector which stores all the residual norms
484 mutable std::vector<ValueType> residual_norms{};
485 // Vector which stores all the iteration numbers
486 mutable std::vector<std::size_t> iterations{};
487 // Whether or not to compute the residual at every iteration,
488 // rather than using the recurrent norm
489 const bool compute_real_residual;
490};
491
492/**
493* This class wraps a Ginkgo Executor for use in MFEM.
494* Note that objects in the Ginkgo namespace intended to work
495* together, e.g. a Ginkgo solver and preconditioner, should use the same
496* GinkgoExecutor object. In general, most users will want to create
497* one GinkgoExecutor object for use with all Ginkgo-related objects.
498* The wrapper can be created to match MFEM's device configuration.
499*/
501{
502public:
503 // Types of Ginkgo Executors.
505 {
506 /// Reference CPU Executor.
508 /// OpenMP CPU Executor.
509 OMP = 1,
510 /// CUDA GPU Executor.
511 CUDA = 2,
512 /// HIP GPU Executor.
513 HIP = 3
514 };
515 /**
516 * Constructor.
517 * Takes an @p GinkgoExecType argument and creates an Executor.
518 * In Ginkgo, GPU Executors must have an associated host Executor.
519 * This routine will select a CPU Executor based on the OpenMP support
520 * for Ginkgo.
521 */
522 GinkgoExecutor(ExecType exec_type);
523
524 /**
525 * Constructor.
526 * Takes an @p GinkgoExecType argument and creates an Executor.
527 * In Ginkgo, GPU Executors must have an associated host Executor.
528 * This routine allows for explicite setting of the CPU Executor
529 * for GPU backends.
530 */
531 GinkgoExecutor(ExecType exec_type, ExecType host_exec_type);
532
533 /**
534 * Constructor.
535 * Takes an MFEM @p Device object and creates an Executor
536 * that "matches" (e.g., if MFEM is using the CPU, Ginkgo
537 * will choose the Reference or OmpExecutor based on MFEM's
538 * configuration and Ginkgo's capabilities; if MFEM is using
539 * CUDA, Ginkgo will choose the CudaExecutor with a default
540 * CPU Executor based on Ginkgo's OpenMP support).
541 */
542 GinkgoExecutor(Device &mfem_device);
543
544 /**
545 * Constructor.
546 * Takes an MFEM @p Device object and creates an Executor
547 * that "matches", but allows the user to specify the host
548 * Executor for GPU backends.
549 */
550 GinkgoExecutor(Device &mfem_device, ExecType host_exec_type);
551
552 /**
553 * Destructor.
554 */
555 virtual ~GinkgoExecutor() = default;
556
557 std::shared_ptr<gko::Executor> GetExecutor() const
558 {
559 return this->executor;
560 };
561
562private:
563 std::shared_ptr<gko::Executor> executor;
564
565};
566
567/**
568* This class forms the base class for all of Ginkgo's preconditioners. The
569* various derived classes only take the additional data that is specific to them.
570* The entire collection of preconditioners that Ginkgo implements is available
571* at the Ginkgo documentation and manual pages,
572* https://ginkgo-project.github.io/ginkgo/doc/develop.
573*
574* @ingroup Ginkgo
575*/
577{
578public:
579 /**
580 * Constructor.
581 *
582 * The @p exec defines the paradigm where the solution is computed.
583 * Ginkgo currently supports four different executor types:
584 *
585 * + OmpExecutor specifies that the data should be stored and the
586 * associated operations executed on an OpenMP-supporting device (e.g.
587 * host CPU);
588 * + CudaExecutor specifies that the data should be stored and the
589 * operations executed on the NVIDIA GPU accelerator;
590 * + HipExecutor specifies that the data should be stored and the
591 * operations executed on the GPU accelerator using HIP;
592 * + ReferenceExecutor executes a non-optimized reference implementation,
593 * which can be used to debug the library.
594 */
596
597 /**
598 * Destructor.
599 */
600 virtual ~GinkgoPreconditioner() = default;
601
602 /**
603 * Generate the preconditioner for the given matrix @p op,
604 * which must be of MFEM SparseMatrix type.
605 * Calling this function is only required when creating a
606 * preconditioner for use with another MFEM solver; to use with
607 * a Ginkgo solver, get the LinOpFactory pointer through @p GetFactory()
608 * and pass to the Ginkgo solver constructor.
609 */
610 virtual void SetOperator(const Operator &op);
611
612 /**
613 * Apply the preconditioner to input vector @p x, with out @p y.
614 */
615 virtual void Mult(const Vector &x, Vector &y) const;
616
617 /**
618 * Return a pointer to the LinOpFactory that will generate the preconditioner
619 * with the parameters set through the specific constructor.
620 */
621 const std::shared_ptr<gko::LinOpFactory> GetFactory() const
622 {
623 return this->precond_gen;
624 };
625
626 /**
627 * Return a pointer to the generated preconditioner for a specific matrix
628 * (that has previously been set with @p SetOperator).
629 */
630 const std::shared_ptr<gko::LinOp> GetGeneratedPreconditioner() const
631 {
632 return this->generated_precond;
633 };
634
635 /**
636 * Return whether this GinkgoPreconditioner object has an explicitly-
637 * generated preconditioner, built for a specific matrix.
638 */
640 {
641 return this->has_generated_precond;
642 };
643
644protected:
645 /**
646 * The Ginkgo generated solver factory object.
647 */
648 std::shared_ptr<gko::LinOpFactory> precond_gen;
649
650 /**
651 * Generated Ginkgo preconditioner for a specific matrix, created through
652 * @p SetOperator(), or a wrapped MFEM preconditioner.
653 * Must exist to use @p Mult().
654 */
655 std::shared_ptr<gko::LinOp> generated_precond;
656
657 /**
658 * The execution paradigm in Ginkgo. The choices are between
659 * `gko::OmpExecutor`, `gko::CudaExecutor` and `gko::ReferenceExecutor`
660 * and more details can be found in Ginkgo's documentation.
661 */
662 std::shared_ptr<gko::Executor> executor;
663
664 /**
665 * Whether or not we have generated a specific preconditioner for
666 * a matrix.
667 */
669
670};
671
672/**
673* This class forms the base class for all of Ginkgo's iterative solvers.
674* It is not intended to be used directly by MFEM applications. The various
675* derived classes only take the additional data that is specific to them
676* and solve the given linear system. The entire collection of solvers that
677* Ginkgo implements is available at the Ginkgo documentation and manual pages,
678* https://ginkgo-project.github.io/ginkgo/doc/develop.
679*
680* @ingroup Ginkgo
681*/
683{
684public:
685 /**
686 * Return a pointer to the LinOpFactory that will generate the solver
687 * with the parameters set through the specific constructor.
688 */
689 const std::shared_ptr<gko::LinOpFactory> GetFactory() const
690 {
691 return this->solver_gen;
692 };
693
694 void SetPrintLevel(int print_lvl) { print_level = print_lvl; }
695
696 int GetNumIterations() const { return final_iter; }
697 int GetConverged() const { return converged; }
698 real_t GetFinalNorm() const { return final_norm; }
699
700 /**
701 * If the Operator is a SparseMatrix, set up a Ginkgo Csr matrix
702 * to use its data directly. If the Operator is not a matrix,
703 * create an OperatorWrapper for it and store.
704 */
705 virtual void SetOperator(const Operator &op);
706
707 /**
708 * Solve the linear system <tt>Ax=y</tt>. Dependent on the information
709 * provided by derived classes one of Ginkgo's linear solvers is chosen.
710 */
711 virtual void Mult(const Vector &x, Vector &y) const;
712
713 /**
714 * Return whether this GinkgoIterativeSolver object will use
715 * VectorWrapper types for input and output vectors.
716 * Note that Mult() will automatically create these wrappers if needed.
717 */
719 {
720 return this->needs_wrapped_vecs;
721 };
722
723 /**
724 * Destructor.
725 */
726 virtual ~GinkgoIterativeSolver() = default;
727
728protected:
729 /**
730 * Constructor.
731 *
732 * The @p exec defines the paradigm where the solution is computed.
733 * @p use_implicit_res_norm is for internal use by the derived classes
734 * for specific Ginkgo solvers; it indicates whether the solver makes
735 * an implicit residual norm estimate available for convergence checking.
736 * Each derived class automatically sets the correct value when calling this
737 * base class constructor.
738 *
739 */
742
749 mutable int final_iter;
750 mutable int converged;
751
752 /**
753 * The Ginkgo solver factory object, to generate specific solvers.
754 */
755 std::shared_ptr<gko::LinOpFactory> solver_gen;
756
757 /**
758 * The Ginkgo solver object, generated for a specific operator.
759 */
760 std::shared_ptr<gko::LinOp> solver;
761
762 /**
763 * The residual criterion object that controls the reduction of the residual
764 * relative to the initial residual.
765 */
766 std::shared_ptr<gko::stop::ResidualNorm<>::Factory>
768
769 /**
770 * The residual criterion object that controls the reduction of the residual
771 * based on an absolute tolerance.
772 */
773 std::shared_ptr<gko::stop::ResidualNorm<>::Factory>
775
776 /**
777 * The implicit residual criterion object that controls the reduction of the residual
778 * relative to the initial residual, based on an implicit residual norm value.
779 */
780 std::shared_ptr<gko::stop::ImplicitResidualNorm<>::Factory>
782
783 /**
784 * The implicit residual criterion object that controls the reduction of the residual
785 * based on an absolute tolerance, based on an implicit residual norm value.
786 */
787 std::shared_ptr<gko::stop::ImplicitResidualNorm<>::Factory>
789
790 /**
791 * The Ginkgo convergence logger used to check for convergence and other
792 * solver data if needed.
793 */
794 mutable std::shared_ptr<gko::log::Convergence<>> convergence_logger;
795
796 /**
797 * The residual logger object used to check for convergence and other solver
798 * data if needed.
799 */
800 mutable std::shared_ptr<ResidualLogger<>> residual_logger;
801
802 /**
803 * The Ginkgo combined factory object is used to create a combined stopping
804 * criterion to be passed to the solver.
805 */
806 std::shared_ptr<gko::stop::Combined::Factory> combined_factory;
807
808 /**
809 * The execution paradigm in Ginkgo. The choices are between
810 * `gko::OmpExecutor`, `gko::CudaExecutor` and `gko::ReferenceExecutor`
811 * and more details can be found in Ginkgo's documentation.
812 */
813 std::shared_ptr<gko::Executor> executor;
814
815 /**
816 * Whether or not we need to use VectorWrapper types with this solver.
817 */
819
820 /**
821 * Whether or not we need to use VectorWrapper types with the preconditioner
822 * or an inner solver. This value is set upon creation of the
823 * GinkgoIterativeSolver object and should never change.
824 */
826
827 /** Rebuild the Ginkgo stopping criterion factory with the latest values
828 * of rel_tol, abs_tol, and max_iter.
829 */
830 void update_stop_factory();
831
832private:
833 /**
834 * Initialize the Ginkgo logger object with event masks. Refer to the logging
835 * event masks in Ginkgo's .../include/ginkgo/core/log/logger.hpp.
836 */
837 void
838 initialize_ginkgo_log(gko::matrix::Dense<real_t>* b) const;
839
840 /**
841 * Pointer to either a Ginkgo CSR matrix or an OperatorWrapper wrapping
842 * an MFEM Operator (for matrix-free evaluation).
843 */
844 std::shared_ptr<gko::LinOp> system_oper;
845
846};
847
848/**
849 * This class adds helper functions for updating Ginkgo factories
850 * and solvers, when the full class type is needed. The derived classes
851 * should inherit from this class, rather than from GinkgoIterativeSolver
852 * directly.
853 */
854template<typename SolverType>
856{
857public:
860
861 void SetRelTol(real_t rtol)
862 {
863 rel_tol = rtol;
864 this->update_stop_factory();
865 auto current_params = gko::as<typename SolverType::Factory>
866 (solver_gen)->get_parameters();
867 this->solver_gen = current_params.with_criteria(this->combined_factory)
868 .on(this->executor);
869 if (solver)
870 {
871 gko::as<SolverType>(solver)->set_stop_criterion_factory(combined_factory);
872 }
873 }
874
875 void SetAbsTol(real_t atol)
876 {
877 abs_tol = atol;
878 this->update_stop_factory();
879 auto current_params = gko::as<typename SolverType::Factory>
880 (solver_gen)->get_parameters();
881 this->solver_gen = current_params.with_criteria(this->combined_factory)
882 .on(this->executor);
883 if (solver)
884 {
885 gko::as<SolverType>(solver)->set_stop_criterion_factory(combined_factory);
886 }
887 }
888
889 void SetMaxIter(int max_it)
890 {
891 max_iter = max_it;
892 this->update_stop_factory();
893 auto current_params = gko::as<typename SolverType::Factory>
894 (solver_gen)->get_parameters();
895 this->solver_gen = current_params.with_criteria(this->combined_factory)
896 .on(this->executor);
897 if (solver)
898 {
899 gko::as<SolverType>(solver)->set_stop_criterion_factory(combined_factory);
900 }
901 }
902};
903
904
905/**
906 * An implementation of the solver interface using the Ginkgo CG solver.
907 *
908 * @ingroup Ginkgo
909 */
910class CGSolver : public EnableGinkgoSolver<gko::solver::Cg<real_t>>
911{
912public:
913 /**
914 * Constructor.
915 *
916 * @param[in] exec The execution paradigm for the solver.
917 */
919
920 /**
921 * Constructor.
922 *
923 * @param[in] exec The execution paradigm for the solver.
924 * @param[in] preconditioner The preconditioner for the solver.
925 */
927 const GinkgoPreconditioner &preconditioner);
928};
929
930
931/**
932 * An implementation of the solver interface using the Ginkgo BiCGStab solver.
933 *
934 * @ingroup Ginkgo
935 */
936class BICGSTABSolver : public EnableGinkgoSolver<gko::solver::Bicgstab<real_t>>
937{
938public:
939 /**
940 * Constructor.
941 *
942 * @param[in] exec The execution paradigm for the solver.
943 */
945
946 /**
947 * Constructor.
948 *
949 * @param[in] exec The execution paradigm for the solver.
950 * @param[in] preconditioner The preconditioner for the solver.
951 */
953 const GinkgoPreconditioner &preconditioner);
954
955};
956
957/**
958 * An implementation of the solver interface using the Ginkgo CGS solver.
959 *
960 * CGS or the conjugate gradient square method is an iterative type Krylov
961 * subspace method which is suitable for general systems.
962 *
963 * @ingroup Ginkgo
964 */
965class CGSSolver : public EnableGinkgoSolver<gko::solver::Cgs<real_t>>
966{
967public:
968 /**
969 * Constructor.
970 *
971 * @param[in] exec The execution paradigm for the solver.
972 */
974
975 /**
976 * Constructor.
977 *
978 * @param[in] exec The execution paradigm for the solver.
979 * @param[in] preconditioner The preconditioner for the solver.
980 */
982 const GinkgoPreconditioner &preconditioner);
983};
984
985/**
986 * An implementation of the solver interface using the Ginkgo FCG solver.
987 *
988 * FCG or the flexible conjugate gradient method is an iterative type Krylov
989 * subspace method which is suitable for symmetric positive definite methods.
990 *
991 * Though this method performs very well for symmetric positive definite
992 * matrices, it is in general not suitable for general matrices.
993 *
994 * In contrast to the standard CG based on the Polack-Ribiere formula, the
995 * flexible CG uses the Fletcher-Reeves formula for creating the orthonormal
996 * vectors spanning the Krylov subspace. This increases the computational cost
997 * of every Krylov solver iteration but allows for non-constant preconditioners.
998 *
999 * @ingroup Ginkgo
1000 */
1001class FCGSolver : public EnableGinkgoSolver<gko::solver::Fcg<real_t>>
1002{
1003public:
1004 /**
1005 * Constructor.
1006 *
1007 * @param[in] exec The execution paradigm for the solver.
1008 */
1010
1011 /**
1012 * Constructor.
1013 *
1014 * @param[in] exec The execution paradigm for the solver.
1015 * @param[in] preconditioner The preconditioner for the solver.
1016 */
1018 const GinkgoPreconditioner &preconditioner);
1019};
1020
1021/**
1022 * An implementation of the solver interface using the Ginkgo GMRES solver.
1023 *
1024 * @ingroup Ginkgo
1025 */
1026class GMRESSolver : public EnableGinkgoSolver<gko::solver::Gmres<real_t>>
1027{
1028public:
1029 /**
1030 * Constructor.
1031 *
1032 * @param[in] exec The execution paradigm for the solver.
1033 * @param[in] dim The Krylov dimension of the solver. Value of 0 will
1034 * let Ginkgo use its own internal default value.
1035 */
1036 GMRESSolver(GinkgoExecutor &exec, int dim = 0);
1037
1038 /**
1039 * Constructor.
1040 *
1041 * @param[in] exec The execution paradigm for the solver.
1042 * @param[in] preconditioner The preconditioner for the solver.
1043 * @param[in] dim The Krylov dimension of the solver. Value of 0 will
1044 * let Ginkgo use its own internal default value.
1045 */
1047 const GinkgoPreconditioner &preconditioner,
1048 int dim = 0);
1049
1050 /**
1051 * Change the Krylov dimension of the solver.
1052 */
1053 void SetKDim(int dim);
1054
1055protected:
1056 int m; // Dimension of Krylov subspace
1057};
1058
1059using gko::solver::cb_gmres::storage_precision;
1060/**
1061 * An implementation of the solver interface using the Ginkgo
1062 * Compressed Basis GMRES solver. With CB-GMRES, the Krylov basis
1063 * is "compressed" by storing in a lower precision. Currently, computations
1064 * are always performed in the MFEM-defined `real_t` precision, when using this
1065 * MFEM integration.
1066 * The Ginkgo storage precision options are accessed
1067 * through Ginkgo::storage_precision::*. The default choice
1068 * is Ginkgo::storage_precision::reduce1, i.e., store in float
1069 * instead of double or half instead of float.
1070 *
1071 * @ingroup Ginkgo
1072 */
1073class CBGMRESSolver : public EnableGinkgoSolver<gko::solver::CbGmres<real_t>>
1074{
1075public:
1076 /**
1077 * Constructor.
1078 *
1079 * @param[in] exec The execution paradigm for the solver.
1080 * @param[in] dim The Krylov dimension of the solver. Value of 0 will
1081 * let Ginkgo use its own internal default value.
1082 * @param[in] prec The storage precision used in the CB-GMRES. Options
1083 * are: keep (keep `real_t` precision), reduce1 (double -> float
1084 * or float -> half), reduce2 (double -> half or float -> half),
1085 * integer (`real_t` -> int64), ireduce1 (double -> int32 or
1086 * float -> int16), ireduce2 (double -> int16 or float -> int16).
1087 * See Ginkgo documentation for more about CB-GMRES.
1088 */
1089 CBGMRESSolver(GinkgoExecutor &exec, int dim = 0,
1090 storage_precision prec = storage_precision::reduce1);
1091
1092 /**
1093 * Constructor.
1094 *
1095 * @param[in] exec The execution paradigm for the solver.
1096 * @param[in] preconditioner The preconditioner for the solver.
1097 * @param[in] dim The Krylov dimension of the solver. Value of 0 will
1098 * let Ginkgo use its own internal default value.
1099 * @param[in] prec The storage precision used in the CB-GMRES. Options
1100 * are: keep (keep `real_t` precision), reduce1 (double -> float
1101 * or float -> half), reduce2 (double -> half or float -> half),
1102 * integer (`real_t` -> int64), ireduce1 (double -> int32 or
1103 * float -> int16), ireduce2 (double -> int16 or float -> int16).
1104 * See Ginkgo documentation for more about CB-GMRES.
1105 */
1107 const GinkgoPreconditioner &preconditioner,
1108 int dim = 0,
1109 storage_precision prec = storage_precision::reduce1);
1110
1111 /**
1112 * Change the Krylov dimension of the solver.
1113 */
1114 void SetKDim(int dim);
1115
1116protected:
1117 int m; // Dimension of Krylov subspace
1118};
1119
1120/**
1121 * An implementation of the solver interface using the Ginkgo IR solver.
1122 *
1123 * Iterative refinement (IR) is an iterative method that uses another coarse
1124 * method to approximate the error of the current solution via the current
1125 * residual.
1126 *
1127 * @ingroup Ginkgo
1128 */
1129class IRSolver : public EnableGinkgoSolver<gko::solver::Ir<real_t>>
1130{
1131public:
1132 /**
1133 * Constructor.
1134 *
1135 * @param[in] exec The execution paradigm for the solver.
1136 */
1137 IRSolver(GinkgoExecutor &exec);
1138
1139 /**
1140 * Constructor.
1141 *
1142 * @param[in] exec The execution paradigm for the solver.
1143 * @param[in] inner_solver The inner solver for the main solver.
1144 */
1146 const GinkgoIterativeSolver &inner_solver);
1147
1148};
1149
1150/**
1151 * An implementation of the preconditioner interface using the Ginkgo Jacobi
1152 * preconditioner.
1153 *
1154 * @ingroup Ginkgo
1155 */
1157{
1158public:
1159 /**
1160 * Constructor.
1161 *
1162 * @param[in] exec The execution paradigm for the preconditioner.
1163 * @param[in] storage_opt The storage optimization parameter.
1164 * @param[in] accuracy The relative accuracy for the adaptive version.
1165 * @param[in] max_block_size Maximum block size.
1166 * See the Ginkgo documentation for more information on these parameters.
1167 */
1169 GinkgoExecutor &exec,
1170 const std::string &storage_opt = "none",
1171 const real_t accuracy = 1.e-1,
1172 const int max_block_size = 32
1173 );
1174};
1175
1176/**
1177 * An implementation of the preconditioner interface using the Ginkgo
1178 * Incomplete LU preconditioner (ILU(0)).
1179 *
1180 * @ingroup Ginkgo
1181 */
1183{
1184public:
1185 /**
1186 * Constructor.
1187 *
1188 * @param[in] exec The execution paradigm for the preconditioner.
1189 * @param[in] factorization_type The factorization type: "exact" or
1190 * "parilu".
1191 * @param[in] sweeps The number of sweeps to do in the ParIlu
1192 * factorization algorithm. A value of 0 tells Ginkgo to use its
1193 * internal default value. This parameter is ignored in the case
1194 * of an exact factorization.
1195 * @param[in] skip_sort Only set this to true if the input matrix
1196 * that will be used to generate this preconditioner is guaranteed
1197 * to be sorted by column.
1198 *
1199 * Note: The use of this preconditioner will sort any input matrix
1200 * given to it, potentially changing the order of the stored values.
1201 */
1203 GinkgoExecutor &exec,
1204 const std::string &factorization_type = "exact",
1205 const int sweeps = 0,
1206 const bool skip_sort = false
1207 );
1208};
1209
1210/**
1211 * An implementation of the preconditioner interface using the Ginkgo
1212 * Incomplete LU-Incomplete Sparse Approximate Inverse preconditioner.
1213 * The Ilu-ISAI preconditioner differs from the Ilu preconditioner in
1214 * that Incomplete Sparse Approximate Inverses (ISAIs) are formed
1215 * to approximate solving the triangular systems defined by L and U.
1216 * When the preconditioner is applied, these ISAI matrices are applied
1217 * through matrix-vector multiplication.
1218 *
1219 * @ingroup Ginkgo
1220 */
1222{
1223public:
1224 /**
1225 * Constructor.
1226 *
1227 * @param[in] exec The execution paradigm for the preconditioner.
1228 * @param[in] factorization_type The factorization type: "exact" or
1229 * "parilu".
1230 * @param[in] sweeps The number of sweeps to do in the ParIlu
1231 * factorization algorithm. A value of 0 tells Ginkgo to use its
1232 * internal default value. This parameter is ignored in the case
1233 * of an exact factorization.
1234 * @param[in] sparsity_power Parameter determining the sparsity pattern of
1235 * the ISAI approximations.
1236 * @param[in] skip_sort Only set this to true if the input matrix
1237 * that will be used to generate this preconditioner is guaranteed
1238 * to be sorted by column.
1239 * See the Ginkgo documentation for more information on these parameters.
1240 *
1241 * Note: The use of this preconditioner will sort any input matrix
1242 * given to it, potentially changing the order of the stored values.
1243 */
1245 GinkgoExecutor &exec,
1246 const std::string &factorization_type = "exact",
1247 const int sweeps = 0,
1248 const int sparsity_power = 1,
1249 const bool skip_sort = false
1250 );
1251};
1252
1253/**
1254 * An implementation of the preconditioner interface using the Ginkgo
1255 * Incomplete Cholesky preconditioner (IC(0)).
1256 *
1257 * @ingroup Ginkgo
1258 */
1260{
1261public:
1262 /**
1263 * Constructor.
1264 *
1265 * @param[in] exec The execution paradigm for the preconditioner.
1266 * @param[in] factorization_type The factorization type: "exact" or
1267 * "paric".
1268 * @param[in] sweeps The number of sweeps to do in the ParIc
1269 * factorization algorithm. A value of 0 tells Ginkgo to use its
1270 * internal default value. This parameter is ignored in the case
1271 * of an exact factorization.
1272 * @param[in] skip_sort Only set this to true if the input matrix
1273 * that will be used to generate this preconditioner is guaranteed
1274 * to be sorted by column.
1275 *
1276 * Note: The use of this preconditioner will sort any input matrix
1277 * given to it, potentially changing the order of the stored values.
1278 */
1280 GinkgoExecutor &exec,
1281 const std::string &factorization_type = "exact",
1282 const int sweeps = 0,
1283 const bool skip_sort = false
1284 );
1285};
1286
1287/**
1288 * An implementation of the preconditioner interface using the Ginkgo
1289 * Incomplete Cholesky-Incomplete Sparse Approximate Inverse preconditioner.
1290 * The Ic-ISAI preconditioner differs from the Ic preconditioner in
1291 * that Incomplete Sparse Approximate Inverses (ISAIs) are formed
1292 * to approximate solving the triangular systems defined by L and L^T.
1293 * When the preconditioner is applied, these ISAI matrices are applied
1294 * through matrix-vector multiplication.
1295 *
1296 * @ingroup Ginkgo
1297 */
1299{
1300public:
1301 /**
1302 * Constructor.
1303 *
1304 * @param[in] exec The execution paradigm for the preconditioner.
1305 * @param[in] factorization_type The factorization type: "exact" or
1306 * "paric".
1307 * @param[in] sweeps The number of sweeps to do in the ParIc
1308 * factorization algorithm. A value of 0 tells Ginkgo to use its
1309 * internal default value. This parameter is ignored in the case
1310 * of an exact factorization.
1311 * @param[in] sparsity_power Parameter determining the sparsity pattern of
1312 * the ISAI approximations.
1313 * @param[in] skip_sort Only set this to true if the input matrix
1314 * that will be used to generate this preconditioner is guaranteed
1315 * to be sorted by column.
1316 * See the Ginkgo documentation for more information on these parameters.
1317 *
1318 * Note: The use of this preconditioner will sort any input matrix
1319 * given to it, potentially changing the order of the stored values.
1320 */
1322 GinkgoExecutor &exec,
1323 const std::string &factorization_type = "exact",
1324 const int sweeps = 0,
1325 const int sparsity_power = 1,
1326 const bool skip_sort = false
1327 );
1328};
1329
1330/**
1331 * A wrapper that allows Ginkgo to use MFEM preconditioners.
1332 *
1333 * @ingroup Ginkgo
1334 */
1336{
1337public:
1338 /**
1339 * Constructor.
1340 *
1341 * @param[in] exec The execution paradigm for the preconditioner.
1342 * @param[in] mfem_precond The MFEM Preconditioner to wrap.
1343 */
1345 GinkgoExecutor &exec,
1346 const Solver &mfem_precond
1347 );
1348
1349 /**
1350 * SetOperator is not allowed for this type of preconditioner;
1351 * this function overrides the base class in order to give an
1352 * error if SetOperator() is called for this class.
1353 */
1354 virtual void SetOperator(const Operator &op)
1355 {
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.");
1359 };
1360};
1361} // namespace Ginkgo
1362
1363} // namespace mfem
1364
1365#endif // MFEM_USE_GINKGO
1366
1367#endif // MFEM_GINKGO
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
BICGSTABSolver(GinkgoExecutor &exec)
Definition ginkgo.cpp:743
CBGMRESSolver(GinkgoExecutor &exec, int dim=0, storage_precision prec=storage_precision::reduce1)
Definition ginkgo.cpp:952
CGSSolver(GinkgoExecutor &exec)
Definition ginkgo.cpp:782
CGSolver(GinkgoExecutor &exec)
Definition ginkgo.cpp:704
void SetRelTol(real_t rtol)
Definition ginkgo.hpp:861
void SetAbsTol(real_t atol)
Definition ginkgo.hpp:875
EnableGinkgoSolver(GinkgoExecutor &exec, bool use_implicit_res_norm)
Definition ginkgo.hpp:858
FCGSolver(GinkgoExecutor &exec)
Definition ginkgo.cpp:820
GMRESSolver(GinkgoExecutor &exec, int dim=0)
Definition ginkgo.cpp:858
void SetKDim(int dim)
Definition ginkgo.cpp:937
std::shared_ptr< gko::Executor > GetExecutor() const
Definition ginkgo.hpp:557
GinkgoExecutor(ExecType exec_type)
Definition ginkgo.cpp:33
virtual ~GinkgoExecutor()=default
@ CUDA
CUDA GPU Executor.
Definition ginkgo.hpp:511
@ HIP
HIP GPU Executor.
Definition ginkgo.hpp:513
@ REFERENCE
Reference CPU Executor.
Definition ginkgo.hpp:507
@ OMP
OpenMP CPU Executor.
Definition ginkgo.hpp:509
std::shared_ptr< gko::stop::ImplicitResidualNorm<>::Factory > imp_abs_criterion
Definition ginkgo.hpp:788
virtual void Mult(const Vector &x, Vector &y) const
Definition ginkgo.cpp:508
std::shared_ptr< ResidualLogger<> > residual_logger
Definition ginkgo.hpp:800
std::shared_ptr< gko::stop::ResidualNorm<>::Factory > abs_criterion
Definition ginkgo.hpp:774
std::shared_ptr< gko::stop::Combined::Factory > combined_factory
Definition ginkgo.hpp:806
const std::shared_ptr< gko::LinOpFactory > GetFactory() const
Definition ginkgo.hpp:689
std::shared_ptr< gko::log::Convergence<> > convergence_logger
Definition ginkgo.hpp:794
GinkgoIterativeSolver(GinkgoExecutor &exec, bool use_implicit_res_norm)
Definition ginkgo.cpp:342
void SetPrintLevel(int print_lvl)
Definition ginkgo.hpp:694
std::shared_ptr< gko::LinOpFactory > solver_gen
Definition ginkgo.hpp:755
std::shared_ptr< gko::stop::ResidualNorm<>::Factory > rel_criterion
Definition ginkgo.hpp:767
virtual ~GinkgoIterativeSolver()=default
std::shared_ptr< gko::stop::ImplicitResidualNorm<>::Factory > imp_rel_criterion
Definition ginkgo.hpp:781
std::shared_ptr< gko::Executor > executor
Definition ginkgo.hpp:813
virtual void SetOperator(const Operator &op)
Definition ginkgo.cpp:640
std::shared_ptr< gko::LinOp > solver
Definition ginkgo.hpp:760
std::shared_ptr< gko::Executor > executor
Definition ginkgo.hpp:662
const std::shared_ptr< gko::LinOpFactory > GetFactory() const
Definition ginkgo.hpp:621
virtual ~GinkgoPreconditioner()=default
std::shared_ptr< gko::LinOp > generated_precond
Definition ginkgo.hpp:655
virtual void SetOperator(const Operator &op)
Definition ginkgo.cpp:1123
virtual void Mult(const Vector &x, Vector &y) const
Definition ginkgo.cpp:1089
GinkgoPreconditioner(GinkgoExecutor &exec)
Definition ginkgo.cpp:1080
std::shared_ptr< gko::LinOpFactory > precond_gen
Definition ginkgo.hpp:648
const std::shared_ptr< gko::LinOp > GetGeneratedPreconditioner() const
Definition ginkgo.hpp:630
IRSolver(GinkgoExecutor &exec)
Definition ginkgo.cpp:1054
IcIsaiPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const int sparsity_power=1, const bool skip_sort=false)
Definition ginkgo.cpp:1353
IcPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const bool skip_sort=false)
Definition ginkgo.cpp:1309
IluIsaiPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const int sparsity_power=1, const bool skip_sort=false)
Definition ginkgo.cpp:1240
IluPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const bool skip_sort=false)
Definition ginkgo.cpp:1198
JacobiPreconditioner(GinkgoExecutor &exec, const std::string &storage_opt="none", const real_t accuracy=1.e-1, const int max_block_size=32)
Definition ginkgo.cpp:1167
MFEMPreconditioner(GinkgoExecutor &exec, const Solver &mfem_precond)
Definition ginkgo.cpp:1408
virtual void SetOperator(const Operator &op)
Definition ginkgo.hpp:1354
OperatorWrapper(std::shared_ptr< const gko::Executor > exec, gko::size_type size=0, const Operator *oper=NULL)
Definition ginkgo.hpp:229
void apply_impl(const gko::LinOp *b, gko::LinOp *x) const override
Definition ginkgo.cpp:423
static std::unique_ptr< VectorWrapper > create(std::shared_ptr< const gko::Executor > exec, gko::size_type size, Vector *mfem_vec, bool ownership=false)
Definition ginkgo.hpp:108
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
Definition ginkgo.hpp:148
VectorWrapper(std::shared_ptr< const gko::Executor > exec, gko::size_type size, Vector *mfem_vec, bool ownership=false)
Definition ginkgo.hpp:75
virtual std::unique_ptr< gko::matrix::Dense< real_t > > create_with_same_config() const override
Definition ginkgo.hpp:127
const Vector & get_mfem_vec_const_ref() const
Definition ginkgo.hpp:122
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
Definition ginkgo.hpp:178
void operator()(pointer ptr) const noexcept
Definition ginkgo.hpp:60
Abstract operator.
Definition operator.hpp:25
Base class for solvers.
Definition operator.hpp:683
Vector data type.
Definition vector.hpp:80
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:136
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:602
Vector beta
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
real_t get_norm(const gko::matrix::Dense< ValueType > *norm)
Definition ginkgo.hpp:250
real_t compute_norm(const gko::matrix::Dense< ValueType > *b)
Definition ginkgo.hpp:261
gko::Array< T > gko_array
Definition ginkgo.hpp:43
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
float real_t
Definition config.hpp:43
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
Definition ginkgo.hpp:327
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
Definition ginkgo.hpp:349
ResidualLogger(std::shared_ptr< const gko::Executor > exec, const gko::LinOp *matrix, const gko_dense *b, bool compute_real_residual=false)
Definition ginkgo.hpp:365
gko::matrix::Dense< ValueType > gko_dense
Definition ginkgo.hpp:324
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
Definition ginkgo.hpp:337