]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Hide creation of ginkgo executor from the user.
authorPratik Nayak <pratik.nayak4@gmail.com>
Tue, 5 Feb 2019 14:16:25 +0000 (15:16 +0100)
committerPratik Nayak <pratik.nayak4@gmail.com>
Thu, 25 Apr 2019 08:57:01 +0000 (10:57 +0200)
include/deal.II/lac/ginkgo_solver.h
source/lac/ginkgo_solver.cc
tests/ginkgo/solver.cc

index 3b327511ed8615a14ff016b1ba0b48ab918a3c4c..861cb047f5965b7bc4de0f91acf6bcfb0fb5005a 100644 (file)
@@ -100,7 +100,7 @@ namespace GinkgoWrappers
      * deal.II iterative solvers.
      */
     SolverBase(SolverControl &                solver_control,
-               std::shared_ptr<gko::Executor> executor);
+               std::string exec_type);
 
     /**
      * Destructor.
@@ -194,6 +194,13 @@ namespace GinkgoWrappers
      * @todo Templatize based on Matrix type.
      */
     std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>> 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<gko::Executor> 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<gko::Executor> executor,
+           std::string exec_type,
            std::shared_ptr<gko::LinOpFactory> preconditioner,
            const AdditionalData &         data = AdditionalData());
   protected:
index 2aa25a55ee8a42a88433e6f87e29177df2fd69e9..4dba46b52b1c10d6f1f811aa33ceb800b0500b28 100644 (file)
@@ -33,10 +33,22 @@ namespace GinkgoWrappers
   template <typename ValueType, typename IndexType>
   SolverBase<ValueType, IndexType>::SolverBase(
     SolverControl &                solver_control,
-    std::shared_ptr<gko::Executor> 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 <typename ValueType, typename IndexType>
   SolverCG<ValueType, IndexType>::SolverCG(
     SolverControl &                solver_control,
-    std::shared_ptr<gko::Executor> executor,
+    std::string exec_type,
     const AdditionalData &         data)
-    : SolverBase<ValueType, IndexType>(solver_control, executor)
+    : SolverBase<ValueType, IndexType>(solver_control, exec_type)
     , additional_data(data)
   {
     using cg = gko::solver::Cg<ValueType>;
     this->solver_gen =
-      cg::build().with_criteria(this->combined_factory).on(executor);
+      cg::build().with_criteria(this->combined_factory).on(this->executor);
   }
 
   template <typename ValueType, typename IndexType>
   SolverCG<ValueType, IndexType>::SolverCG(
                                            SolverControl &                solver_control,
-                                           std::shared_ptr<gko::Executor> executor,
+                                           std::string exec_type,
                                            std::shared_ptr<gko::LinOpFactory> preconditioner,
                                            const AdditionalData &         data)
-    : SolverBase<ValueType, IndexType>(solver_control, executor)
+    : SolverBase<ValueType, IndexType>(solver_control, exec_type)
     , additional_data(data)
   {
     using cg = gko::solver::Cg<ValueType>;
     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) \
index 46deda45c81454926dc2c5fdb9d96094937e6ef3..fc0059ce1d91bcc6a4d1721d5dad9aaa76fa427f 100644 (file)
@@ -56,10 +56,11 @@ main(int argc, char **argv)
     // std::shared_ptr<gko::Executor> exec = gko::CudaExecutor::create(0,
     //                                                               gko::OmpExecutor::create());
     // std::shared_ptr<gko::Executor> exec = gko::OmpExecutor::create();
+    auto executor = "reference";
     std::shared_ptr<gko::Executor> exec = gko::ReferenceExecutor::create();
     std::shared_ptr<gko::LinOpFactory> 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(),

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.