MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ginkgo.cpp
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#include "../config/config.hpp"
13
14#ifdef MFEM_USE_GINKGO
15
16#include "ginkgo.hpp"
17#include "sparsemat.hpp"
19#include "../general/error.hpp"
20#include <iostream>
21#include <iomanip>
22#include <algorithm>
23#include <cmath>
24#include <cstring>
25
26namespace mfem
27{
28
29namespace Ginkgo
30{
31
32// Create a GinkgoExecutor of type exec_type.
34{
35#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
36 gko::version_info gko_version = gko::version_info::get();
37 bool gko_with_omp_support = (strcmp(gko_version.omp_version.tag,
38 "not compiled") != 0);
39#endif
40 switch (exec_type)
41 {
43 {
44 executor = gko::ReferenceExecutor::create();
45 break;
46 }
48 {
49 executor = gko::OmpExecutor::create();
50 break;
51 }
53 {
54 if (gko::CudaExecutor::get_num_devices() > 0)
55 {
56#ifdef MFEM_USE_CUDA
57 int current_device = 0;
58 MFEM_GPU_CHECK(cudaGetDevice(&current_device));
59 if (gko_with_omp_support)
60 {
61 executor = gko::CudaExecutor::create(current_device,
62 gko::OmpExecutor::create());
63 }
64 else
65 {
66 executor = gko::CudaExecutor::create(current_device,
67 gko::ReferenceExecutor::create());
68 }
69#endif
70 }
71 else
72 {
73 MFEM_ABORT("gko::CudaExecutor::get_num_devices() did not report "
74 "any valid devices.");
75 }
76 break;
77 }
79 {
80 if (gko::HipExecutor::get_num_devices() > 0)
81 {
82#ifdef MFEM_USE_HIP
83 int current_device = 0;
84 MFEM_GPU_CHECK(hipGetDevice(&current_device));
85 if (gko_with_omp_support)
86 {
87 executor = gko::HipExecutor::create(current_device,
88 gko::OmpExecutor::create());
89 }
90 else
91 {
92 executor = gko::HipExecutor::create(current_device,
93 gko::ReferenceExecutor::create());
94 }
95#endif
96 }
97 else
98 {
99 MFEM_ABORT("gko::HipExecutor::get_num_devices() did not report "
100 "any valid devices.");
101 }
102 break;
103 }
104 default:
105 MFEM_ABORT("Invalid ExecType specified");
106 }
107}
108
109// Create a GinkgoExecutor of type exec_type, with host_exec_type for the
110// related CPU Executor (only applicable to GPU backends).
112{
113 switch (exec_type)
114 {
116 {
117 MFEM_WARNING("Parameter host_exec_type ignored for CPU GinkgoExecutor.");
118 executor = gko::ReferenceExecutor::create();
119 break;
120 }
122 {
123 MFEM_WARNING("Parameter host_exec_type ignored for CPU GinkgoExecutor.");
124 executor = gko::OmpExecutor::create();
125 break;
126 }
128 {
129 if (gko::CudaExecutor::get_num_devices() > 0)
130 {
131#ifdef MFEM_USE_CUDA
132 int current_device = 0;
133 MFEM_GPU_CHECK(cudaGetDevice(&current_device));
134 if (host_exec_type == GinkgoExecutor::OMP)
135 {
136 executor = gko::CudaExecutor::create(current_device,
137 gko::OmpExecutor::create());
138 }
139 else
140 {
141 executor = gko::CudaExecutor::create(current_device,
142 gko::ReferenceExecutor::create());
143 }
144#endif
145 }
146 else
147 {
148 MFEM_ABORT("gko::CudaExecutor::get_num_devices() did not report "
149 "any valid devices.");
150 }
151 break;
152 }
154 {
155 if (gko::HipExecutor::get_num_devices() > 0)
156 {
157#ifdef MFEM_USE_HIP
158 int current_device = 0;
159 MFEM_GPU_CHECK(hipGetDevice(&current_device));
160 if (host_exec_type == GinkgoExecutor::OMP)
161 {
162 executor = gko::HipExecutor::create(current_device,
163 gko::OmpExecutor::create());
164 }
165 else
166 {
167 executor = gko::HipExecutor::create(current_device,
168 gko::ReferenceExecutor::create());
169 }
170#endif
171 }
172 else
173 {
174 MFEM_ABORT("gko::HipExecutor::get_num_devices() did not report "
175 "any valid devices.");
176 }
177 break;
178 }
179 default:
180 MFEM_ABORT("Invalid ExecType specified");
181 }
182}
183
184// Create a GinkgoExecutor to match MFEM's device configuration.
186{
187 gko::version_info gko_version = gko::version_info::get();
188 bool gko_with_omp_support = (strcmp(gko_version.omp_version.tag,
189 "not compiled") != 0);
190 if (mfem_device.Allows(Backend::CUDA_MASK))
191 {
192 if (gko::CudaExecutor::get_num_devices() > 0)
193 {
194#ifdef MFEM_USE_CUDA
195 int current_device = 0;
196 MFEM_GPU_CHECK(cudaGetDevice(&current_device));
197 if (gko_with_omp_support)
198 {
199 executor = gko::CudaExecutor::create(current_device,
200 gko::OmpExecutor::create());
201 }
202 else
203 {
204 executor = gko::CudaExecutor::create(current_device,
205 gko::ReferenceExecutor::create());
206 }
207#endif
208 }
209 else
210 {
211 MFEM_ABORT("gko::CudaExecutor::get_num_devices() did not report "
212 "any valid devices.");
213 }
214 }
215 else if (mfem_device.Allows(Backend::HIP_MASK))
216 {
217 if (gko::HipExecutor::get_num_devices() > 0)
218 {
219#ifdef MFEM_USE_HIP
220 int current_device = 0;
221 MFEM_GPU_CHECK(hipGetDevice(&current_device));
222 if (gko_with_omp_support)
223 {
224 executor = gko::HipExecutor::create(current_device,
225 gko::OmpExecutor::create());
226 }
227 else
228 {
229 executor = gko::HipExecutor::create(current_device,
230 gko::ReferenceExecutor::create());
231 }
232#endif
233 }
234 else
235 {
236 MFEM_ABORT("gko::HipExecutor::get_num_devices() did not report "
237 "any valid devices.");
238 }
239 }
240 else
241 {
242 if (mfem_device.Allows(Backend::OMP_MASK))
243 {
244 // Also use OpenMP for Ginkgo, if Ginkgo supports it
245 if (gko_with_omp_support)
246 {
247 executor = gko::OmpExecutor::create();
248 }
249 else
250 {
251 executor = gko::ReferenceExecutor::create();
252 }
253 }
254 else
255 {
256 executor = gko::ReferenceExecutor::create();
257 }
258 }
259}
260
261// Create a GinkgoExecutor to match MFEM's device configuration, with
262// a specific host_exec_type for the associated CPU Executor (only
263// applicable to GPU backends).
264GinkgoExecutor::GinkgoExecutor(Device &mfem_device, ExecType host_exec_type)
265{
266
267 if (mfem_device.Allows(Backend::CUDA_MASK))
268 {
269 if (gko::CudaExecutor::get_num_devices() > 0)
270 {
271#ifdef MFEM_USE_CUDA
272 int current_device = 0;
273 MFEM_GPU_CHECK(cudaGetDevice(&current_device));
274 if (host_exec_type == GinkgoExecutor::OMP)
275 {
276 executor = gko::CudaExecutor::create(current_device,
277 gko::OmpExecutor::create());
278 }
279 else
280 {
281 executor = gko::CudaExecutor::create(current_device,
282 gko::ReferenceExecutor::create());
283 }
284#endif
285 }
286 else
287 {
288 MFEM_ABORT("gko::CudaExecutor::get_num_devices() did not report "
289 "any valid devices.");
290 }
291 }
292 else if (mfem_device.Allows(Backend::HIP_MASK))
293 {
294 if (gko::HipExecutor::get_num_devices() > 0)
295 {
296#ifdef MFEM_USE_HIP
297 int current_device = 0;
298 MFEM_GPU_CHECK(hipGetDevice(&current_device));
299 if (host_exec_type == GinkgoExecutor::OMP)
300 {
301 executor = gko::HipExecutor::create(current_device,
302 gko::OmpExecutor::create());
303 }
304 else
305 {
306 executor = gko::HipExecutor::create(current_device,
307 gko::ReferenceExecutor::create());
308 }
309#endif
310 }
311 else
312 {
313 MFEM_ABORT("gko::HipExecutor::get_num_devices() did not report "
314 "any valid devices.");
315 }
316 }
317 else
318 {
319 MFEM_WARNING("Parameter host_exec_type ignored for CPU GinkgoExecutor.");
320 if (mfem_device.Allows(Backend::OMP_MASK))
321 {
322 // Also use OpenMP for Ginkgo, if Ginkgo supports it
323 gko::version_info gko_version = gko::version_info::get();
324 bool gko_with_omp_support = (strcmp(gko_version.omp_version.tag,
325 "not compiled") != 0);
326 if (gko_with_omp_support)
327 {
328 executor = gko::OmpExecutor::create();
329 }
330 else
331 {
332 executor = gko::ReferenceExecutor::create();
333 }
334 }
335 else
336 {
337 executor = gko::ReferenceExecutor::create();
338 }
339 }
340}
341
343 bool use_implicit_res_norm)
344 : Solver(),
345 use_implicit_res_norm(use_implicit_res_norm)
346{
347 executor = exec.GetExecutor();
348 print_level = -1;
349
350 // Build default stopping criterion factory
351 max_iter = 10;
352 rel_tol = 0.0;
353 abs_tol = 0.0;
354 this->update_stop_factory();
355
356 needs_wrapped_vecs = false;
358}
359
361{
362 using ResidualCriterionFactory = gko::stop::ResidualNorm<>;
363 using ImplicitResidualCriterionFactory = gko::stop::ImplicitResidualNorm<>;
364
366 {
367 imp_rel_criterion = ImplicitResidualCriterionFactory::build()
368 .with_reduction_factor(sqrt(rel_tol))
369 .with_baseline(gko::stop::mode::initial_resnorm)
370 .on(executor);
371 imp_abs_criterion = ImplicitResidualCriterionFactory::build()
372 .with_reduction_factor(sqrt(abs_tol))
373 .with_baseline(gko::stop::mode::absolute)
374 .on(executor);
376 gko::stop::Combined::build()
377 .with_criteria(imp_rel_criterion,
379 gko::stop::Iteration::build()
380 .with_max_iters(static_cast<unsigned long>(max_iter))
381 .on(executor))
382 .on(executor);
383 }
384 else
385 {
386 rel_criterion = ResidualCriterionFactory::build()
387 .with_reduction_factor(rel_tol)
388 .with_baseline(gko::stop::mode::initial_resnorm)
389 .on(executor);
390 abs_criterion = ResidualCriterionFactory::build()
391 .with_reduction_factor(abs_tol)
392 .with_baseline(gko::stop::mode::absolute)
393 .on(executor);
395 gko::stop::Combined::build()
396 .with_criteria(rel_criterion,
398 gko::stop::Iteration::build()
399 .with_max_iters(static_cast<unsigned long>(max_iter))
400 .on(executor))
401 .on(executor);
402 }
403}
404
405void
406GinkgoIterativeSolver::initialize_ginkgo_log(gko::matrix::Dense<real_t>* b)
407const
408{
409 // Add the logger object. See the different masks available in Ginkgo's
410 // documentation
411#if MFEM_GINKGO_VERSION < 10500
412 convergence_logger = gko::log::Convergence<>::create(
413 executor, gko::log::Logger::criterion_check_completed_mask);
414#else
415 convergence_logger = gko::log::Convergence<>::create(
416 gko::log::Logger::criterion_check_completed_mask);
417#endif
418 residual_logger = std::make_shared<ResidualLogger<>>(executor,
419 system_oper.get(),b);
420
421}
422
423void OperatorWrapper::apply_impl(const gko::LinOp *b, gko::LinOp *x) const
424{
425
426 // Cast to VectorWrapper; only accept this type for this impl
427 const VectorWrapper *mfem_b = gko::as<const VectorWrapper>(b);
428 VectorWrapper *mfem_x = gko::as<VectorWrapper>(x);
429
430 this->wrapped_oper->Mult(mfem_b->get_mfem_vec_const_ref(),
431 mfem_x->get_mfem_vec_ref());
432}
433void OperatorWrapper::apply_impl(const gko::LinOp *alpha,
434 const gko::LinOp *b,
435 const gko::LinOp *beta,
436 gko::LinOp *x) const
437{
438 // x = alpha * op (b) + beta * x
439 // Cast to VectorWrapper; only accept this type for this impl
440 const VectorWrapper *mfem_b = gko::as<const VectorWrapper>(b);
441 VectorWrapper *mfem_x = gko::as<VectorWrapper>(x);
442
443 // Check that alpha and beta are Dense<real_t> of size (1,1):
444 if (alpha->get_size()[0] > 1 || alpha->get_size()[1] > 1)
445 {
446 throw gko::BadDimension(
447 __FILE__, __LINE__, __func__, "alpha", alpha->get_size()[0],
448 alpha->get_size()[1],
449 "Expected an object of size [1 x 1] for scaling "
450 " in this operator's apply_impl");
451 }
452 if (beta->get_size()[0] > 1 || beta->get_size()[1] > 1)
453 {
454 throw gko::BadDimension(
455 __FILE__, __LINE__, __func__, "beta", beta->get_size()[0],
456 beta->get_size()[1],
457 "Expected an object of size [1 x 1] for scaling "
458 " in this operator's apply_impl");
459 }
460 real_t alpha_f;
461 real_t beta_f;
462
463 if (alpha->get_executor() == alpha->get_executor()->get_master())
464 {
465 // Access value directly
466 alpha_f = gko::as<gko::matrix::Dense<real_t>>(alpha)->at(0, 0);
467 }
468 else
469 {
470 // Copy from device to host
471 this->get_executor()->get_master().get()->copy_from(
472 this->get_executor().get(),
473 1, gko::as<gko::matrix::Dense<real_t>>(alpha)->get_const_values(),
474 &alpha_f);
475 }
476 if (beta->get_executor() == beta->get_executor()->get_master())
477 {
478 // Access value directly
479 beta_f = gko::as<gko::matrix::Dense<real_t>>(beta)->at(0, 0);
480 }
481 else
482 {
483 // Copy from device to host
484 this->get_executor()->get_master().get()->copy_from(
485 this->get_executor().get(),
486 1, gko::as<gko::matrix::Dense<real_t>>(beta)->get_const_values(),
487 &beta_f);
488 }
489 // Scale x by beta
490 mfem_x->get_mfem_vec_ref() *= beta_f;
491 // Multiply operator with b and store in tmp
492 mfem::Vector mfem_tmp =
493 mfem::Vector(mfem_x->get_size()[0],
495 // Set UseDevice flag to match mfem_x (not automatically done through
496 // MemoryType)
497 mfem_tmp.UseDevice(mfem_x->get_mfem_vec_ref().UseDevice());
498
499 // Apply the operator
500 this->wrapped_oper->Mult(mfem_b->get_mfem_vec_const_ref(), mfem_tmp);
501 // Scale tmp by alpha and add
502 mfem_x->get_mfem_vec_ref().Add(alpha_f, mfem_tmp);
503
504 mfem_tmp.Destroy();
505}
506
507void
509{
510
511 MFEM_VERIFY(system_oper, "System matrix or operator not initialized");
512 MFEM_VERIFY(executor, "executor is not initialized");
513 MFEM_VERIFY(y.Size() == x.Size(),
514 "Mismatching sizes for rhs and solution");
515
516 using vec = gko::matrix::Dense<real_t>;
517 if (!iterative_mode)
518 {
519 y = 0.0;
520 }
521
522 // Create x and y vectors in Ginkgo's format. Wrap MFEM's data directly,
523 // on CPU or GPU.
524 bool on_device = false;
525 if (executor != executor->get_master())
526 {
527 on_device = true;
528 }
529 std::unique_ptr<vec> gko_x;
530 std::unique_ptr<vec> gko_y;
531
532 // If we do not have an OperatorWrapper for the system operator or
533 // preconditioner, or have an inner solver using VectorWrappers (as
534 // for IR), then directly create Ginkgo vectors from MFEM's data.
536 {
537 gko_x = vec::create(executor, gko::dim<2>(x.Size(), 1),
539 x.Size(), const_cast<real_t *>(
540 x.Read(on_device))), 1);
541 gko_y = vec::create(executor, gko::dim<2>(y.Size(), 1),
543 y.Size(),
544 y.ReadWrite(on_device)), 1);
545 }
546 else // We have at least one wrapped MFEM operator; need wrapped vectors
547 {
548 gko_x = std::unique_ptr<vec>(
549 new VectorWrapper(executor, x.Size(),
550 const_cast<Vector *>(&x), false));
551 gko_y = std::unique_ptr<vec>(
552 new VectorWrapper(executor, y.Size(), &y,
553 false));
554 }
555
556 // Create the logger object to log some data from the solvers to confirm
557 // convergence.
558 initialize_ginkgo_log(gko_x.get());
559
560 MFEM_VERIFY(convergence_logger, "convergence logger not initialized" );
561 if (print_level==1)
562 {
563 MFEM_VERIFY(residual_logger, "residual logger not initialized" );
564 solver->clear_loggers(); // Clear any loggers from previous Mult() calls
565 solver->add_logger(residual_logger);
566 }
567
568 // Add the convergence logger object to the combined factory to retrieve the
569 // solver and other data
570 combined_factory->clear_loggers();
572
573 // Finally, apply the solver to x and get the solution in y.
574#if MFEM_GINKGO_VERSION < 10600
575 solver->apply(gko::lend(gko_x), gko::lend(gko_y));
576#else
577 solver->apply(gko_x, gko_y);
578#endif
579
580 // Get the number of iterations taken to converge to the solution.
581 final_iter = convergence_logger->get_num_iterations();
582
583 // Some residual norm and convergence print outs.
584 real_t final_res_norm = 0.0;
585
586 // The convergence_logger object contains the residual vector after the
587 // solver has returned. use this vector to compute the residual norm of the
588 // solution. Get the residual norm from the logger. As the convergence logger
589 // returns a `linop`, it is necessary to convert it to a Dense matrix.
590 // Additionally, if the logger is logging on the gpu, it is necessary to copy
591 // the data to the host and hence the `residual_norm_d_master`
592 auto residual_norm = convergence_logger->get_residual_norm();
593 auto imp_residual_norm = convergence_logger->get_implicit_sq_resnorm();
594
596 {
597 auto imp_residual_norm_d =
598 gko::as<gko::matrix::Dense<real_t>>(imp_residual_norm);
599 auto imp_residual_norm_d_master =
600 gko::matrix::Dense<real_t>::create(executor->get_master(),
601 gko::dim<2> {1, 1});
602 imp_residual_norm_d_master->copy_from(imp_residual_norm_d);
603
604 final_res_norm = imp_residual_norm_d_master->at(0,0);
605 }
606 else
607 {
608 auto residual_norm_d =
609 gko::as<gko::matrix::Dense<real_t>>(residual_norm);
610 auto residual_norm_d_master =
611 gko::matrix::Dense<real_t>::create(executor->get_master(),
612 gko::dim<2> {1, 1});
613 residual_norm_d_master->copy_from(residual_norm_d);
614
615 final_res_norm = residual_norm_d_master->at(0,0);
616 }
617
618 converged = 0;
619 if (convergence_logger->has_converged())
620 {
621 converged = 1;
622 }
623
624 if (print_level == 1)
625 {
626 residual_logger->write();
627 }
628 if (converged == 0)
629 {
630 mfem::err << "No convergence!" << '\n';
631 }
632 if (print_level >=2 && converged==1 )
633 {
634 mfem::out << "Converged in " << final_iter <<
635 " iterations with final residual norm "
636 << final_res_norm << '\n';
637 }
638}
639
641{
642
643 if (system_oper)
644 {
645 // If the solver currently needs VectorWrappers, but not due to a
646 // "sub-operator" (preconditioner or inner solver), then it's
647 // because the current system_oper needs them. Reset the property
648 // to false in case the new op is a SparseMatrix.
649 if (needs_wrapped_vecs == true && sub_op_needs_wrapped_vecs == false)
650 {
651 needs_wrapped_vecs = false;
652 }
653 // Reset the pointer
654 system_oper.reset();
655 // Reset the solver generated for the previous operator
656 solver.reset();
657 }
658
659 // Check for SparseMatrix:
660 SparseMatrix *op_mat = const_cast<SparseMatrix*>(
661 dynamic_cast<const SparseMatrix*>(&op));
662 if (op_mat != NULL)
663 {
664 // Needs to be a square matrix
665 MFEM_VERIFY(op_mat->Height() == op_mat->Width(),
666 "System matrix is not square");
667
668 bool on_device = false;
669 if (executor != executor->get_master())
670 {
671 on_device = true;
672 }
673
674 using mtx = gko::matrix::Csr<real_t, int>;
675 const int nnz = op_mat->GetMemoryData().Capacity();
676 system_oper = mtx::create(
677 executor, gko::dim<2>(op_mat->Height(), op_mat->Width()),
679 nnz,
680 op_mat->ReadWriteData(on_device)),
682 nnz,
683 op_mat->ReadWriteJ(on_device)),
684 gko_array<int>::view(executor, op_mat->Height() + 1,
685 op_mat->ReadWriteI(on_device)));
686
687 }
688 else
689 {
690 needs_wrapped_vecs = true;
691 system_oper = std::shared_ptr<OperatorWrapper>(
692 new OperatorWrapper(executor, op.Height(), &op));
693 }
694
695 // Set MFEM Solver size values
696 height = op.Height();
697 width = op.Width();
698
699 // Generate the solver from the solver using the system matrix or operator.
700 solver = solver_gen->generate(system_oper);
701}
702
703/* ---------------------- CGSolver ------------------------ */
705 : EnableGinkgoSolver(exec, true)
706{
707 using cg = gko::solver::Cg<real_t>;
708 this->solver_gen =
709 cg::build().with_criteria(this->combined_factory).on(this->executor);
710}
711
713 const GinkgoPreconditioner &preconditioner)
714 : EnableGinkgoSolver(exec, true)
715{
716 using cg = gko::solver::Cg<real_t>;
717 // Check for a previously-generated preconditioner (for a specific matrix)
718 if (preconditioner.HasGeneratedPreconditioner())
719 {
720 this->solver_gen = cg::build()
721 .with_criteria(this->combined_factory)
722 .with_generated_preconditioner(
723 preconditioner.GetGeneratedPreconditioner())
724 .on(this->executor);
725 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
726 GetGeneratedPreconditioner().get()))
727 {
728 this->sub_op_needs_wrapped_vecs = true;
729 this->needs_wrapped_vecs = true;
730 }
731 }
732 else // Pass a preconditioner factory (will use same matrix as the solver)
733 {
734 this->solver_gen = cg::build()
735 .with_criteria(this->combined_factory)
736 .with_preconditioner(preconditioner.GetFactory())
737 .on(this->executor);
738 }
739}
740
741
742/* ---------------------- BICGSTABSolver ------------------------ */
744 : EnableGinkgoSolver(exec, true)
745{
746 using bicgstab = gko::solver::Bicgstab<real_t>;
747 this->solver_gen = bicgstab::build()
748 .with_criteria(this->combined_factory)
749 .on(this->executor);
750}
751
753 const GinkgoPreconditioner &preconditioner)
754 : EnableGinkgoSolver(exec, true)
755{
756 using bicgstab = gko::solver::Bicgstab<real_t>;
757 if (preconditioner.HasGeneratedPreconditioner())
758 {
759 this->solver_gen = bicgstab::build()
760 .with_criteria(this->combined_factory)
761 .with_generated_preconditioner(
762 preconditioner.GetGeneratedPreconditioner())
763 .on(this->executor);
764 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
765 GetGeneratedPreconditioner().get()))
766 {
767 this->sub_op_needs_wrapped_vecs = true;
768 this->needs_wrapped_vecs = true;
769 }
770 }
771 else
772 {
773 this->solver_gen = bicgstab::build()
774 .with_criteria(this->combined_factory)
775 .with_preconditioner(preconditioner.GetFactory())
776 .on(this->executor);
777 }
778}
779
780
781/* ---------------------- CGSSolver ------------------------ */
783 : EnableGinkgoSolver(exec, true)
784{
785 using cgs = gko::solver::Cgs<real_t>;
786 this->solver_gen =
787 cgs::build().with_criteria(this->combined_factory).on(this->executor);
788}
789
791 const GinkgoPreconditioner &preconditioner)
792 : EnableGinkgoSolver(exec, true)
793{
794 using cgs = gko::solver::Cgs<real_t>;
795 if (preconditioner.HasGeneratedPreconditioner())
796 {
797 this->solver_gen = cgs::build()
798 .with_criteria(this->combined_factory)
799 .with_generated_preconditioner(
800 preconditioner.GetGeneratedPreconditioner())
801 .on(this->executor);
802 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
803 GetGeneratedPreconditioner().get()))
804 {
805 this->sub_op_needs_wrapped_vecs = true;
806 this->needs_wrapped_vecs = true;
807 }
808 }
809 else
810 {
811 this->solver_gen = cgs::build()
812 .with_criteria(this->combined_factory)
813 .with_preconditioner(preconditioner.GetFactory())
814 .on(this->executor);
815 }
816}
817
818
819/* ---------------------- FCGSolver ------------------------ */
821 : EnableGinkgoSolver(exec, true)
822{
823 using fcg = gko::solver::Fcg<real_t>;
824 this->solver_gen =
825 fcg::build().with_criteria(this->combined_factory).on(this->executor);
826}
827
829 const GinkgoPreconditioner &preconditioner)
830 : EnableGinkgoSolver(exec, true)
831{
832 using fcg = gko::solver::Fcg<real_t>;
833 if (preconditioner.HasGeneratedPreconditioner())
834 {
835 this->solver_gen = fcg::build()
836 .with_criteria(this->combined_factory)
837 .with_generated_preconditioner(
838 preconditioner.GetGeneratedPreconditioner())
839 .on(this->executor);
840 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
841 GetGeneratedPreconditioner().get()))
842 {
843 this->sub_op_needs_wrapped_vecs = true;
844 this->needs_wrapped_vecs = true;
845 }
846 }
847 else
848 {
849 this->solver_gen = fcg::build()
850 .with_criteria(this->combined_factory)
851 .with_preconditioner(preconditioner.GetFactory())
852 .on(this->executor);
853 }
854}
855
856
857/* ---------------------- GMRESSolver ------------------------ */
859 : EnableGinkgoSolver(exec, false),
860 m{dim}
861{
862 using gmres = gko::solver::Gmres<real_t>;
863 if (this->m == 0) // Don't set a dimension, but let Ginkgo use its default
864 {
865 this->solver_gen = gmres::build()
866 .with_criteria(this->combined_factory)
867 .on(this->executor);
868 }
869 else
870 {
871 this->solver_gen = gmres::build()
872 .with_krylov_dim(static_cast<unsigned long>(m))
873 .with_criteria(this->combined_factory)
874 .on(this->executor);
875 }
876}
877
879 const GinkgoPreconditioner &preconditioner, int dim)
880 : EnableGinkgoSolver(exec, false),
881 m{dim}
882{
883 using gmres = gko::solver::Gmres<real_t>;
884 // Check for a previously-generated preconditioner (for a specific matrix)
885 if (this->m == 0) // Don't set a dimension, but let Ginkgo use its default
886 {
887 if (preconditioner.HasGeneratedPreconditioner())
888 {
889 this->solver_gen = gmres::build()
890 .with_criteria(this->combined_factory)
891 .with_generated_preconditioner(
892 preconditioner.GetGeneratedPreconditioner())
893 .on(this->executor);
894 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
895 GetGeneratedPreconditioner().get()))
896 {
897 this->sub_op_needs_wrapped_vecs = true;
898 this->needs_wrapped_vecs = true;
899 }
900 }
901 else
902 {
903 this->solver_gen = gmres::build()
904 .with_criteria(this->combined_factory)
905 .with_preconditioner(preconditioner.GetFactory())
906 .on(this->executor);
907 }
908 }
909 else
910 {
911 if (preconditioner.HasGeneratedPreconditioner())
912 {
913 this->solver_gen = gmres::build()
914 .with_krylov_dim(static_cast<unsigned long>(m))
915 .with_criteria(this->combined_factory)
916 .with_generated_preconditioner(
917 preconditioner.GetGeneratedPreconditioner())
918 .on(this->executor);
919 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
920 GetGeneratedPreconditioner().get()))
921 {
922 this->sub_op_needs_wrapped_vecs = true;
923 this->needs_wrapped_vecs = true;
924 }
925 }
926 else
927 {
928 this->solver_gen = gmres::build()
929 .with_krylov_dim(static_cast<unsigned long>(m))
930 .with_criteria(this->combined_factory)
931 .with_preconditioner(preconditioner.GetFactory())
932 .on(this->executor);
933 }
934 }
935}
936
938{
939 m = dim;
940 using gmres = gko::solver::Gmres<real_t>;
941 // Create new solver factory with other parameters the same, but new value for krylov_dim
942 auto current_params = gko::as<gmres::Factory>(solver_gen)->get_parameters();
943 this->solver_gen = current_params.with_krylov_dim(static_cast<unsigned long>(m))
944 .on(this->executor);
945 if (solver)
946 {
947 gko::as<gmres>(solver)->set_krylov_dim(static_cast<unsigned long>(m));
948 }
949}
950
951/* ---------------------- CBGMRESSolver ------------------------ */
953 storage_precision prec)
954 : EnableGinkgoSolver(exec, false),
955 m{dim}
956{
957 using gmres = gko::solver::CbGmres<real_t>;
958 if (this->m == 0) // Don't set a dimension, but let Ginkgo use its default
959 {
960 this->solver_gen = gmres::build()
961 .with_criteria(this->combined_factory)
962 .with_storage_precision(prec)
963 .on(this->executor);
964 }
965 else
966 {
967 this->solver_gen = gmres::build()
968 .with_krylov_dim(static_cast<unsigned long>(m))
969 .with_criteria(this->combined_factory)
970 .with_storage_precision(prec)
971 .on(this->executor);
972 }
973}
974
976 const GinkgoPreconditioner &preconditioner,
977 int dim, storage_precision prec)
978 : EnableGinkgoSolver(exec, false),
979 m{dim}
980{
981 using gmres = gko::solver::CbGmres<real_t>;
982 // Check for a previously-generated preconditioner (for a specific matrix)
983 if (this->m == 0) // Don't set a dimension, but let Ginkgo use its default
984 {
985 if (preconditioner.HasGeneratedPreconditioner())
986 {
987 this->solver_gen = gmres::build()
988 .with_criteria(this->combined_factory)
989 .with_storage_precision(prec)
990 .with_generated_preconditioner(
991 preconditioner.GetGeneratedPreconditioner())
992 .on(this->executor);
993 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
994 GetGeneratedPreconditioner().get()))
995 {
996 this->sub_op_needs_wrapped_vecs = true;
997 this->needs_wrapped_vecs = true;
998 }
999 }
1000 else
1001 {
1002 this->solver_gen = gmres::build()
1003 .with_criteria(this->combined_factory)
1004 .with_storage_precision(prec)
1005 .with_preconditioner(preconditioner.GetFactory())
1006 .on(this->executor);
1007 }
1008 }
1009 else
1010 {
1011 if (preconditioner.HasGeneratedPreconditioner())
1012 {
1013 this->solver_gen = gmres::build()
1014 .with_krylov_dim(static_cast<unsigned long>(m))
1015 .with_criteria(this->combined_factory)
1016 .with_storage_precision(prec)
1017 .with_generated_preconditioner(
1018 preconditioner.GetGeneratedPreconditioner())
1019 .on(this->executor);
1020 if (dynamic_cast<const OperatorWrapper*>(preconditioner.
1021 GetGeneratedPreconditioner().get()))
1022 {
1023 this->sub_op_needs_wrapped_vecs = true;
1024 this->needs_wrapped_vecs = true;
1025 }
1026 }
1027 else
1028 {
1029 this->solver_gen = gmres::build()
1030 .with_krylov_dim(static_cast<unsigned long>(m))
1031 .with_criteria(this->combined_factory)
1032 .with_storage_precision(prec)
1033 .with_preconditioner(preconditioner.GetFactory())
1034 .on(this->executor);
1035 }
1036 }
1037}
1038
1040{
1041 m = dim;
1042 using gmres = gko::solver::CbGmres<real_t>;
1043 // Create new solver factory with other parameters the same, but new value for krylov_dim
1044 auto current_params = gko::as<gmres::Factory>(solver_gen)->get_parameters();
1045 this->solver_gen = current_params.with_krylov_dim(static_cast<unsigned long>(m))
1046 .on(this->executor);
1047 if (solver)
1048 {
1049 gko::as<gmres>(solver)->set_krylov_dim(static_cast<unsigned long>(m));
1050 }
1051}
1052
1053/* ---------------------- IRSolver ------------------------ */
1055 : EnableGinkgoSolver(exec, false)
1056{
1057 using ir = gko::solver::Ir<real_t>;
1058 this->solver_gen =
1059 ir::build().with_criteria(this->combined_factory).on(this->executor);
1060}
1061
1063 const GinkgoIterativeSolver &inner_solver)
1064 : EnableGinkgoSolver(exec, false)
1065{
1066 using ir = gko::solver::Ir<real_t>;
1067 this->solver_gen = ir::build()
1068 .with_criteria(this->combined_factory)
1069 .with_solver(inner_solver.GetFactory())
1070 .on(this->executor);
1071 if (inner_solver.UsesVectorWrappers())
1072 {
1073 this->sub_op_needs_wrapped_vecs = true;
1074 this->needs_wrapped_vecs = true;
1075 }
1076}
1077
1078/* --------------------------------------------------------------- */
1079/* ---------------------- Preconditioners ------------------------ */
1087
1088void
1090{
1091
1092 MFEM_VERIFY(generated_precond, "Preconditioner not initialized");
1093 MFEM_VERIFY(executor, "executor is not initialized");
1094
1095 using vec = gko::matrix::Dense<real_t>;
1096 if (!iterative_mode)
1097 {
1098 y = 0.0;
1099 }
1100
1101 // Create x and y vectors in Ginkgo's format. Wrap MFEM's data directly,
1102 // on CPU or GPU.
1103 bool on_device = false;
1104 if (executor != executor->get_master())
1105 {
1106 on_device = true;
1107 }
1108 auto gko_x = vec::create(executor, gko::dim<2>(x.Size(), 1),
1110 x.Size(), const_cast<real_t *>(
1111 x.Read(on_device))), 1);
1112 auto gko_y = vec::create(executor, gko::dim<2>(y.Size(), 1),
1114 y.Size(),
1115 y.ReadWrite(on_device)), 1);
1116#if MFEM_GINKGO_VERSION < 10600
1117 generated_precond.get()->apply(gko::lend(gko_x), gko::lend(gko_y));
1118#else
1119 generated_precond.get()->apply(gko_x, gko_y);
1120#endif
1121}
1122
1124{
1125
1127 {
1128 generated_precond.reset();
1129 has_generated_precond = false;
1130 }
1131
1132 // Only accept SparseMatrix for this type.
1133 SparseMatrix *op_mat = const_cast<SparseMatrix*>(
1134 dynamic_cast<const SparseMatrix*>(&op));
1135 MFEM_VERIFY(op_mat != NULL,
1136 "GinkgoPreconditioner::SetOperator : not a SparseMatrix!");
1137
1138 bool on_device = false;
1139 if (executor != executor->get_master())
1140 {
1141 on_device = true;
1142 }
1143
1144 using mtx = gko::matrix::Csr<real_t, int>;
1145 const int nnz = op_mat->GetMemoryData().Capacity();
1146 auto gko_matrix = mtx::create(
1147 executor, gko::dim<2>(op_mat->Height(), op_mat->Width()),
1149 nnz,
1150 op_mat->ReadWriteData(on_device)),
1152 nnz,
1153 op_mat->ReadWriteJ(on_device)),
1154 gko_array<int>::view(executor, op_mat->Height() + 1,
1155 op_mat->ReadWriteI(on_device)));
1156
1157 generated_precond = precond_gen->generate(gko::give(gko_matrix));
1158 has_generated_precond = true;
1159
1160 // Set MFEM Solver size values
1161 height = op.Height();
1162 width = op.Width();
1163}
1164
1165
1166/* ---------------------- JacobiPreconditioner ------------------------ */
1168 GinkgoExecutor &exec,
1169 const std::string &storage_opt,
1170 const real_t accuracy,
1171 const int max_block_size
1172)
1173 : GinkgoPreconditioner(exec)
1174{
1175
1176 if (storage_opt == "auto")
1177 {
1178 precond_gen = gko::preconditioner::Jacobi<real_t, int>::build()
1179 .with_storage_optimization(
1180 gko::precision_reduction::autodetect())
1181 .with_accuracy(accuracy)
1182 .with_max_block_size(static_cast<unsigned int>(max_block_size))
1183 .on(executor);
1184 }
1185 else
1186 {
1187 precond_gen = gko::preconditioner::Jacobi<real_t, int>::build()
1188 .with_storage_optimization(
1189 gko::precision_reduction(0, 0))
1190 .with_accuracy(accuracy)
1191 .with_max_block_size(static_cast<unsigned int>(max_block_size))
1192 .on(executor);
1193 }
1194
1195}
1196
1197/* ---------------------- Ilu/IluIsaiPreconditioner ------------------------ */
1199 GinkgoExecutor &exec,
1200 const std::string &factorization_type,
1201 const int sweeps,
1202 const bool skip_sort
1203)
1204 : GinkgoPreconditioner(exec)
1205{
1206 if (factorization_type == "exact")
1207 {
1208 using ilu_fact_type = gko::factorization::Ilu<real_t, int>;
1209 std::shared_ptr<ilu_fact_type::Factory> fact_factory =
1210 ilu_fact_type::build()
1211 .with_skip_sorting(skip_sort)
1212 .on(executor);
1213 precond_gen = gko::preconditioner::Ilu<>::build()
1214#if MFEM_GINKGO_VERSION < 10700
1215 .with_factorization_factory(fact_factory)
1216#else
1217 .with_factorization(fact_factory)
1218#endif
1219 .on(executor);
1220 }
1221 else
1222 {
1223 using ilu_fact_type = gko::factorization::ParIlu<real_t, int>;
1224 std::shared_ptr<ilu_fact_type::Factory> fact_factory =
1225 ilu_fact_type::build()
1226 .with_iterations(static_cast<unsigned long>(sweeps))
1227 .with_skip_sorting(skip_sort)
1228 .on(executor);
1229 precond_gen = gko::preconditioner::Ilu<>::build()
1230#if MFEM_GINKGO_VERSION < 10700
1231 .with_factorization_factory(fact_factory)
1232#else
1233 .with_factorization(fact_factory)
1234#endif
1235 .on(executor);
1236 }
1237
1238}
1239
1241 GinkgoExecutor &exec,
1242 const std::string &factorization_type,
1243 const int sweeps,
1244 const int sparsity_power,
1245 const bool skip_sort
1246)
1247 : GinkgoPreconditioner(exec)
1248{
1249 using l_solver_type = gko::preconditioner::LowerIsai<>;
1250 using u_solver_type = gko::preconditioner::UpperIsai<>;
1251
1252 std::shared_ptr<l_solver_type::Factory> l_solver_factory =
1253 l_solver_type::build()
1254 .with_sparsity_power(sparsity_power)
1255 .on(executor);
1256 std::shared_ptr<u_solver_type::Factory> u_solver_factory =
1257 u_solver_type::build()
1258 .with_sparsity_power(sparsity_power)
1259 .on(executor);
1260
1261
1262
1263 if (factorization_type == "exact")
1264 {
1265 using ilu_fact_type = gko::factorization::Ilu<real_t, int>;
1266 std::shared_ptr<ilu_fact_type::Factory> fact_factory =
1267 ilu_fact_type::build()
1268 .with_skip_sorting(skip_sort)
1269 .on(executor);
1270 precond_gen = gko::preconditioner::Ilu<l_solver_type,
1271 u_solver_type>::build()
1272#if MFEM_GINKGO_VERSION < 10700
1273 .with_factorization_factory(fact_factory)
1274 .with_l_solver_factory(l_solver_factory)
1275 .with_u_solver_factory(u_solver_factory)
1276#else
1277 .with_factorization(fact_factory)
1278 .with_l_solver(l_solver_factory)
1279 .with_u_solver(u_solver_factory)
1280#endif
1281 .on(executor);
1282
1283 }
1284 else
1285 {
1286 using ilu_fact_type = gko::factorization::ParIlu<real_t, int>;
1287 std::shared_ptr<ilu_fact_type::Factory> fact_factory =
1288 ilu_fact_type::build()
1289 .with_iterations(static_cast<unsigned long>(sweeps))
1290 .with_skip_sorting(skip_sort)
1291 .on(executor);
1292 precond_gen = gko::preconditioner::Ilu<l_solver_type,
1293 u_solver_type>::build()
1294#if MFEM_GINKGO_VERSION < 10700
1295 .with_factorization_factory(fact_factory)
1296 .with_l_solver_factory(l_solver_factory)
1297 .with_u_solver_factory(u_solver_factory)
1298#else
1299 .with_factorization(fact_factory)
1300 .with_l_solver(l_solver_factory)
1301 .with_u_solver(u_solver_factory)
1302#endif
1303 .on(executor);
1304 }
1305}
1306
1307
1308/* ---------------------- Ic/IcIsaiPreconditioner ------------------------ */
1310 GinkgoExecutor &exec,
1311 const std::string &factorization_type,
1312 const int sweeps,
1313 const bool skip_sort
1314)
1315 : GinkgoPreconditioner(exec)
1316{
1317
1318 if (factorization_type == "exact")
1319 {
1320 using ic_fact_type = gko::factorization::Ic<real_t, int>;
1321 std::shared_ptr<ic_fact_type::Factory> fact_factory =
1322 ic_fact_type::build()
1323 .with_both_factors(false)
1324 .with_skip_sorting(skip_sort)
1325 .on(executor);
1326 precond_gen = gko::preconditioner::Ic<>::build()
1327#if MFEM_GINKGO_VERSION < 10700
1328 .with_factorization_factory(fact_factory)
1329#else
1330 .with_factorization(fact_factory)
1331#endif
1332 .on(executor);
1333 }
1334 else
1335 {
1336 using ic_fact_type = gko::factorization::ParIc<real_t, int>;
1337 std::shared_ptr<ic_fact_type::Factory> fact_factory =
1338 ic_fact_type::build()
1339 .with_both_factors(false)
1340 .with_iterations(static_cast<unsigned long>(sweeps))
1341 .with_skip_sorting(skip_sort)
1342 .on(executor);
1343 precond_gen = gko::preconditioner::Ic<>::build()
1344#if MFEM_GINKGO_VERSION < 10700
1345 .with_factorization_factory(fact_factory)
1346#else
1347 .with_factorization(fact_factory)
1348#endif
1349 .on(executor);
1350 }
1351}
1352
1354 GinkgoExecutor &exec,
1355 const std::string &factorization_type,
1356 const int sweeps,
1357 const int sparsity_power,
1358 const bool skip_sort
1359)
1360 : GinkgoPreconditioner(exec)
1361{
1362
1363 using l_solver_type = gko::preconditioner::LowerIsai<>;
1364 std::shared_ptr<l_solver_type::Factory> l_solver_factory =
1365 l_solver_type::build()
1366 .with_sparsity_power(sparsity_power)
1367 .on(executor);
1368 if (factorization_type == "exact")
1369 {
1370 using ic_fact_type = gko::factorization::Ic<real_t, int>;
1371 std::shared_ptr<ic_fact_type::Factory> fact_factory =
1372 ic_fact_type::build()
1373 .with_both_factors(false)
1374 .with_skip_sorting(skip_sort)
1375 .on(executor);
1376 precond_gen = gko::preconditioner::Ic<l_solver_type>::build()
1377#if MFEM_GINKGO_VERSION < 10700
1378 .with_factorization_factory(fact_factory)
1379 .with_l_solver_factory(l_solver_factory)
1380#else
1381 .with_factorization(fact_factory)
1382 .with_l_solver(l_solver_factory)
1383#endif
1384 .on(executor);
1385 }
1386 else
1387 {
1388 using ic_fact_type = gko::factorization::ParIc<real_t, int>;
1389 std::shared_ptr<ic_fact_type::Factory> fact_factory =
1390 ic_fact_type::build()
1391 .with_both_factors(false)
1392 .with_iterations(static_cast<unsigned long>(sweeps))
1393 .with_skip_sorting(skip_sort)
1394 .on(executor);
1395 precond_gen = gko::preconditioner::Ic<l_solver_type>::build()
1396#if MFEM_GINKGO_VERSION < 10700
1397 .with_factorization_factory(fact_factory)
1398 .with_l_solver_factory(l_solver_factory)
1399#else
1400 .with_factorization(fact_factory)
1401 .with_l_solver(l_solver_factory)
1402#endif
1403 .on(executor);
1404 }
1405}
1406
1407/* ---------------------- MFEMPreconditioner ------------------------ */
1409 GinkgoExecutor &exec,
1410 const Solver &mfem_precond
1411)
1412 : GinkgoPreconditioner(exec)
1413{
1414 generated_precond = std::shared_ptr<OperatorWrapper>(
1416 mfem_precond.Height(), &mfem_precond));
1417 has_generated_precond = true;
1418}
1419
1420} // namespace Ginkgo
1421
1422} // namespace mfem
1423
1424#endif // MFEM_USE_GINKGO
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
Definition device.hpp:259
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
FCGSolver(GinkgoExecutor &exec)
Definition ginkgo.cpp:820
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
@ 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
std::shared_ptr< gko::LinOpFactory > solver_gen
Definition ginkgo.hpp:755
std::shared_ptr< gko::stop::ResidualNorm<>::Factory > rel_criterion
Definition ginkgo.hpp:767
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
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
void apply_impl(const gko::LinOp *b, gko::LinOp *x) const override
Definition ginkgo.cpp:423
const Vector & get_mfem_vec_const_ref() const
Definition ginkgo.hpp:122
int Capacity() const
Return the size of the allocated memory.
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
Abstract operator.
Definition operator.hpp:25
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
Base class for solvers.
Definition operator.hpp:683
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
Data type sparse matrix.
Definition sparsemat.hpp:51
int * ReadWriteI(bool on_dev=true)
real_t * ReadWriteData(bool on_dev=true)
int * ReadWriteJ(bool on_dev=true)
Memory< real_t > & GetMemoryData()
Vector data type.
Definition vector.hpp:80
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:474
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:490
Memory< real_t > & GetMemory()
Return a reference to the Memory object used by the Vector.
Definition vector.hpp:249
void Destroy()
Destroy a vector.
Definition vector.hpp:615
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:136
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:247
Vector beta
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
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
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition globals.hpp:71
float real_t
Definition config.hpp:43
@ HIP_MASK
Biwise-OR of all HIP backends.
Definition device.hpp:91
@ CUDA_MASK
Biwise-OR of all CUDA backends.
Definition device.hpp:89
@ OMP_MASK
Biwise-OR of all OpenMP backends.
Definition device.hpp:93