]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Annotate MatrixType with linear operator concepts. 17010/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 10 May 2024 09:30:04 +0000 (15:00 +0530)
committerWolfgang Bangerth <bangerth@colostate.edu>
Sun, 12 May 2024 16:17:27 +0000 (10:17 -0600)
include/deal.II/lac/solver_bicgstab.h
include/deal.II/lac/solver_cg.h
include/deal.II/lac/solver_fire.h
include/deal.II/lac/solver_gmres.h
include/deal.II/lac/solver_idr.h
include/deal.II/lac/solver_minres.h
include/deal.II/lac/solver_qmrs.h
include/deal.II/lac/solver_relaxation.h
include/deal.II/lac/solver_richardson.h
include/deal.II/lac/solver_selector.h

index 79b3d67974293826fa0c70b2d98fe85211736686..2695c61f311994a2b2d42aa0b63edf4118dc9417 100644 (file)
@@ -139,11 +139,13 @@ public:
    * Solve primal problem only.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
 protected:
   /**
@@ -405,6 +407,9 @@ typename SolverBicgstab<VectorType>::IterationResult
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverBicgstab<VectorType>::solve(const MatrixType         &A,
                                        VectorType               &x,
                                        const VectorType         &b,
index 498df053c42b2ae138c7260b0030d061e6451a4c..dca6a25b68a7dbdf1859f3a5aef44111a0f1d8ec 100644 (file)
@@ -213,11 +213,13 @@ public:
    * Solve the linear system $Ax=b$ for x.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
   /**
    * Connect a slot to retrieve the CG coefficients. The slot will be called
@@ -1265,6 +1267,9 @@ namespace internal
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverCG<VectorType>::solve(const MatrixType         &A,
                                  VectorType               &x,
                                  const VectorType         &b,
index 4efed16e8592f8b774084e01bddbd37c156821ea..279905ae9d407a2fd9868c995b48a5f5fc1acf1b 100644 (file)
@@ -157,11 +157,13 @@ public:
    * = \frac{1}{2} \mathbf x^{T} \mathbf A \mathbf x - \mathbf x^{T} \mathbf b$.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
 protected:
   /**
@@ -353,6 +355,9 @@ void SolverFIRE<VectorType>::solve(
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverFIRE<VectorType>::solve(const MatrixType         &A,
                                    VectorType               &x,
                                    const VectorType         &b,
index 71466995b107a921fc057bee3aa1b0c58ddbec09..82f306f40489988c9b60184b80a2430ecd6c4504 100644 (file)
@@ -457,11 +457,13 @@ public:
    * Solve the linear system $Ax=b$ for x.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
   /**
    * Connect a slot to retrieve the estimated condition number. Called on each
@@ -691,11 +693,13 @@ public:
    * Solve the linear system $Ax=b$ for x.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
 private:
   /**
@@ -1777,6 +1781,9 @@ inline void SolverGMRES<VectorType>::compute_eigs_and_cond(
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverGMRES<VectorType>::solve(const MatrixType         &A,
                                     VectorType               &x,
                                     const VectorType         &b,
@@ -2164,6 +2171,9 @@ SolverFGMRES<VectorType>::SolverFGMRES(SolverControl        &cn,
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverFGMRES<VectorType>::solve(const MatrixType         &A,
                                      VectorType               &x,
                                      const VectorType         &b,
index bb9e059b4a2bd4b14abbfd2b6019125ce2956bd1..693868ca2b35e19bf0a20e7e873b1c7145249909 100644 (file)
@@ -158,11 +158,13 @@ public:
    * Solve the linear system <code>Ax=b</code> for x.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
 protected:
   /**
@@ -312,6 +314,9 @@ void SolverIDR<VectorType>::print_vectors(const unsigned int,
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverIDR<VectorType>::solve(const MatrixType         &A,
                                   VectorType               &x,
                                   const VectorType         &b,
index 50782660878e6120fc8169efa3ec53f52ab25eaa..6c2b62eb56cfbb2d47b6ee1705d55757ec3bcf47 100644 (file)
@@ -101,11 +101,13 @@ public:
    * Solve the linear system $Ax=b$ for x.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
   /**
    * @addtogroup Exceptions
@@ -192,6 +194,9 @@ void SolverMinRes<VectorType>::print_vectors(const unsigned int,
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverMinRes<VectorType>::solve(const MatrixType         &A,
                                      VectorType               &x,
                                      const VectorType         &b,
index 75e29c73b695ebb0facec0f6f153422b4c73dce4..1eab3ca3e65518cf3236941364d95cd934f004ea 100644 (file)
@@ -176,11 +176,13 @@ public:
    * Solve the linear system $Ax=b$ for x.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
   /**
    * Interface for derived class. This function gets the current iteration
@@ -288,6 +290,9 @@ void SolverQMRS<VectorType>::print_vectors(const unsigned int,
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverQMRS<VectorType>::solve(const MatrixType         &A,
                                    VectorType               &x,
                                    const VectorType         &b,
index cf7df9167e314d6ba9d9eb0eb9887a306843f99e..f240a798a0b1f2a48c8be641d1a5329f5ea8c1b8 100644 (file)
@@ -77,11 +77,15 @@ public:
    * residual.
    */
   template <typename MatrixType, typename RelaxationType>
