]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce make_rcp(...) 16270/head
authorSebastian Kinnewig <kinnewig@ifam.uni-hannover.de>
Wed, 15 Nov 2023 12:08:20 +0000 (13:08 +0100)
committerSebastian Kinnewig <sebastian@kinnewig.org>
Fri, 17 Nov 2023 16:29:01 +0000 (17:29 +0100)
doc/news/changes/minor/20231115SebastianKinnewig [new file with mode: 0644]
include/deal.II/base/index_set.h
include/deal.II/base/trilinos_utilities.h
include/deal.II/lac/trilinos_tpetra_vector.templates.h
source/base/index_set.cc
source/base/trilinos_utilities.cc

diff --git a/doc/news/changes/minor/20231115SebastianKinnewig b/doc/news/changes/minor/20231115SebastianKinnewig
new file mode 100644 (file)
index 0000000..881ab68
--- /dev/null
@@ -0,0 +1,4 @@
+New: Introduce a new function, Utilities::Trilinos::internal::make_rcp(),
+to create Teuchos::RCP objects, avoiding the usage of 'operator new'.
+<br>
+(Sebastian Kinnewig, 2023/11/15)
index 28868eac534d505578aa143b38dbf5cb84fe0496..8cc267ebbae0413132bdf2a40fa4371e223081ed 100644 (file)
@@ -21,6 +21,7 @@
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/mpi_stub.h>
 #include <deal.II/base/mutex.h>
+#include <deal.II/base/trilinos_utilities.h>
 
 #include <boost/container/small_vector.hpp>
 
index 761a1683d0e5c15e7808c55845130f8cf9e5cc2f..bb3e1e5958adebd4607c0f7da4fbe23bc0d24081 100644 (file)
 #  endif
 #endif
 
+#ifdef DEAL_II_TRILINOS_WITH_TPETRA
+#  include <Teuchos_RCPDecl.hpp>
+#endif // DEAL_II_TRILINOS_WITH_TPETRA
+
 DEAL_II_NAMESPACE_OPEN
 
 /**
@@ -180,6 +184,35 @@ namespace Utilities
     duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm);
   } // namespace Trilinos
 #endif
+
+
+#ifdef DEAL_II_TRILINOS_WITH_TPETRA
+  namespace Trilinos
+  {
+    namespace internal
+    {
+      /**
+       * Creates and returns a
+       * <a
+       * href="https://docs.trilinos.org/dev/packages/teuchos/doc/html/classTeuchos_1_1RCP.html">Teuchos::RCP</a>
+       * object for type T.
+       *
+       * @note In Trilinos 14.0.0, the function
+       * <a
+       * href="https://docs.trilinos.org/dev/packages/teuchos/doc/html/namespaceTeuchos.html#a280c0ab8c9ee8d0481114d4edf5a3393">Teuchos::make_rcp()</a>
+       * was introduced, which should be preferred to this function.
+       */
+#  if defined(DOXYGEN) || !DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0)
+      template <class T, class... Args>
+      Teuchos::RCP<T>
+      make_rcp(Args &&...args);
+#  else
+      using Teuchos::make_rcp;
+#  endif // defined DOXYGEN || !DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0)
+    }    // namespace internal
+  }      // namespace Trilinos
+#endif   // DEAL_II_TRILINOS_WITH_TPETRA
+
 } // namespace Utilities
 
 DEAL_II_NAMESPACE_CLOSE
index d5c4745fb2b7e45b32b9ea00571ca9296b812f04..06f32439f4a759196128925f78e38a54542d582c 100644 (file)
@@ -47,8 +47,11 @@ namespace LinearAlgebra
     template <typename Number>
     Vector<Number>::Vector()
       : Subscriptor()
-      , vector(new VectorType(Teuchos::rcp(
-          new MapType(0, 0, Utilities::Trilinos::tpetra_comm_self()))))
+      , vector(Utilities::Trilinos::internal::make_rcp<VectorType>(
+          Utilities::Trilinos::internal::make_rcp<MapType>(
+            0,
+            0,
+            Utilities::Trilinos::tpetra_comm_self())))
     {}
 
 
@@ -56,7 +59,9 @@ namespace LinearAlgebra
     template <typename Number>
     Vector<Number>::Vector(const Vector<Number> &V)
       : Subscriptor()
-      , vector(new VectorType(V.trilinos_vector(), Teuchos::Copy))
+      , vector(Utilities::Trilinos::internal::make_rcp<VectorType>(
+          V.trilinos_vector(),
+          Teuchos::Copy))
     {}
 
 
@@ -73,7 +78,7 @@ namespace LinearAlgebra
     Vector<Number>::Vector(const IndexSet &parallel_partitioner,
                            const MPI_Comm  communicator)
       : Subscriptor()
