]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Permit TrilinosPayload to use itself as an exemplar operator
authorJean-Paul Pelteret <jppelteret@gmail.com>
Sun, 12 Sep 2021 20:37:48 +0000 (22:37 +0200)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Tue, 14 Sep 2021 19:32:34 +0000 (21:32 +0200)
include/deal.II/lac/trilinos_sparse_matrix.h
source/lac/trilinos_sparse_matrix.cc

index 33f40bcdeefa4887d0b2bc9b3ee4e1c3161e83a3..d9c33137f31f4022eac762ae15a169ce8eaaafea 100644 (file)
@@ -2080,6 +2080,12 @@ namespace TrilinosWrappers
         TrilinosPayload(const TrilinosWrappers::SparseMatrix &matrix_exemplar,
                         const TrilinosWrappers::SparseMatrix &matrix);
 
+        /**
+         * Constructor for a sparse matrix based on an exemplary payload
+         */
+        TrilinosPayload(const TrilinosPayload &               payload_exemplar,
+                        const TrilinosWrappers::SparseMatrix &matrix);
+
         /**
          * Constructor for a preconditioner based on an exemplary matrix
          */
@@ -2094,6 +2100,13 @@ namespace TrilinosWrappers
           const TrilinosWrappers::PreconditionBase &preconditioner_exemplar,
           const TrilinosWrappers::PreconditionBase &preconditioner);
 
+        /**
+         * Constructor for a preconditioner based on an exemplary payload
+         */
+        TrilinosPayload(
+          const TrilinosPayload &                   payload_exemplar,
+          const TrilinosWrappers::PreconditionBase &preconditioner);
+
         /**
          * Default copy constructor
          */
index 35e67ac9e78cfe419d6061f2917d4cdc5ea53b6f..ad7864f12a2b3cb3acc487fe868c8839fb382003 100644 (file)
@@ -2428,6 +2428,81 @@ namespace TrilinosWrappers
 
 
 
+      TrilinosPayload::TrilinosPayload(
+        const TrilinosPayload &               payload_exemplar,
+        const TrilinosWrappers::SparseMatrix &matrix)
+        : use_transpose(payload_exemplar.UseTranspose())
+        , communicator(payload_exemplar.get_mpi_communicator())
+        , domain_map(
+            payload_exemplar.locally_owned_domain_indices().make_trilinos_map(
+              communicator.Comm()))
+        , range_map(
+            payload_exemplar.locally_owned_range_indices().make_trilinos_map(
+              communicator.Comm()))
+      {
+        vmult = [&payload_exemplar, &matrix](Range &       tril_dst,
+                                             const Domain &tril_src) {
+          // Duplicated from TrilinosWrappers::SparseMatrix::vmult
+          Assert(&tril_src != &tril_dst,
+                 TrilinosWrappers::SparseMatrix::ExcSourceEqualsDestination());
+          Assert(matrix.trilinos_matrix().Filled(),
+                 TrilinosWrappers::SparseMatrix::ExcMatrixNotCompressed());
+          internal::check_vector_map_equality(payload_exemplar,
+                                              tril_src,
+                                              tril_dst,
+                                              payload_exemplar.UseTranspose());
+          internal::check_vector_map_equality(
+            matrix.trilinos_matrix(),
+            tril_src,
+            tril_dst,
+            matrix.trilinos_matrix().UseTranspose());
+
+          const int ierr = matrix.trilinos_matrix().Apply(tril_src, tril_dst);
+          AssertThrow(ierr == 0, ExcTrilinosError(ierr));
+        };
+
+        Tvmult = [&payload_exemplar, &matrix](Domain &     tril_dst,
+                                              const Range &tril_src) {
+          // Duplicated from TrilinosWrappers::SparseMatrix::Tvmult
+          Assert(&tril_src != &tril_dst,
+                 TrilinosWrappers::SparseMatrix::ExcSourceEqualsDestination());
+          Assert(matrix.trilinos_matrix().Filled(),
+                 TrilinosWrappers::SparseMatrix::ExcMatrixNotCompressed());
+          internal::check_vector_map_equality(payload_exemplar,
+                                              tril_src,
+                                              tril_dst,
+                                              !payload_exemplar.UseTranspose());
+          internal::check_vector_map_equality(
+            matrix.trilinos_matrix(),
+            tril_src,
+            tril_dst,
+            !matrix.trilinos_matrix().UseTranspose());
+
+          Epetra_CrsMatrix &tril_mtrx_non_const =
+            const_cast<Epetra_CrsMatrix &>(matrix.trilinos_matrix());
+          tril_mtrx_non_const.SetUseTranspose(
+            !matrix.trilinos_matrix().UseTranspose());
+          const int ierr = matrix.trilinos_matrix().Apply(tril_src, tril_dst);
+          AssertThrow(ierr == 0, ExcTrilinosError(ierr));
+          tril_mtrx_non_const.SetUseTranspose(
+            !matrix.trilinos_matrix().UseTranspose());
+        };
+
+        inv_vmult = [](Domain &, const Range &) {
+          Assert(false,
+                 ExcMessage("Uninitialized TrilinosPayload::inv_vmult called "
+                            "(Matrix constructor with matrix exemplar)"));
+        };
+
+        inv_Tvmult = [](Range &, const Domain &) {
+          Assert(false,
+                 ExcMessage("Uninitialized TrilinosPayload::inv_Tvmult called "
+                            "(Matrix constructor with matrix exemplar)"));
+        };
+      }
+
+
+
       TrilinosPayload::TrilinosPayload(
         const TrilinosWrappers::SparseMatrix &    matrix_exemplar,
         const TrilinosWrappers::PreconditionBase &preconditioner)
