]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Only enable a function if a template type is right.
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 15 Feb 2024 22:47:34 +0000 (15:47 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 15 Feb 2024 22:49:17 +0000 (15:49 -0700)
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 48de7bd615105dfe438e2f0a405ab5e1a9a730f2..60b6bcbf27f210ceeee614bbe6e17b7c92b494de 100644 (file)
@@ -319,7 +319,8 @@ namespace LinearAlgebra
        * processors in order to avoid a dead lock.
        */
       template <typename SparsityPatternType>
-      void
+      std::enable_if_t<
+        !std::is_same_v<SparsityPatternType, dealii::SparseMatrix<double>>>
       reinit(const IndexSet            &parallel_partitioning,
              const SparsityPatternType &sparsity_pattern,
              const MPI_Comm             communicator  = MPI_COMM_WORLD,
@@ -338,7 +339,8 @@ namespace LinearAlgebra
        * processors in order to avoid a dead lock.
        */
       template <typename SparsityPatternType>
-      void
+      std::enable_if_t<
+        !std::is_same_v<SparsityPatternType, dealii::SparseMatrix<double>>>
       reinit(const IndexSet            &row_parallel_partitioning,
              const IndexSet            &col_parallel_partitioning,
              const SparsityPatternType &sparsity_pattern,
index de1d09db5053e49a6ec6b051e7af6073e5bdc4c2..b703deaf66f7e7444bd0e46beb690b83257968f4 100644 (file)
@@ -430,7 +430,8 @@ namespace LinearAlgebra
 
     template <typename Number, typename MemorySpace>
     template <typename SparsityPatternType>
-    inline void
+    inline std::enable_if_t<
+      !std::is_same_v<SparsityPatternType, dealii::SparseMatrix<double>>>
     SparseMatrix<Number, MemorySpace>::reinit(
       const IndexSet            &parallel_partitioning,
       const SparsityPatternType &sparsity_pattern,
@@ -448,10 +449,10 @@ namespace LinearAlgebra
 
     template <typename Number, typename MemorySpace>
     template <typename SparsityPatternType>
-    void
+    std::enable_if_t<
+      !std::is_same_v<SparsityPatternType, dealii::SparseMatrix<double>>>
     SparseMatrix<Number, MemorySpace>::reinit(
-      const IndexSet &row_parallel_partitioning,
-
+      const IndexSet            &row_parallel_partitioning,
       const IndexSet            &col_parallel_partitioning,
       const SparsityPatternType &sparsity_pattern,
       const MPI_Comm             communicator,
index 63f970f3edd38f9b8685d18ae4ebd14b08159286..36945b1843cd35560f83dfadc3d21068bcb8c0b2 100644 (file)
@@ -29,6 +29,13 @@ namespace LinearAlgebra
   {
     template class SparseMatrix<double>;
 
+    template void
+    SparseMatrix<double>::reinit(
+      const IndexSet                       &parallel_partitioning,
+      const dealii::DynamicSparsityPattern &sparsity_pattern,
+      const MPI_Comm                        communicator,
+      const bool                            exchange_data);
+
     template void
     SparseMatrix<double>::reinit(
       const IndexSet                       &row_parallel_partitioning,

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.