]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Template TpetraWrappers::SparseMatrix::vmult on the vector type. 16614/head
authorSebastian Kinnewig <kinnewig@ifam.uni-hannover.de>
Fri, 9 Feb 2024 10:40:23 +0000 (11:40 +0100)
committerSebastian Kinnewig <kinnewig@ifam.uni-hannover.de>
Wed, 28 Feb 2024 08:31:23 +0000 (09:31 +0100)
include/deal.II/lac/trilinos_tpetra_sparse_matrix.h
include/deal.II/lac/trilinos_tpetra_sparse_matrix.templates.h
source/lac/trilinos_tpetra_sparse_matrix.cc

index f29536400ed0d31fd8eaac277dee986aa70ad3ad..e466f8985177e3d72e20b821479955d3de39b26b 100644 (file)
 #  include <deal.II/lac/sparsity_pattern.h>
 #  include <deal.II/lac/trilinos_tpetra_sparsity_pattern.h>
 #  include <deal.II/lac/trilinos_tpetra_vector.h>
+#  include <deal.II/lac/vector.h>
 
 // Tpetra includes
 #  include <Tpetra_Core.hpp>
 #  include <Tpetra_CrsMatrix.hpp>
 
+#  include <type_traits>
+
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -162,6 +165,12 @@ namespace LinearAlgebra
       using GraphType =
         Tpetra::CrsGraph<int, dealii::types::signed_global_dof_index, NodeType>;
 
+      /**
+       * Typedef for Tpetra::Vector
+       */
+      using VectorType = Tpetra::
+        Vector<Number, int, dealii::types::signed_global_dof_index, NodeType>;
+
       /**
        * @name Constructors and initialization.
        */
@@ -799,9 +808,9 @@ namespace LinearAlgebra
        * initialized with the same IndexSet that was used for the column indices
        * of the matrix.
        */
+      template <typename InputVectorType>
       void
-      vmult(Vector<Number, MemorySpace>       &dst,
-            const Vector<Number, MemorySpace> &src) const;
+      vmult(InputVectorType &dst, const InputVectorType &src) const;
 
       /*
        * Matrix-vector multiplication: let <i>dst = M<sup>T</sup>*src</i> with
@@ -810,9 +819,9 @@ namespace LinearAlgebra
        *
        * Source and destination must not be the same vector.
        */
+      template <typename InputVectorType>
       void
-      Tvmult(Vector<Number, MemorySpace>       &dst,
-             const Vector<Number, MemorySpace> &src) const;
+      Tvmult(InputVectorType &dst, const InputVectorType &src) const;
 
       /**
        * Adding matrix-vector multiplication. Add <i>M*src</i> on <i>dst</i>
@@ -820,10 +829,9 @@ namespace LinearAlgebra
        *
        * Source and destination must not be the same vector.
        */
+      template <typename InputVectorType>
       void
-      vmult_add(Vector<Number, MemorySpace>       &dst,
-                const Vector<Number, MemorySpace> &src) const;
-
+      vmult_add(InputVectorType &dst, const InputVectorType &src) const;
 
       /**
        * Adding matrix-vector multiplication. Add <i>M<sup>T</sup>*src</i> to
@@ -832,9 +840,9 @@ namespace LinearAlgebra
        *
        * Source and destination must not be the same vector.
        */
+      template <typename InputVectorType>
       void