@@ -2651,6 +2726,114 @@ namespace TrilinosWrappers
 
 
 
+      TrilinosPayload::TrilinosPayload(
+        const TrilinosPayload &                   payload_exemplar,
+        const TrilinosWrappers::PreconditionBase &preconditioner)
+        : use_transpose(payload_exemplar.UseTranspose())
+        , communicator(payload_exemplar.get_mpi_communicator())
+        , domain_map(
+            payload_exemplar.locally_owned_domain_indices().make_trilinos_map(
+              communicator.Comm()))
+        , range_map(
+            payload_exemplar.locally_owned_range_indices().make_trilinos_map(
+              communicator.Comm()))
+      {
+        vmult = [&payload_exemplar, &preconditioner](Range &       tril_dst,
+                                                     const Domain &tril_src) {
+          // Duplicated from TrilinosWrappers::PreconditionBase::vmult
+          // as well as from TrilinosWrappers::SparseMatrix::Tvmult
+          Assert(&tril_src != &tril_dst,
+                 TrilinosWrappers::SparseMatrix::ExcSourceEqualsDestination());
+          internal::check_vector_map_equality(payload_exemplar,
+                                              tril_src,
+                                              tril_dst,
+                                              payload_exemplar.UseTranspose());
+          internal::check_vector_map_equality(
+            preconditioner.trilinos_operator(),
+            tril_src,
+            tril_dst,
+            preconditioner.trilinos_operator().UseTranspose());
+
+          const int ierr =
+            preconditioner.trilinos_operator().ApplyInverse(tril_src, tril_dst);
+          AssertThrow(ierr == 0, ExcTrilinosError(ierr));
+        };
+
+        Tvmult = [&payload_exemplar, &preconditioner](Domain &     tril_dst,
+                                                      const Range &tril_src) {
+          // Duplicated from TrilinosWrappers::PreconditionBase::vmult
+          // as well as from TrilinosWrappers::SparseMatrix::Tvmult
+          Assert(&tril_src != &tril_dst,
+                 TrilinosWrappers::SparseMatrix::ExcSourceEqualsDestination());
+          internal::check_vector_map_equality(payload_exemplar,
+                                              tril_src,
+                                              tril_dst,
+                                              !payload_exemplar.UseTranspose());
+          internal::check_vector_map_equality(
+            preconditioner.trilinos_operator(),
+            tril_src,
+            tril_dst,
+            !preconditioner.trilinos_operator().UseTranspose());
+
+          preconditioner.trilinos_operator().SetUseTranspose(
+            !preconditioner.trilinos_operator().UseTranspose());
+          const int ierr =
+            preconditioner.trilinos_operator().ApplyInverse(tril_src, tril_dst);
+          AssertThrow(ierr == 0, ExcTrilinosError(ierr));
+          preconditioner.trilinos_operator().SetUseTranspose(
+            !preconditioner.trilinos_operator().UseTranspose());
+        };
+
+        inv_vmult = [&payload_exemplar,
+                     &preconditioner](Domain &tril_dst, const Range &tril_src) {
+          // Duplicated from TrilinosWrappers::PreconditionBase::vmult
+          // as well as from TrilinosWrappers::SparseMatrix::Tvmult
+          Assert(&tril_src != &tril_dst,
+                 TrilinosWrappers::SparseMatrix::ExcSourceEqualsDestination());
+          internal::check_vector_map_equality(payload_exemplar,
+                                              tril_src,
+                                              tril_dst,
+                                              !payload_exemplar.UseTranspose());
+          internal::check_vector_map_equality(
+            preconditioner.trilinos_operator(),
+            tril_src,
+            tril_dst,
+            !preconditioner.trilinos_operator().UseTranspose());
+
+          const int ierr =
+            preconditioner.trilinos_operator().ApplyInverse(tril_src, tril_dst);
+          AssertThrow(ierr == 0, ExcTrilinosError(ierr));
+        };
+
+        inv_Tvmult = [&payload_exemplar,
+                      &preconditioner](Range &       tril_dst,
+                                       const Domain &tril_src) {
+          // Duplicated from TrilinosWrappers::PreconditionBase::vmult
+          // as well as from TrilinosWrappers::SparseMatrix::Tvmult
+          Assert(&tril_src != &tril_dst,
+                 TrilinosWrappers::SparseMatrix::ExcSourceEqualsDestination());
+          internal::check_vector_map_equality(payload_exemplar,
+                                              tril_src,
+                                              tril_dst,
+                                              payload_exemplar.UseTranspose());
+          internal::check_vector_map_equality(
+            preconditioner.trilinos_operator(),
+            tril_src,
+            tril_dst,
+            preconditioner.trilinos_operator().UseTranspose());
+
+          preconditioner.trilinos_operator().SetUseTranspose(
+            !preconditioner.trilinos_operator().UseTranspose());
+          const int ierr =
+            preconditioner.trilinos_operator().ApplyInverse(tril_src, tril_dst);
+          AssertThrow(ierr == 0, ExcTrilinosError(ierr));
+          preconditioner.trilinos_operator().SetUseTranspose(
+            !preconditioner.trilinos_operator().UseTranspose());
+        };
+      }
+
+
+
       TrilinosPayload::TrilinosPayload(const TrilinosPayload &payload)
         : vmult(payload.vmult)
         , Tvmult(payload.Tvmult)

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.