]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Update ginkgo html docs, clang-format.
authorPratik Nayak <pratik.nayak4@gmail.com>
Fri, 26 Apr 2019 10:04:00 +0000 (12:04 +0200)
committerPratik Nayak <pratik.nayak4@gmail.com>
Fri, 26 Apr 2019 10:04:00 +0000 (12:04 +0200)
doc/external-libs/ginkgo.html
include/deal.II/lac/ginkgo_solver.h
tests/ginkgo/solver.cc

index d918e17f7357e64f24feedba3c4efae9c4ba6f34..e89e3730c096eabd80a718ea33da642bcbc08c83 100644 (file)
                                       target="_top">Ginkgo</a>.
         The different flags that can be added are:
         <ul>
-           <li> <code> -DBUILD_REFERENCE={ON,OFF}</code>: Builds the reference single thread implementation
+           <li> <code> -DGINKGO_BUILD_REFERENCE={ON,OFF}</code>: Builds the reference single thread implementation
         used to check the more sophisticated implementations of the GPU/CPU.
-           <li> <code> -DBUILD_CUDA={ON,OFF}</code>: Builds the GPU implementation, specifically the NVIDIA
+           <li> <code> -DGINKGO_BUILD_CUDA={ON,OFF}</code>: Builds the GPU implementation, specifically the NVIDIA
         CUDA implementations. This needs CUDA to be installed on the machine.
-           <li> <code> -DBUILD_OMP={ON,OFF}</code>: Builds the multithreaded OpenMP implementations.
+           <li> <code> -DGINKGO_BUILD_OMP={ON,OFF}</code>: Builds the multithreaded OpenMP implementations.
         </ul>
         To install Ginkgo, you would need to:
       <pre>
        git clone https://github.com/ginkgo-project/ginkgo.git
-  mkdir build; cd build 
-  cmake -DBUILD_REFERENCE=on/off -DBUILD_CUDA=on/off -DBUILD_OMP=on/off -DCMAKE_INSTALL_PREFIX=/path/to/install/ -DCMAKE_BUILD_TYPE=Debug
+  mkdir build; cd build
+  cmake -DGINKGO_BUILD_REFERENCE=on/off -DGINKGO_BUILD_CUDA=on/off -DGINKGO_BUILD_OMP=on/off -DCMAKE_INSTALL_PREFIX=/path/to/install/ -DCMAKE_BUILD_TYPE=Release/Debug
        make install
       </pre>
     </p>
index eaca09895ef10de253f38f2dcb14a907ac046e9c..9307385b9a7accc0da4415ece47890a3517746b4 100644 (file)
@@ -99,8 +99,7 @@ namespace GinkgoWrappers
      * The @p solver_control object is the same as for other
      * deal.II iterative solvers.
      */
-    SolverBase(SolverControl &                solver_control,
-               std::string exec_type);
+    SolverBase(SolverControl &solver_control, std::string exec_type);
 
     /**
      * Destructor.
@@ -228,25 +227,26 @@ namespace GinkgoWrappers
      *
      * @p exec_type The execution paradigm for the CG solver.
      */
-    SolverCG(SolverControl &                solver_control,
-             std::string exec_type,
-             const AdditionalData &         data = AdditionalData());
+    SolverCG(SolverControl &       solver_control,
+             std::string           exec_type,
+             const AdditionalData &data = AdditionalData());
+
+    /**
+     * Constructor.
+     *
+     * @p solver_control The solver control object is then used to set the
+     * parameters and setup the CG solver from the CG factory which solves the
+     * linear system.
+     *
+     * @p exec_type The execution paradigm for the CG solver.
+     *
+     * @p preconditioner The preconditioner for the solver.
+     */
+    SolverCG(SolverControl &                    solver_control,
+             std::string                        exec_type,
+             std::shared_ptr<gko::LinOpFactory> preconditioner,
+             const AdditionalData &             data = AdditionalData());
 
-  /**
-   * Constructor.
-   *
-   * @p solver_control The solver control object is then used to set the
-   * parameters and setup the CG solver from the CG factory which solves the
-   * linear system.
-   *
-   * @p exec_type The execution paradigm for the CG solver.
-   *
-   * @p preconditioner The preconditioner for the solver.
-   */
-  SolverCG(SolverControl &                solver_control,
-           std::string exec_type,
-           std::shared_ptr<gko::LinOpFactory> preconditioner,
-           const AdditionalData &         data = AdditionalData());
   protected:
     /**
      * Store a copy of the settings for this particular solver.
@@ -274,30 +274,31 @@ namespace GinkgoWrappers
      * Constructor.
      *
      * @p solver_control The solver control object is then used to set the
-     * parameters and setup the BICGSTAB solver from the BICGSTAB factory which solves the
-     * linear system.
+     * parameters and setup the BICGSTAB solver from the BICGSTAB factory which
+     * solves the linear system.
      *
      * @p exec_type The execution paradigm for the BICGSTAB solver.
      */
