22 #ifdef DEAL_II_WITH_GINKGO
33 template <
typename ValueType,
typename IndexType>
35 const std::string &exec_type)
36 : solver_control(solver_control)
37 , exec_type(exec_type)
39 if (exec_type ==
"reference")
41 executor = gko::ReferenceExecutor::create();
43 else if (exec_type ==
"omp")
45 executor = gko::OmpExecutor::create();
47 else if (exec_type ==
"cuda" && gko::CudaExecutor::get_num_devices() > 0)
49 executor = gko::CudaExecutor::create(0, gko::OmpExecutor::create());
56 " exec_type needs to be one of the three strings: \"reference\", \"cuda\" or \"omp\" "));
58 using ResidualCriterionFactory = gko::stop::ResidualNormReduction<>;
59 residual_criterion = ResidualCriterionFactory::build()
60 .with_reduction_factor(solver_control.
tolerance())
64 gko::stop::Combined::build()
65 .with_criteria(residual_criterion,
66 gko::stop::Iteration::build()
67 .with_max_iters(solver_control.
max_steps())
74 template <
typename ValueType,
typename IndexType>
80 convergence_logger = gko::log::Convergence<>::create(
81 executor, gko::log::Logger::criterion_check_completed_mask);
86 template <
typename ValueType,
typename IndexType>
89 const Vector<ValueType> &rhs)
92 using val_array = gko::Array<ValueType>;
93 using vec = gko::matrix::Dense<ValueType>;
97 Assert(rhs.size() == solution.size(),
101 auto solver = solver_gen->generate(system_matrix);
104 std::vector<ValueType> f(rhs.size());
105 std::copy(rhs.begin(), rhs.begin() + rhs.size(), f.begin());
107 vec::create(executor,
108 gko::dim<2>(rhs.size(), 1),
109 val_array::view(executor->get_master(), rhs.size(), f.data()),
113 std::vector<ValueType> u(solution.size());
114 std::copy(solution.begin(), solution.begin() + solution.size(), u.begin());
115 auto x = vec::create(executor,
116 gko::dim<2>(solution.size(), 1),
117 val_array::view(executor->get_master(),
124 initialize_ginkgo_log();
129 combined_factory->add_logger(convergence_logger);
132 solver->apply(gko::lend(
b), gko::lend(x));
141 auto residual_norm = convergence_logger->get_residual_norm();
142 auto residual_norm_d =
143 gko::as<gko::matrix::Dense<ValueType>>(residual_norm);
144 auto residual_norm_d_master =
145 gko::matrix::Dense<ValueType>::create(executor->get_master(),
147 residual_norm_d_master->copy_from(residual_norm_d);
150 auto num_iteration = convergence_logger->get_num_iterations();
155 auto b_norm = gko::matrix::Dense<ValueType>::create(executor->get_master(),
157 if (executor != executor->get_master())
159 auto b_master = vec::create(executor->get_master(),
160 gko::dim<2>(rhs.size(), 1),
161 val_array::view(executor->get_master(),
165 b_master->compute_norm2(b_norm.get());
169 b->compute_norm2(b_norm.get());
178 solver_control.check(num_iteration,
179 residual_norm_d_master->at(0, 0) / b_norm->at(0, 0));
185 solver_control.last_value()));
189 if (executor != executor->get_master())
191 auto x_master = vec::create(executor->get_master(),
192 gko::dim<2>(solution.size(), 1),
193 val_array::view(executor,
197 x.reset(x_master.release());
201 x->get_values() + solution.size(),
207 template <
typename ValueType,
typename IndexType>
211 return solver_control;
216 template <
typename ValueType,
typename IndexType>
227 using mtx = gko::matrix::Csr<ValueType, IndexType>;
228 std::shared_ptr<mtx> system_matrix_compute;
229 system_matrix_compute = mtx::create(executor->get_master(),
231 matrix.n_nonzero_elements());
232 ValueType *mat_values = system_matrix_compute->get_values();
233 IndexType *mat_row_ptrs = system_matrix_compute->get_row_ptrs();
234 IndexType *mat_col_idxs = system_matrix_compute->get_col_idxs();
247 mat_row_ptrs[row - 1] +
matrix.get_row_length(row - 1);
256 std::vector<IndexType> row_pointers(
N + 1);
257 std::copy(system_matrix_compute->get_row_ptrs(),
258 system_matrix_compute->get_row_ptrs() +
N + 1,
259 row_pointers.begin());
271 mat_col_idxs[row_pointers[row]] = p->column();
272 mat_values[row_pointers[row]] = p->value();
284 mtx::create(executor, gko::dim<2>(
N),
matrix.n_nonzero_elements());
285 system_matrix->copy_from(system_matrix_compute.get());
290 template <
typename ValueType,
typename IndexType>
293 Vector<ValueType> & solution,
294 const Vector<ValueType> &rhs)
297 apply(solution, rhs);
303 template <
typename ValueType,
typename IndexType>
305 const std::string & exec_type,
307 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
308 , additional_data(data)
310 using cg = gko::solver::Cg<ValueType>;
317 template <
typename ValueType,
typename IndexType>
320 const std::string & exec_type,
321 const std::shared_ptr<gko::LinOpFactory> &preconditioner,
323 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
324 , additional_data(data)
326 using cg = gko::solver::Cg<ValueType>;
329 .with_preconditioner(preconditioner)
336 template <
typename ValueType,
typename IndexType>
339 const std::string & exec_type,
341 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
342 , additional_data(data)
344 using bicgstab = gko::solver::Bicgstab<ValueType>;
352 template <
typename ValueType,
typename IndexType>
355 const std::string & exec_type,
356 const std::shared_ptr<gko::LinOpFactory> &preconditioner,
358 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
359 , additional_data(data)
361 using bicgstab = gko::solver::Bicgstab<ValueType>;
364 .with_preconditioner(preconditioner)
371 template <
typename ValueType,
typename IndexType>
373 const std::string &exec_type,
375 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
376 , additional_data(data)
378 using cgs = gko::solver::Cgs<ValueType>;
385 template <
typename ValueType,
typename IndexType>
388 const std::string & exec_type,
389 const std::shared_ptr<gko::LinOpFactory> &preconditioner,
391 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
392 , additional_data(data)
394 using cgs = gko::solver::Cgs<ValueType>;
397 .with_preconditioner(preconditioner)
404 template <
typename ValueType,
typename IndexType>
406 const std::string &exec_type,
408 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
409 , additional_data(data)
411 using fcg = gko::solver::Fcg<ValueType>;
418 template <
typename ValueType,
typename IndexType>
421 const std::string & exec_type,
422 const std::shared_ptr<gko::LinOpFactory> &preconditioner,
424 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
425 , additional_data(data)
427 using fcg = gko::solver::Fcg<ValueType>;
430 .with_preconditioner(preconditioner)
437 template <
typename ValueType,
typename IndexType>
439 const unsigned int restart_parameter)
440 : restart_parameter(restart_parameter)
445 template <
typename ValueType,
typename IndexType>
452 using gmres = gko::solver::Gmres<ValueType>;
461 template <
typename ValueType,
typename IndexType>
464 const std::string & exec_type,
465 const std::shared_ptr<gko::LinOpFactory> &preconditioner,
467 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
468 , additional_data(data)
470 using gmres = gko::solver::Gmres<ValueType>;
474 .with_preconditioner(preconditioner)
481 template <
typename ValueType,
typename IndexType>
483 const std::string & exec_type,
485 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
486 , additional_data(data)
488 using ir = gko::solver::Ir<ValueType>;
495 template <
typename ValueType,
typename IndexType>
498 const std::string & exec_type,
499 const std::shared_ptr<gko::LinOpFactory> &inner_solver,
501 :
SolverBase<ValueType, IndexType>(solver_control, exec_type)
502 , additional_data(data)
504 using ir = gko::solver::Ir<ValueType>;
507 .with_solver(inner_solver)
514 # define DEALII_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(_macro) \
515 template _macro(float, int32_t); \
516 template _macro(double, int32_t); \
517 template _macro(float, int64_t); \
518 template _macro(double, int64_t);
520 # define DECLARE_SOLVER_BASE(ValueType, IndexType) \
521 class SolverBase<ValueType, IndexType>
523 # undef DECLARE_SOLVER_BASE
525 # define DECLARE_SOLVER_CG(ValueType, IndexType) \
526 class SolverCG<ValueType, IndexType>
528 # undef DECLARE_SOLVER_CG
530 # define DECLARE_SOLVER_Bicgstab(ValueType, IndexType) \
531 class SolverBicgstab<ValueType, IndexType>
533 # undef DECLARE_SOLVER_Bicgstab
535 # define DECLARE_SOLVER_CGS(ValueType, IndexType) \
536 class SolverCGS<ValueType, IndexType>
538 # undef DECLARE_SOLVER_CGS
540 # define DECLARE_SOLVER_FCG(ValueType, IndexType) \
541 class SolverFCG<ValueType, IndexType>
543 # undef DECLARE_SOLVER_FCG
545 # define DECLARE_SOLVER_GMRES(ValueType, IndexType) \
546 class SolverGMRES<ValueType, IndexType>
548 # undef DECLARE_SOLVER_GMRES
550 # define DECLARE_SOLVER_IR(ValueType, IndexType) \
551 class SolverIR<ValueType, IndexType>
553 # undef DECLARE_SOLVER_IR
560 #endif // DEAL_II_WITH_GINKGO