-      , vector(new VectorType(
+      , vector(Utilities::Trilinos::internal::make_rcp<VectorType>(
           parallel_partitioner.make_tpetra_map_rcp(communicator, false)))
     {}
 
@@ -89,7 +94,7 @@ namespace LinearAlgebra
         parallel_partitioner.make_tpetra_map_rcp(communicator, false);
 
       if (vector->getMap()->isSameAs(*input_map) == false)
-        vector = Teuchos::rcp(new VectorType(input_map));
+        Utilities::Trilinos::internal::make_rcp<VectorType>(input_map);
       else if (omit_zeroing_entries == false)
         {
           vector->putScalar(0.);
@@ -110,7 +115,7 @@ namespace LinearAlgebra
       Teuchos::RCP<MapType> input_map =
         parallel_partitioner.make_tpetra_map_rcp(communicator, true);
 
-      vector = Teuchos::rcp(new VectorType(input_map));
+      Utilities::Trilinos::internal::make_rcp<VectorType>(input_map);
     }
 
 
@@ -188,7 +193,8 @@ namespace LinearAlgebra
                                Tpetra::REPLACE);
             }
           else
-            vector = Teuchos::rcp(new VectorType(V.trilinos_vector()));
+            vector = Utilities::Trilinos::internal::make_rcp<VectorType>(
+              V.trilinos_vector());
         }
 
       return *this;
@@ -763,10 +769,10 @@ namespace LinearAlgebra
                                                const MPI_Comm  mpi_comm)
     {
       source_stored_elements = source_index_set;
-
-      tpetra_comm_pattern =
-        Teuchos::rcp(new TpetraWrappers::CommunicationPattern(
-          locally_owned_elements(), source_index_set, mpi_comm));
+      tpetra_comm_pattern    = Utilities::Trilinos::internal::make_rcp<
+        TpetraWrappers::CommunicationPattern>(locally_owned_elements(),
+                                              source_index_set,
+                                              mpi_comm);
     }
   } // namespace TpetraWrappers
 } // namespace LinearAlgebra
index 5dd4b6d251f90374ef4231ee0d32d439c502d6ab..590a814b4a31d5a76e01ceb4940d0f0a8712b943 100644 (file)
@@ -995,17 +995,18 @@ IndexSet::make_tpetra_map_rcp(const MPI_Comm communicator,
   const bool linear =
     overlapping ? false : is_ascending_and_one_to_one(communicator);
   if (linear)
-    return Teuchos::RCP<Tpetra::Map<int, types::signed_global_dof_index>>(
-      new Tpetra::Map<int, types::signed_global_dof_index>(
-        size(),
-        n_elements(),
-        0,
+    return Utilities::Trilinos::internal::make_rcp<
+      Tpetra::Map<int, types::signed_global_dof_index>>(
+      size(),
+      n_elements(),
+      0,
 #    ifdef DEAL_II_WITH_MPI
-        Teuchos::rcp(new Teuchos::MpiComm<int>(communicator))
+      Utilities::Trilinos::internal::make_rcp<Teuchos::MpiComm<int>>(
+        communicator)
 #    else
-        Teuchos::rcp(new Teuchos::Comm<int>())
-#    endif
-          ));
+      Utilities::Trilinos::internal::make_rcp<Teuchos::Comm<int>>()
+#    endif // DEAL_WITH_MPI
+    );
   else
     {
       const std::vector<size_type>                indices = get_index_vector();
@@ -1013,17 +1014,19 @@ IndexSet::make_tpetra_map_rcp(const MPI_Comm communicator,
       std::copy(indices.begin(), indices.end(), int_indices.begin());
       const Teuchos::ArrayView<types::signed_global_dof_index> arr_view(
         int_indices);
-      return Teuchos::RCP<Tpetra::Map<int, types::signed_global_dof_index>>(
-        new Tpetra::Map<int, types::signed_global_dof_index>(
-          size(),
-          arr_view,
-          0,
+
+      return Utilities::Trilinos::internal::make_rcp<
+        Tpetra::Map<int, types::signed_global_dof_index>>(
+        size(),
+        arr_view,
+        0,
 #    ifdef DEAL_II_WITH_MPI
-          Teuchos::rcp(new Teuchos::MpiComm<int>(communicator))
+        Utilities::Trilinos::internal::make_rcp<Teuchos::MpiComm<int>>(
+          communicator)
 #    else
-          Teuchos::rcp(new Teuchos::Comm<int>())
-#    endif
-            ));
+        Utilities::Trilinos::internal::make_rcp<Teuchos::Comm<int>>()
+#    endif // DEAL_II_WITH_MPI
+      );
     }
 }
 #  endif
index e9f2b105e5a3bb32f4b6f0de92a390b54ccfbe8b..9fb095e3454fe6ddcaae25c2aab7c7ce078c5785 100644 (file)
@@ -169,6 +169,23 @@ namespace Utilities
     }
   } // namespace Trilinos
 #endif
+
+
+#ifdef DEAL_II_TRILINOS_WITH_TPETRA
+  namespace Trilinos
+  {
+#  if !DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0)
+    template <class T, class... Args>
+    Teuchos::RCP<T>
+    make_rcp(Args &&...args)
+    {
+      return Teuchos::RCP<T>(new T(std::forward<Args>(args)...));
+    }
+#  endif // !DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0)
+
+  }    // namespace Trilinos
+#endif // DEAL_II_TRILINOS_WITH_TPETRA
+
 } // namespace Utilities
 
 DEAL_II_NAMESPACE_CLOSE

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.