From e27db51a695ef3cb57fd34ec829d199c80c5a629 Mon Sep 17 00:00:00 2001 From: Pratik Nayak Date: Tue, 5 Feb 2019 15:16:25 +0100 Subject: [PATCH] Hide creation of ginkgo executor from the user. --- include/deal.II/lac/ginkgo_solver.h | 13 ++++++++++--- source/lac/ginkgo_solver.cc | 28 ++++++++++++++++++++-------- tests/ginkgo/solver.cc | 5 +++-- 3 files changed, 33 insertions(+), 13 deletions(-) diff --git a/include/deal.II/lac/ginkgo_solver.h b/include/deal.II/lac/ginkgo_solver.h index 3b327511ed..861cb047f5 100644 --- a/include/deal.II/lac/ginkgo_solver.h +++ b/include/deal.II/lac/ginkgo_solver.h @@ -100,7 +100,7 @@ namespace GinkgoWrappers * deal.II iterative solvers. */ SolverBase(SolverControl & solver_control, - std::shared_ptr executor); + std::string exec_type); /** * Destructor. @@ -194,6 +194,13 @@ namespace GinkgoWrappers * @todo Templatize based on Matrix type. */ std::shared_ptr> system_matrix; + + /** + * The execution paradigm to be set by the user. The choices are between + * `omp`, `cuda` and `reference` + * and more details can be found in Ginkgo's documentation. + */ + std::string exec_type; }; @@ -222,7 +229,7 @@ namespace GinkgoWrappers * @p executor The execution paradigm for the CG solver. */ SolverCG(SolverControl & solver_control, - std::shared_ptr executor, + std::string exec_type, const AdditionalData & data = AdditionalData()); /** @@ -235,7 +242,7 @@ namespace GinkgoWrappers * @p executor The execution paradigm for the CG solver. */ SolverCG(SolverControl & solver_control, - std::shared_ptr executor, + std::string exec_type, std::shared_ptr preconditioner, const AdditionalData & data = AdditionalData()); protected: diff --git a/source/lac/ginkgo_solver.cc b/source/lac/ginkgo_solver.cc index 2aa25a55ee..4dba46b52b 100644 --- a/source/lac/ginkgo_solver.cc +++ b/source/lac/ginkgo_solver.cc @@ -33,10 +33,22 @@ namespace GinkgoWrappers template SolverBase::SolverBase( SolverControl & solver_control, - std::shared_ptr executor) + std::string exec_type) : solver_control(solver_control) - , executor(executor) + , exec_type(exec_type) { + if (exec_type == "reference") { + executor = gko::ReferenceExecutor::create(); + } else if (exec_type == "omp") { + executor = gko::OmpExecutor::create(); + } else if (exec_type == "cuda" && + gko::CudaExecutor::get_num_devices() > 0) { + executor = gko::CudaExecutor::create(0, gko::OmpExecutor::create()); + } else { + std::cerr << "exec_type needs to be one of the three strings: reference, cuda or omp, but provided with" + << exec_type<< std::endl; + std::exit(-1); + } using ResidualCriterionFactory = gko::stop::ResidualNormReduction<>; residual_criterion = ResidualCriterionFactory::build() .with_reduction_factor(solver_control.tolerance()) @@ -278,28 +290,28 @@ namespace GinkgoWrappers template SolverCG::SolverCG( SolverControl & solver_control, - std::shared_ptr executor, + std::string exec_type, const AdditionalData & data) - : SolverBase(solver_control, executor) + : SolverBase(solver_control, exec_type) , additional_data(data) { using cg = gko::solver::Cg; this->solver_gen = - cg::build().with_criteria(this->combined_factory).on(executor); + cg::build().with_criteria(this->combined_factory).on(this->executor); } template SolverCG::SolverCG( SolverControl & solver_control, - std::shared_ptr executor, + std::string exec_type, std::shared_ptr preconditioner, const AdditionalData & data) - : SolverBase(solver_control, executor) + : SolverBase(solver_control, exec_type) , additional_data(data) { using cg = gko::solver::Cg; this->solver_gen = - cg::build().with_criteria(this->combined_factory).with_preconditioner(preconditioner).on(executor); + cg::build().with_criteria(this->combined_factory).with_preconditioner(preconditioner).on(this->executor); } // Explicit instantiations in GinkgoWrappers # define DEALII_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(_macro) \ diff --git a/tests/ginkgo/solver.cc b/tests/ginkgo/solver.cc index 46deda45c8..fc0059ce1d 100644 --- a/tests/ginkgo/solver.cc +++ b/tests/ginkgo/solver.cc @@ -56,10 +56,11 @@ main(int argc, char **argv) // std::shared_ptr exec = gko::CudaExecutor::create(0, // gko::OmpExecutor::create()); // std::shared_ptr exec = gko::OmpExecutor::create(); + auto executor = "reference"; std::shared_ptr exec = gko::ReferenceExecutor::create(); std::shared_ptr jacobi = gko::preconditioner::Jacobi<>::build().on(exec); - GinkgoWrappers::SolverCG<> solver(control, exec); - GinkgoWrappers::SolverCG<> solver_with_jacobi_precond(control, exec, jacobi); + GinkgoWrappers::SolverCG<> solver(control, executor); + GinkgoWrappers::SolverCG<> solver_with_jacobi_precond(control, executor, jacobi); check_solver_within_range(solver.solve(A, u, f), control.last_step(), -- 2.39.5