MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ginkgo.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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"
18 #include "../general/globals.hpp"
19 #include "../general/error.hpp"
20 #include <iostream>
21 #include <iomanip>
22 #include <algorithm>
23 #include <cmath>
24 
25 namespace mfem
26 {
27 
28 namespace GinkgoWrappers
29 {
30 
32  const std::string &exec_type, int print_iter, int max_num_iter,
33  double RTOLERANCE, double ATOLERANCE)
34  : exec_type(exec_type),
35  print_lvl(print_iter),
36  max_iter(max_num_iter),
37  rel_tol(sqrt(RTOLERANCE)),
38  abs_tol(sqrt(ATOLERANCE))
39 {
40  if (exec_type == "reference")
41  {
42  executor = gko::ReferenceExecutor::create();
43  }
44  else if (exec_type == "omp")
45  {
46  executor = gko::OmpExecutor::create();
47  }
48  else if (exec_type == "cuda" && gko::CudaExecutor::get_num_devices() > 0)
49  {
50  executor = gko::CudaExecutor::create(0, gko::OmpExecutor::create());
51  }
52  else
53  {
54  mfem::err <<
55  " exec_type needs to be one of the three strings: \"reference\", \"cuda\" or \"omp\" "
56  << std::endl;
57  }
58  using ResidualCriterionFactory = gko::stop::ResidualNormReduction<>;
59  residual_criterion = ResidualCriterionFactory::build()
60  .with_reduction_factor(rel_tol)
61  .on(executor);
62 
64  gko::stop::Combined::build()
65  .with_criteria(residual_criterion,
66  gko::stop::Iteration::build()
67  .with_max_iters(max_iter)
68  .on(executor))
69  .on(executor);
70 }
71 
72 void
73 GinkgoIterativeSolverBase::initialize_ginkgo_log(gko::matrix::Dense<double>* b)
74 {
75  // Add the logger object. See the different masks available in Ginkgo's
76  // documentation
77  convergence_logger = gko::log::Convergence<>::create(
78  executor, gko::log::Logger::criterion_check_completed_mask);
79  residual_logger = std::make_shared<ResidualLogger<>>(executor,
80  gko::lend(system_matrix),b);
81 
82 }
83 
84 void
86  const Vector &rhs)
87 {
88  // some shortcuts.
89  using val_array = gko::Array<double>;
90  using vec = gko::matrix::Dense<double>;
91 
92  MFEM_VERIFY(system_matrix, "System matrix not initialized");
93  MFEM_VERIFY(executor, "executor is not initialized");
94  MFEM_VERIFY(rhs.Size() == solution.Size(),
95  "Mismatching sizes for rhs and solution");
96  // Create the rhs vector in Ginkgo's format.
97  std::vector<double> f(rhs.Size());
98  std::copy(rhs.GetData(), rhs.GetData() + rhs.Size(), f.begin());
99  auto b =
100  vec::create(executor,
101  gko::dim<2>(rhs.Size(), 1),
102  val_array::view(executor->get_master(), rhs.Size(), f.data()),
103  1);
104 
105  // Create the solution vector in Ginkgo's format.
106  std::vector<double> u(solution.Size());
107  std::copy(solution.GetData(), solution.GetData() + solution.Size(), u.begin());
108  auto x = vec::create(executor,
109  gko::dim<2>(solution.Size(), 1),
110  val_array::view(executor->get_master(),
111  solution.Size(),
112  u.data()),
113  1);
114 
115  // Create the logger object to log some data from the solvers to confirm
116  // convergence.
117  initialize_ginkgo_log(gko::lend(b));
118 
119  MFEM_VERIFY(convergence_logger, "convergence logger not initialized" );
120  if (print_lvl==1)
121  {
122  MFEM_VERIFY(residual_logger, "residual logger not initialized" );
123  solver_gen->add_logger(residual_logger);
124  }
125 
126  // Generate the solver from the solver using the system matrix.
127  auto solver = solver_gen->generate(system_matrix);
128 
129  // Add the convergence logger object to the combined factory to retrieve the
130  // solver and other data
132 
133  // Finally, apply the solver to b and get the solution in x.
134  solver->apply(gko::lend(b), gko::lend(x));
135 
136  // The convergence_logger object contains the residual vector after the
137  // solver has returned. use this vector to compute the residual norm of the
138  // solution. Get the residual norm from the logger. As the convergence logger
139  // returns a `linop`, it is necessary to convert it to a Dense matrix.
140  // Additionally, if the logger is logging on the gpu, it is necessary to copy
141  // the data to the host and hence the `residual_norm_d_master`
142  auto residual_norm = convergence_logger->get_residual_norm();
143  auto residual_norm_d =
144  gko::as<gko::matrix::Dense<double>>(residual_norm);
145  auto residual_norm_d_master =
146  gko::matrix::Dense<double>::create(executor->get_master(),
147  gko::dim<2> {1, 1});
148  residual_norm_d_master->copy_from(residual_norm_d);
149 
150  // Get the number of iterations taken to converge to the solution.
151  auto num_iteration = convergence_logger->get_num_iterations();
152 
153  // Ginkgo works with a relative residual norm through its
154  // ResidualNormReduction criterion. Therefore, to get the normalized
155  // residual, we divide by the norm of the rhs.
156  auto b_norm = gko::matrix::Dense<double>::create(executor->get_master(),
157  gko::dim<2> {1, 1});
158  if (executor != executor->get_master())
159  {
160  auto b_master = vec::create(executor->get_master(),
161  gko::dim<2>(rhs.Size(), 1),
162  val_array::view(executor->get_master(),
163  rhs.Size(),
164  f.data()),
165  1);
166  b_master->compute_norm2(b_norm.get());
167  }
168  else
169  {
170  b->compute_norm2(b_norm.get());
171  }
172 
173  MFEM_VERIFY(b_norm.get()->at(0, 0) != 0.0, " rhs norm is zero");
174  // Some residual norm and convergence print outs. As both
175  // `residual_norm_d_master` and `b_norm` are seen as Dense matrices, we use
176  // the `at` function to get the first value here. In case of multiple right
177  // hand sides, this will need to be modified.
178  auto fin_res_norm = std::pow(residual_norm_d_master->at(0,0) / b_norm->at(0,0),
179  2);
180  if (num_iteration==max_iter &&
181  fin_res_norm > rel_tol )
182  {
183  converged = 1;
184  }
185  if (fin_res_norm < rel_tol)
186  {
187  converged =0;
188  }
189  if (print_lvl ==1)
190  {
191  residual_logger->write();
192  }
193  if (converged!=0)
194  {
195  mfem::err << "No convergence!" << '\n';
196  mfem::out << "(B r_N, r_N) = " << fin_res_norm << '\n'
197  << "Number of iterations: " << num_iteration << '\n';
198  }
199  if (print_lvl >=2 && converged==0 )
200  {
201  mfem::out << "Converged in " << num_iteration <<
202  " iterations with final residual norm "
203  << fin_res_norm << '\n';
204  }
205 
206  // Check if the solution is on a CUDA device, if so, copy it over to the
207  // host.
208  if (executor != executor->get_master())
209  {
210  auto x_master = vec::create(executor->get_master(),
211  gko::dim<2>(solution.Size(), 1),
212  val_array::view(executor,
213  solution.Size(),
214  x->get_values()),
215  1);
216  x.reset(x_master.release());
217  }
218  // Finally copy over the solution vector to mfem's solution vector.
219  std::copy(x->get_values(),
220  x->get_values() + solution.Size(),
221  solution.GetData());
222 }
223 
224 void
226  const SparseMatrix *matrix)
227 {
228  // Needs to be a square matrix
229  MFEM_VERIFY(matrix->Height() == matrix->Width(), "System matrix is not square");
230 
231  const int N = matrix->Size();
232  using mtx = gko::matrix::Csr<double, int>;
233  std::shared_ptr<mtx> system_matrix_compute;
234  system_matrix_compute = mtx::create(executor->get_master(),
235  gko::dim<2>(N),
236  matrix->NumNonZeroElems());
237  double *mat_values = system_matrix_compute->get_values();
238  int *mat_row_ptrs = system_matrix_compute->get_row_ptrs();
239  int *mat_col_idxs = system_matrix_compute->get_col_idxs();
240  mat_row_ptrs[0] =0;
241  for (int r=0; r< N; ++r)
242  {
243  const int* col = matrix->GetRowColumns(r);
244  const double * val = matrix->GetRowEntries(r);
245  mat_row_ptrs[r+1] = mat_row_ptrs[r] + matrix->RowSize(r);
246  for (int cj=0; cj < matrix->RowSize(r); cj++ )
247  {
248  mat_values[mat_row_ptrs[r]+cj] = val[cj];
249  mat_col_idxs[mat_row_ptrs[r]+cj] = col[cj];
250  }
251  }
252  system_matrix =
253  mtx::create(executor, gko::dim<2>(N), matrix->NumNonZeroElems());
254  system_matrix->copy_from(system_matrix_compute.get());
255 }
256 
257 void
259  Vector &solution,
260  const Vector &rhs)
261 {
262  initialize(matrix);
263  apply(solution, rhs);
264 }
265 
266 
267 /* ---------------------- CGSolver ------------------------ */
269  const std::string &exec_type,
270  int print_iter,
271  int max_num_iter,
272  double RTOLERANCE,
273  double ATOLERANCE
274 )
275  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
276  ATOLERANCE)
277 {
278  using cg = gko::solver::Cg<double>;
279  this->solver_gen =
280  cg::build().with_criteria(this->combined_factory).on(this->executor);
281 }
282 
284  const std::string &exec_type,
285  int print_iter,
286  int max_num_iter,
287  double RTOLERANCE,
288  double ATOLERANCE,
289  const gko::LinOpFactory* preconditioner
290 )
291  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
292  ATOLERANCE)
293 {
294  using cg = gko::solver::Cg<double>;
295  this->solver_gen = cg::build()
296  .with_criteria(this->combined_factory)
297  .with_preconditioner(preconditioner)
298  .on(this->executor);
299 }
300 
301 
302 /* ---------------------- BICGSTABSolver ------------------------ */
304  const std::string &exec_type,
305  int print_iter,
306  int max_num_iter,
307  double RTOLERANCE,
308  double ATOLERANCE
309 )
310  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
311  ATOLERANCE)
312 {
313  using bicgstab = gko::solver::Bicgstab<double>;
314  this->solver_gen = bicgstab::build()
315  .with_criteria(this->combined_factory)
316  .on(this->executor);
317 }
318 
320  const std::string &exec_type,
321  int print_iter,
322  int max_num_iter,
323  double RTOLERANCE,
324  double ATOLERANCE,
325  const gko::LinOpFactory* preconditioner
326 )
327  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
328  ATOLERANCE)
329 {
330  using bicgstab = gko::solver::Bicgstab<double>;
331  this->solver_gen = bicgstab::build()
332  .with_criteria(this->combined_factory)
333  .with_preconditioner(preconditioner)
334  .on(this->executor);
335 }
336 
337 
338 /* ---------------------- CGSSolver ------------------------ */
340  const std::string &exec_type,
341  int print_iter,
342  int max_num_iter,
343  double RTOLERANCE,
344  double ATOLERANCE
345 )
346  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
347  ATOLERANCE)
348 {
349  using cgs = gko::solver::Cgs<double>;
350  this->solver_gen =
351  cgs::build().with_criteria(this->combined_factory).on(this->executor);
352 }
353 
355  const std::string &exec_type,
356  int print_iter,
357  int max_num_iter,
358  double RTOLERANCE,
359  double ATOLERANCE,
360  const gko::LinOpFactory* preconditioner
361 )
362  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
363  ATOLERANCE)
364 {
365  using cgs = gko::solver::Cgs<double>;
366  this->solver_gen = cgs::build()
367  .with_criteria(this->combined_factory)
368  .with_preconditioner(preconditioner)
369  .on(this->executor);
370 }
371 
372 
373 /* ---------------------- FCGSolver ------------------------ */
375  const std::string &exec_type,
376  int print_iter,
377  int max_num_iter,
378  double RTOLERANCE,
379  double ATOLERANCE
380 )
381  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
382  ATOLERANCE)
383 {
384  using fcg = gko::solver::Fcg<double>;
385  this->solver_gen =
386  fcg::build().with_criteria(this->combined_factory).on(this->executor);
387 }
388 
390  const std::string &exec_type,
391  int print_iter,
392  int max_num_iter,
393  double RTOLERANCE,
394  double ATOLERANCE,
395  const gko::LinOpFactory* preconditioner
396 )
397  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
398  ATOLERANCE)
399 {
400  using fcg = gko::solver::Fcg<double>;
401  this->solver_gen = fcg::build()
402  .with_criteria(this->combined_factory)
403  .with_preconditioner(preconditioner)
404  .on(this->executor);
405 }
406 
407 
408 /* ---------------------- GMRESSolver ------------------------ */
410  const std::string &exec_type,
411  int print_iter,
412  int max_num_iter,
413  double RTOLERANCE,
414  double ATOLERANCE
415 )
416  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
417  ATOLERANCE)
418 {
419  using gmres = gko::solver::Gmres<double>;
420  this->solver_gen = gmres::build()
421  .with_krylov_dim(m)
422  .with_criteria(this->combined_factory)
423  .on(this->executor);
424 }
425 
427  const std::string &exec_type,
428  int print_iter,
429  int max_num_iter,
430  double RTOLERANCE,
431  double ATOLERANCE,
432  const gko::LinOpFactory* preconditioner
433 )
434  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
435  ATOLERANCE)
436 {
437  using gmres = gko::solver::Gmres<double>;
438  this->solver_gen = gmres::build()
439  .with_krylov_dim(m)
440  .with_criteria(this->combined_factory)
441  .with_preconditioner(preconditioner)
442  .on(this->executor);
443 }
444 
445 
446 /* ---------------------- IRSolver ------------------------ */
448  const std::string & exec_type,
449  int print_iter,
450  int max_num_iter,
451  double RTOLERANCE,
452  double ATOLERANCE
453 )
454  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
455  ATOLERANCE)
456 {
457  using ir = gko::solver::Ir<double>;
458  this->solver_gen =
459  ir::build().with_criteria(this->combined_factory).on(this->executor);
460 }
461 
463  const std::string &exec_type,
464  int print_iter,
465  int max_num_iter,
466  double RTOLERANCE,
467  double ATOLERANCE,
468  const gko::LinOpFactory* inner_solver
469 )
470  : GinkgoIterativeSolverBase(exec_type, print_iter, max_num_iter, RTOLERANCE,
471  ATOLERANCE)
472 {
473  using ir = gko::solver::Ir<double>;
474  this->solver_gen = ir::build()
475  .with_criteria(this->combined_factory)
476  .with_solver(inner_solver)
477  .on(this->executor);
478 }
479 
480 
481 } // namespace GinkgoWrappers
482 
483 } // namespace mfem
484 
485 #endif // MFEM_USE_GINKGO
CGSSolver(const std::string &exec_type, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Definition: ginkgo.cpp:339
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:296
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1359
FCGSolver(const std::string &exec_type, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Definition: ginkgo.cpp:374
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:343
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
std::shared_ptr< ResidualLogger<> > residual_logger
Definition: ginkgo.hpp:305
void solve(const SparseMatrix *matrix, Vector &solution, const Vector &rhs)
Definition: ginkgo.cpp:258
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
GMRESSolver(const std::string &exec_type, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Definition: ginkgo.cpp:409
std::shared_ptr< gko::LinOpFactory > solver_gen
Definition: ginkgo.hpp:286
std::shared_ptr< gko::stop::Combined::Factory > combined_factory
Definition: ginkgo.hpp:311
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:357
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
CGSolver(const std::string &exec_type, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Definition: ginkgo.cpp:268
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
std::shared_ptr< gko::stop::ResidualNormReduction<>::Factory > residual_criterion
Definition: ginkgo.hpp:293
double b
Definition: lissajous.cpp:42
void initialize(const SparseMatrix *matrix)
Definition: ginkgo.cpp:225
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
void apply(Vector &solution, const Vector &rhs)
Definition: ginkgo.cpp:85
int Size() const
For backward compatibility, define Size() to be synonym of Height().
Definition: sparsemat.hpp:164
BICGSTABSolver(const std::string &exec_type, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Definition: ginkgo.cpp:303
std::shared_ptr< gko::log::Convergence<> > convergence_logger
Definition: ginkgo.hpp:299
std::shared_ptr< gko::Executor > executor
Definition: ginkgo.hpp:318
Vector data type.
Definition: vector.hpp:51
IRSolver(const std::string &exec_type, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Definition: ginkgo.cpp:447
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
GinkgoIterativeSolverBase(const std::string &exec_type, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Definition: ginkgo.cpp:31
double f(const Vector &p)