-      Tvmult_add(Vector<Number, MemorySpace>       &dst,
-                 const Vector<Number, MemorySpace> &src) const;
+      Tvmult_add(InputVectorType &dst, const InputVectorType &src) const;
 
       /**
        * Return the square of the norm of the vector $v$ with respect to the
index e0660c7729b1a99b35c587685e8c21dfe7057f78..4a49d36982c5d3808add7c1433e4a2b0cd5ceed6 100644 (file)
@@ -31,6 +31,114 @@ namespace LinearAlgebra
 
   namespace TpetraWrappers
   {
+    namespace internal
+    {
+      template <typename Number, typename MemorySpace>
+      void
+      apply(const SparseMatrix<Number, MemorySpace> &M,
+            const Vector<Number, MemorySpace>       &src,
+            Vector<Number, MemorySpace>             &dst,
+            Teuchos::ETransp                         mode = Teuchos::NO_TRANS,
+            Number alpha = Teuchos::ScalarTraits<Number>::one(),
+            Number beta  = Teuchos::ScalarTraits<Number>::zero())
+      {
+        Assert(&src != &dst,
+               SparseMatrix<double>::ExcSourceEqualsDestination());
+        Assert(M.trilinos_matrix().isFillComplete(),
+               SparseMatrix<double>::ExcMatrixNotCompressed());
+
+        if (mode == Teuchos::NO_TRANS)
+          {
+            Assert(src.trilinos_vector().getMap()->isSameAs(
+                     *M.trilinos_matrix().getDomainMap()),
+                   SparseMatrix<double>::ExcColMapMissmatch());
+            Assert(dst.trilinos_vector().getMap()->isSameAs(
+                     *M.trilinos_matrix().getRangeMap()),
+                   SparseMatrix<double>::ExcDomainMapMissmatch());
+          }
+        else
+          {
+            Assert(dst.trilinos_vector().getMap()->isSameAs(
+                     *M.trilinos_matrix().getDomainMap()),
+                   SparseMatrix<double>::ExcColMapMissmatch());
+            Assert(src.trilinos_vector().getMap()->isSameAs(
+                     *M.trilinos_matrix().getRangeMap()),
+                   SparseMatrix<double>::ExcDomainMapMissmatch());
+          }
+
+        M.trilinos_matrix().apply(
+          src.trilinos_vector(), dst.trilinos_vector(), mode, alpha, beta);
+      }
+
+
+
+      template <typename Number, typename MemorySpace>
+      void
+      apply(const SparseMatrix<Number, MemorySpace> &M,
+            const dealii::Vector<Number>            &src,
+            dealii::Vector<Number>                  &dst,
+            Teuchos::ETransp                         mode = Teuchos::NO_TRANS,
+            Number alpha = Teuchos::ScalarTraits<Number>::one(),
+            Number beta  = Teuchos::ScalarTraits<Number>::zero())
+      {
+        Assert(&src != &dst,
+               SparseMatrix<double>::ExcSourceEqualsDestination());
+        Assert(M.trilinos_matrix().isFillComplete(),
+               SparseMatrix<double>::ExcMatrixNotCompressed());
+
+        // get the size of the input vectors:
+        const size_type dst_local_size = dst.end() - dst.begin();
+        const size_type src_local_size = src.end() - src.begin();
+
+        // For the dst vector:
+        Kokkos::View<Number **, Kokkos::LayoutLeft, Kokkos::HostSpace>
+          kokkos_view_dst(dst.begin(), dst_local_size, 1);
+
+        // get a Kokkos::DualView
+        auto mirror_view_dst = Kokkos::create_mirror_view_and_copy(
+          typename MemorySpace::kokkos_space{}, kokkos_view_dst);
+        typename SparseMatrix<Number, MemorySpace>::VectorType::dual_view_type
+          kokkos_dual_view_dst(mirror_view_dst, kokkos_view_dst);
+
+        // create the Tpetra::Vector
+        typename SparseMatrix<Number, MemorySpace>::VectorType tpetra_dst(
+          M.trilinos_matrix().getRangeMap(), kokkos_dual_view_dst);
+
+        // For the src vector:
+        // create a Kokkos::View from the src vector
+        Kokkos::View<Number **, Kokkos::LayoutLeft, Kokkos::HostSpace>
+          kokkos_view_src(const_cast<Number *>(src.begin()), src_local_size, 1);
+
+        // get a Kokkos::DualView
+        auto mirror_view_src = Kokkos::create_mirror_view_and_copy(
+          typename MemorySpace::kokkos_space{}, kokkos_view_src);
+        typename SparseMatrix<Number, MemorySpace>::VectorType::dual_view_type
+          kokkos_dual_view_src(mirror_view_src, kokkos_view_src);
+
+        // create the Tpetra::Vector
+        typename SparseMatrix<Number, MemorySpace>::VectorType tpetra_src(
+          M.trilinos_matrix().getDomainMap(), kokkos_dual_view_src);
+
+        M.trilinos_matrix().apply(tpetra_src, tpetra_dst, mode, alpha, beta);
+      }
+
+
+
+      template <typename Number, typename MemorySpace, typename VectorType>
+      void
+      apply(SparseMatrix<Number, MemorySpace> &,
+            const VectorType &,
+            VectorType &,
+            Teuchos::ETransp,
+            Number,
+            Number)
+      {
+        DEAL_II_NOT_IMPLEMENTED();
+      }
+    } // namespace internal
+
+
+
     // reinit_matrix():
     namespace
     {
@@ -1170,81 +1278,58 @@ namespace LinearAlgebra
 
 
     // Multiplications
-
     template <typename Number, typename MemorySpace>
+    template <typename InputVectorType>
     void
-    SparseMatrix<Number, MemorySpace>::vmult(
-      Vector<Number, MemorySpace>       &dst,
-      const Vector<Number, MemorySpace> &src) const
+    SparseMatrix<Number, MemorySpace>::vmult(InputVectorType       &dst,
+                                             const InputVectorType &src) const
     {
-      Assert(&src != &dst, ExcSourceEqualsDestination());
-      Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
-      Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
-             ExcColMapMissmatch());
-      Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
-             ExcDomainMapMissmatch());
-      matrix->apply(src.trilinos_vector(), dst.trilinos_vector());
+      internal::apply(*this, src, dst);
     }
 
 
 
     template <typename Number, typename MemorySpace>
+    template <typename InputVectorType>
     void
-    SparseMatrix<Number, MemorySpace>::Tvmult(
-      Vector<Number, MemorySpace>       &dst,
-      const Vector<Number, MemorySpace> &src) const
+    SparseMatrix<Number, MemorySpace>::Tvmult(InputVectorType       &dst,
+                                              const InputVectorType &src) const
     {
-      Assert(&src != &dst, ExcSourceEqualsDestination());
-      Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
-      Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
-             ExcColMapMissmatch());
-      Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
-             ExcDomainMapMissmatch());
-      matrix->apply(src.trilinos_vector(),
-                    dst.trilinos_vector(),
-                    Teuchos::TRANS);
+      internal::apply(*this, src, dst, Teuchos::TRANS);
     }
 
 
 
     template <typename Number, typename MemorySpace>
+    template <typename InputVectorType>
     void
     SparseMatrix<Number, MemorySpace>::vmult_add(
-      Vector<Number, MemorySpace>       &dst,
-      const Vector<Number, MemorySpace> &src) const
+      InputVectorType       &dst,
+      const InputVectorType &src) const
     {
-      Assert(&src != &dst, ExcSourceEqualsDestination());
-      Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
-      Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
-             ExcColMapMissmatch());
-      Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
-             ExcDomainMapMissmatch());
-      matrix->apply(src.trilinos_vector(),
-                    dst.trilinos_vector(),
-                    Teuchos::NO_TRANS,
-                    Teuchos::ScalarTraits<Number>::one(),
-                    Teuchos::ScalarTraits<Number>::one());
+      internal::apply(*this,
+                      src,
+                      dst,
+                      Teuchos::NO_TRANS,
+                      Teuchos::ScalarTraits<Number>::one(),
+                      Teuchos::ScalarTraits<Number>::one());
     }
 
 
 
     template <typename Number, typename MemorySpace>
+    template <typename InputVectorType>
     void
     SparseMatrix<Number, MemorySpace>::Tvmult_add(
-      Vector<Number, MemorySpace>       &dst,
-      const Vector<Number, MemorySpace> &src) const
+      InputVectorType       &dst,
+      const InputVectorType &src) const
     {
-      Assert(&src != &dst, ExcSourceEqualsDestination());
-      Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
-      Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
-             ExcColMapMissmatch());
-      Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
-             ExcDomainMapMissmatch());
-      matrix->apply(src.trilinos_vector(),
-                    dst.trilinos_vector(),
-                    Teuchos::TRANS,
-                    Teuchos::ScalarTraits<Number>::one(),
-                    Teuchos::ScalarTraits<Number>::one());
+      internal::apply(*this,
+                      src,
+                      dst,
+                      Teuchos::TRANS,
+                      Teuchos::ScalarTraits<Number>::one(),
+                      Teuchos::ScalarTraits<Number>::one());
     }
 
 
index 76c5335856a72f26d2d48a3d79bd5d93f7d7df73..3cc43d6226cd7e7973a8a724a526d12f01611e79 100644 (file)
@@ -46,6 +46,37 @@ namespace LinearAlgebra
     template void
     SparseMatrix<double>::reinit(const dealii::DynamicSparsityPattern &);
 
+    template void
+    SparseMatrix<double>::vmult(Vector<double>       &dst,
+                                const Vector<double> &src) const;
+
+    template void
+    SparseMatrix<double>::Tvmult(Vector<double>       &dst,
+                                 const Vector<double> &src) const;
+
+    template void
+    SparseMatrix<double>::vmult_add(Vector<double>       &dst,
+                                    const Vector<double> &src) const;
+
+    template void
+    SparseMatrix<double>::Tvmult_add(Vector<double>       &dst,
+                                     const Vector<double> &src) const;
+
+    template void
+    SparseMatrix<double>::vmult(::dealii::Vector<double>       &dst,
+                                const ::dealii::Vector<double> &src) const;
+
+    template void
+    SparseMatrix<double>::Tvmult(::dealii::Vector<double>       &dst,
+                                 const ::dealii::Vector<double> &src) const;
+
+    template void
+    SparseMatrix<double>::vmult_add(::dealii::Vector<double>       &dst,
+                                    const ::dealii::Vector<double> &src) const;
+
+    template void
+    SparseMatrix<double>::Tvmult_add(::dealii::Vector<double>       &dst,
+                                     const ::dealii::Vector<double> &src) const;
   } // namespace TpetraWrappers
 } // namespace LinearAlgebra
 #  endif // DOXYGEN

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.