]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add 'requires' clauses to PETSc solvers.
authorWolfgang Bangerth <bangerth@colostate.edu>
Sat, 22 Apr 2023 23:23:26 +0000 (17:23 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Sat, 22 Apr 2023 23:23:26 +0000 (17:23 -0600)
include/deal.II/base/template_constraints.h
include/deal.II/lac/petsc_snes.h
include/deal.II/lac/petsc_ts.h

index da9ba1613106620cb88dc9de420e265ff1cc107a..0f735f42f6b2cad6cf04d3ab583b49d1c821e489 100644 (file)
@@ -656,6 +656,7 @@ namespace LinearAlgebra
 #ifdef DEAL_II_WITH_PETSC
 namespace PETScWrappers
 {
+  class VectorBase;
   class Vector;
   class BlockVector;
 
@@ -788,6 +789,40 @@ namespace concepts
       dealii::LinearAlgebra::TpetraWrappers::Vector<Number>> = true;
 #    endif
 #  endif
+
+
+    /**
+     * A template variable that returns whether the template argument is
+     * a valid deal.II vector type that is internally built on PETSc
+     * functionality. Its general definition is `false`, with
+     * specializations dealing with actual vector types for which the
+     * predicate is `true`.
+     */
+    template <typename T>
+    inline constexpr bool is_dealii_petsc_vector_type = false;
+
+#  ifdef DEAL_II_WITH_PETSC
+    template <>
+    inline constexpr bool
+      is_dealii_petsc_vector_type<dealii::PETScWrappers::VectorBase> = true;
+
+    template <>
+    inline constexpr bool
+      is_dealii_petsc_vector_type<dealii::PETScWrappers::Vector> = true;
+
+    template <>
+    inline constexpr bool
+      is_dealii_petsc_vector_type<dealii::PETScWrappers::BlockVector> = true;
+
+    template <>
+    inline constexpr bool
+      is_dealii_petsc_vector_type<dealii::PETScWrappers::MPI::Vector> = true;
+
+    template <>
+    inline constexpr bool
+      is_dealii_petsc_vector_type<dealii::PETScWrappers::MPI::BlockVector> =
+        true;
+#  endif
   } // namespace internal
 
 
@@ -820,6 +855,18 @@ namespace concepts
                                            (std::is_const_v<VectorType> ==
                                             false);
 
+  /**
+   * A concept that tests whether a given template argument is a deal.II
+   * vector type that internally builds on PETSc functionality. This
+   * concept is used to constrain some classes that implement advanced
+   * functionality based on PETSc and that requires that the vector
+   * it works on are PETSc vectors. This includes, for example, the
+   * time stepping and nonlinear solver classes in namespace PETScWrappers.
+   */
+  template <typename VectorType>
+  concept is_dealii_petsc_vector_type =
+    internal::is_dealii_petsc_vector_type<VectorType>;
+
 #endif
 } // namespace concepts
 
index ebae3fa4dabb19d5864ca97c2cc118edcd6cbd3c..11c6db51e5a05772acd52ed80e162a2dd81d5e95 100644 (file)
@@ -224,6 +224,7 @@ namespace PETScWrappers
   template <typename VectorType  = PETScWrappers::VectorBase,
             typename PMatrixType = PETScWrappers::MatrixBase,
             typename AMatrixType = PMatrixType>
+  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_petsc_vector_type<VectorType>)
   class NonlinearSolver
   {
   public:
index 3dc726407f56686bd52af5d0f03134e529bae55e..64e9cb235f77338d8e035426ac4cab9e6959d16d 100644 (file)
@@ -289,6 +289,7 @@ namespace PETScWrappers
   template <typename VectorType  = PETScWrappers::VectorBase,
             typename PMatrixType = PETScWrappers::MatrixBase,
             typename AMatrixType = PMatrixType>
+  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_petsc_vector_type<VectorType>)
   class TimeStepper
   {
   public:

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.