-    SolverBICGSTAB(SolverControl &                solver_control,
-             std::string exec_type,
-             const AdditionalData &         data = AdditionalData());
+    SolverBICGSTAB(SolverControl &       solver_control,
+                   std::string           exec_type,
+                   const AdditionalData &data = AdditionalData());
+
+    /**
+     * Constructor.
+     *
+     * @p solver_control The solver control object is then used to set the
+     * parameters and setup the BICGSTAB solver from the BICGSTAB factory which
+     * solves the linear system.
+     *
+     * @p exec_type The execution paradigm for the BICGSTAB solver.
+     *
+     * @p preconditioner The preconditioner for the solver.
+     */
+    SolverBICGSTAB(SolverControl &                    solver_control,
+                   std::string                        exec_type,
+                   std::shared_ptr<gko::LinOpFactory> preconditioner,
+                   const AdditionalData &             data = AdditionalData());
 
-  /**
-   * Constructor.
-   *
-   * @p solver_control The solver control object is then used to set the
-   * parameters and setup the BICGSTAB solver from the BICGSTAB factory which solves the
-   * linear system.
-   *
-   * @p exec_type The execution paradigm for the BICGSTAB solver.
-   *
-   * @p preconditioner The preconditioner for the solver.
-   */
-  SolverBICGSTAB(SolverControl &                solver_control,
-           std::string exec_type,
-           std::shared_ptr<gko::LinOpFactory> preconditioner,
-           const AdditionalData &         data = AdditionalData());
   protected:
     /**
      * Store a copy of the settings for this particular solver.
@@ -311,7 +312,7 @@ namespace GinkgoWrappers
    * @ingroup GinkgoWrappers
    */
   template <typename ValueType = double, typename IndexType = int32_t>