-  void
-  solve(const MatrixType     &A,
-        VectorType           &x,
-        const VectorType     &b,
-        const RelaxationType &R);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     requires(const RelaxationType &R, VectorType &a, VectorType &b) {
+       R.step(a, b);
+     }))
+  void solve(const MatrixType     &A,
+             VectorType           &x,
+             const VectorType     &b,
+             const RelaxationType &R);
 };
 
 //----------------------------------------------------------------------//
@@ -98,6 +102,11 @@ SolverRelaxation<VectorType>::SolverRelaxation(SolverControl &cn,
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename RelaxationType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   requires(const RelaxationType &R, VectorType &a, VectorType &b) {
+     R.step(a, b);
+   }))
 void SolverRelaxation<VectorType>::solve(const MatrixType     &A,
                                          VectorType           &x,
                                          const VectorType     &b,
index 54bab807bff74b4fc0503e13e8c3184ee2146bc6..db10303de82cd882cc4e5c18a85d653d9753cbe5 100644 (file)
@@ -110,21 +110,25 @@ public:
    * Solve the linear system $Ax=b$ for x.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &preconditioner);
 
   /**
    * Solve $A^Tx=b$ for $x$.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  Tsolve(const MatrixType         &A,
-         VectorType               &x,
-         const VectorType         &b,
-         const PreconditionerType &preconditioner);
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_transpose_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_transpose_linear_operator_on<PreconditionerType, VectorType>))
+  void Tsolve(const MatrixType         &A,
+              VectorType               &x,
+              const VectorType         &b,
+              const PreconditionerType &preconditioner);
 
   /**
    * Set the damping-coefficient. Default is 1., i.e. no damping.
@@ -198,6 +202,9 @@ SolverRichardson<VectorType>::SolverRichardson(SolverControl        &cn,
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverRichardson<VectorType>::solve(
   const MatrixType         &A,
   VectorType               &x,
@@ -256,6 +263,9 @@ void SolverRichardson<VectorType>::solve(
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_transpose_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_transpose_linear_operator_on<PreconditionerType, VectorType>))
 void SolverRichardson<VectorType>::Tsolve(
   const MatrixType         &A,
   VectorType               &x,
index dbf036333eab5656bb0ee0a2f62d4e7ff557db25..d7bba231d81a2259785997076c0829a84e38610b 100644 (file)
@@ -118,11 +118,13 @@ public:
    * SolverName was specified in the constructor.
    */
   template <typename MatrixType, typename PreconditionerType>
-  void
-  solve(const MatrixType         &A,
-        VectorType               &x,
-        const VectorType         &b,
-        const PreconditionerType &precond) const;
+  DEAL_II_CXX20_REQUIRES(
+    (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+     concepts::is_linear_operator_on<PreconditionerType, VectorType>))
+  void solve(const MatrixType         &A,
+             VectorType               &x,
+             const VectorType         &b,
+             const PreconditionerType &precond) const;
 
   /**
    * Select a new solver. Note that all solver names used in this class are
@@ -269,6 +271,9 @@ void SolverSelector<VectorType>::select(const std::string &name)
 template <typename VectorType>
 DEAL_II_CXX20_REQUIRES(concepts::is_vector_space_vector<VectorType>)
 template <typename MatrixType, typename PreconditionerType>
+DEAL_II_CXX20_REQUIRES(
+  (concepts::is_linear_operator_on<MatrixType, VectorType> &&
+   concepts::is_linear_operator_on<PreconditionerType, VectorType>))
 void SolverSelector<VectorType>::solve(const MatrixType         &A,
                                        VectorType               &x,
                                        const VectorType         &b,

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.