-  class SolverCGS: public SolverBase<ValueType, IndexType>
+  class SolverCGS : public SolverBase<ValueType, IndexType>
   {
   public:
     /**
@@ -329,25 +330,26 @@ namespace GinkgoWrappers
      *
      * @p exec_type The execution paradigm for the CGS solver.
      */
-    SolverCGS(SolverControl &                solver_control,
-             std::string exec_type,
-             const AdditionalData &         data = AdditionalData());
+    SolverCGS(SolverControl &       solver_control,
+              std::string           exec_type,
+              const AdditionalData &data = AdditionalData());
+
+    /**
+     * Constructor.
+     *
+     * @p solver_control The solver control object is then used to set the
+     * parameters and setup the CGS solver from the CGS factory which solves the
+     * linear system.
+     *
+     * @p exec_type The execution paradigm for the CGS solver.
+     *
+     * @p preconditioner The preconditioner for the solver.
+     */
+    SolverCGS(SolverControl &                    solver_control,
+              std::string                        exec_type,
+              std::shared_ptr<gko::LinOpFactory> preconditioner,
+              const AdditionalData &             data = AdditionalData());
 
-  /**
-   * Constructor.
-   *
-   * @p solver_control The solver control object is then used to set the
-   * parameters and setup the CGS solver from the CGS factory which solves the
-   * linear system.
-   *
-   * @p exec_type The execution paradigm for the CGS solver.
-   *
-   * @p preconditioner The preconditioner for the solver.
-   */
-  SolverCGS(SolverControl &                solver_control,
-           std::string exec_type,
-           std::shared_ptr<gko::LinOpFactory> preconditioner,
-           const AdditionalData &         data = AdditionalData());
   protected:
     /**
      * Store a copy of the settings for this particular solver.
@@ -361,7 +363,7 @@ namespace GinkgoWrappers
    * @ingroup GinkgoWrappers
    */
   template <typename ValueType = double, typename IndexType = int32_t>
-  class SolverFCG: public SolverBase<ValueType, IndexType>
+  class SolverFCG : public SolverBase<ValueType, IndexType>
   {
   public:
     /**
@@ -379,25 +381,26 @@ namespace GinkgoWrappers
      *
      * @p exec_type The execution paradigm for the FCG solver.
      */
-    SolverFCG(SolverControl &                solver_control,
-             std::string exec_type,
-             const AdditionalData &         data = AdditionalData());
+    SolverFCG(SolverControl &       solver_control,
+              std::string           exec_type,
+              const AdditionalData &data = AdditionalData());
+
+    /**
+     * Constructor.
+     *
+     * @p solver_control The solver control object is then used to set the
+     * parameters and setup the FCG solver from the FCG factory which solves the
+     * linear system.
+     *
+     * @p exec_type The execution paradigm for the FCG solver.
+     *
+     * @p preconditioner The preconditioner for the solver.
+     */
+    SolverFCG(SolverControl &                    solver_control,
+              std::string                        exec_type,
+              std::shared_ptr<gko::LinOpFactory> preconditioner,
+              const AdditionalData &             data = AdditionalData());
 
-  /**
-   * Constructor.
-   *
-   * @p solver_control The solver control object is then used to set the
-   * parameters and setup the FCG solver from the FCG factory which solves the
-   * linear system.
-   *
-   * @p exec_type The execution paradigm for the FCG solver.
-   *
-   * @p preconditioner The preconditioner for the solver.
-   */
-  SolverFCG(SolverControl &                solver_control,
-           std::string exec_type,
-           std::shared_ptr<gko::LinOpFactory> preconditioner,
-           const AdditionalData &         data = AdditionalData());
   protected:
     /**
      * Store a copy of the settings for this particular solver.
@@ -411,55 +414,55 @@ namespace GinkgoWrappers
    * @ingroup GinkgoWrappers
    */
   template <typename ValueType = double, typename IndexType = int32_t>
-  class SolverGMRES: public SolverBase<ValueType, IndexType>
+  class SolverGMRES : public SolverBase<ValueType, IndexType>
   {
   public:
     /**
      * A standardized data struct to pipe additional data to the solver.
      */
     struct AdditionalData
-  {
-    /**
-     * Constructor. By default, set the number of temporary vectors to 30,
-     * i.e. do a restart every 30 iterations.
-     */
-    AdditionalData(const unsigned int restart_parameter     = 30);
+    {
+      /**
+       * Constructor. By default, set the number of temporary vectors to 30,
+       * i.e. do a restart every 30 iterations.
+       */
+      AdditionalData(const unsigned int restart_parameter = 30);
+
+      /**
+       * Maximum number of tmp vectors.
+       */
+      unsigned int restart_parameter;
+    };
 
     /**
-     * Maximum number of tmp vectors.
+     * Constructor.
+     *
+     * @p solver_control The solver control object is then used to set the
+     * parameters and setup the GMRES solver from the GMRES factory which solves
+     * the linear system.
+     *
+     * @p exec_type The execution paradigm for the GMRES solver.
      */
-    unsigned int restart_parameter;
-
-    };
+    SolverGMRES(SolverControl &       solver_control,
+                std::string           exec_type,
+                const AdditionalData &data = AdditionalData());
 
     /**
      * Constructor.
      *
      * @p solver_control The solver control object is then used to set the
-     * parameters and setup the GMRES solver from the GMRES factory which solves the
-     * linear system.
+     * parameters and setup the GMRES solver from the GMRES factory which solves
+     * the linear system.
      *
      * @p exec_type The execution paradigm for the GMRES solver.
+     *
+     * @p preconditioner The preconditioner for the solver.
      */
-    SolverGMRES(SolverControl &                solver_control,
-             std::string exec_type,
-             const AdditionalData &         data = AdditionalData());
+    SolverGMRES(SolverControl &                    solver_control,
+                std::string                        exec_type,
+                std::shared_ptr<gko::LinOpFactory> preconditioner,
+                const AdditionalData &             data = AdditionalData());
 
-  /**
-   * Constructor.
-   *
-   * @p solver_control The solver control object is then used to set the
-   * parameters and setup the GMRES solver from the GMRES factory which solves the
-   * linear system.
-   *
-   * @p exec_type The execution paradigm for the GMRES solver.
-   *
-   * @p preconditioner The preconditioner for the solver.
-   */
-  SolverGMRES(SolverControl &                solver_control,
-           std::string exec_type,
-           std::shared_ptr<gko::LinOpFactory> preconditioner,
-           const AdditionalData &         data = AdditionalData());
   protected:
     /**
      * Store a copy of the settings for this particular solver.
@@ -473,7 +476,7 @@ namespace GinkgoWrappers
    * @ingroup GinkgoWrappers
    */
   template <typename ValueType = double, typename IndexType = int32_t>
-  class SolverIR: public SolverBase<ValueType, IndexType>
+  class SolverIR : public SolverBase<ValueType, IndexType>
   {
   public:
     /**
@@ -491,25 +494,25 @@ namespace GinkgoWrappers
      *
      * @p exec_type The execution paradigm for the IR solver.
      */
-    SolverIR(SolverControl &                solver_control,
-             std::string exec_type,
-             const AdditionalData &         data = AdditionalData());
+    SolverIR(SolverControl &       solver_control,
+             std::string           exec_type,
+             const AdditionalData &data = AdditionalData());
 
-  /**
-   * Constructor.
-   *
-   * @p solver_control The solver control object is then used to set the
-   * parameters and setup the IR solver from the IR factory which solves the
-   * linear system.
-   *
-   * @p exec_type The execution paradigm for the IR solver.
-   *
-   * @p inner_solver The Inner solver for the IR solver.
-   */
-  SolverIR(SolverControl &                solver_control,
-           std::string exec_type,
-           std::shared_ptr<gko::LinOpFactory> inner_solver,
-           const AdditionalData &         data = AdditionalData());
+    /**
+     * Constructor.
+     *
+     * @p solver_control The solver control object is then used to set the
+     * parameters and setup the IR solver from the IR factory which solves the
+     * linear system.
+     *
+     * @p exec_type The execution paradigm for the IR solver.
+     *
+     * @p inner_solver The Inner solver for the IR solver.
+     */
+    SolverIR(SolverControl &                    solver_control,
+             std::string                        exec_type,
+             std::shared_ptr<gko::LinOpFactory> inner_solver,
+             const AdditionalData &             data = AdditionalData());
 
   protected:
     /**
index b68c414d3b00602bfd4052b5c36455c4fc7925f2..12709a450669d085356e0e3041af5d84c71c253d 100644 (file)
@@ -56,23 +56,27 @@ 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);
-    std::shared_ptr<gko::LinOpFactory> inner_cg = gko::solver::Cg<>::build()
-      .with_criteria(
-                     gko::stop::Iteration::build().with_max_iters(45u).on(exec),
-                     gko::stop::ResidualNormReduction<>::build()
-                     .with_reduction_factor(1e-5)
-                     .on(exec))
-      .on(exec);
-    GinkgoWrappers::SolverCG<>     cg_solver(control, executor);
-    GinkgoWrappers::SolverBICGSTAB<>     bicgstab_solver(control, executor);
-    GinkgoWrappers::SolverCGS<>     cgs_solver(control, executor);
-    GinkgoWrappers::SolverFCG<>     fcg_solver(control, executor);
-    GinkgoWrappers::SolverGMRES<>     gmres_solver(control, executor);
-    GinkgoWrappers::SolverIR<>     ir_solver_cg(control, executor, inner_cg);
-    GinkgoWrappers::SolverCG<>     cg_solver_with_jacobi_precond(control, executor, jacobi);
+    auto                               executor = "reference";
+    std::shared_ptr<gko::Executor>     exec = gko::ReferenceExecutor::create();
+    std::shared_ptr<gko::LinOpFactory> jacobi =
+      gko::preconditioner::Jacobi<>::build().on(exec);
+    std::shared_ptr<gko::LinOpFactory> inner_cg =
+      gko::solver::Cg<>::build()
+        .with_criteria(gko::stop::Iteration::build().with_max_iters(45u).on(
+                         exec),
+                       gko::stop::ResidualNormReduction<>::build()
+                         .with_reduction_factor(1e-5)
+                         .on(exec))
+        .on(exec);
+    GinkgoWrappers::SolverCG<>       cg_solver(control, executor);
+    GinkgoWrappers::SolverBICGSTAB<> bicgstab_solver(control, executor);
+    GinkgoWrappers::SolverCGS<>      cgs_solver(control, executor);
+    GinkgoWrappers::SolverFCG<>      fcg_solver(control, executor);
+    GinkgoWrappers::SolverGMRES<>    gmres_solver(control, executor);
+    GinkgoWrappers::SolverIR<>       ir_solver_cg(control, executor, inner_cg);
+    GinkgoWrappers::SolverCG<>       cg_solver_with_jacobi_precond(control,
+                                                             executor,
+                                                             jacobi);
 
     check_solver_within_range(cg_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.