--- /dev/null
+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)
#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>
# endif
#endif
+#ifdef DEAL_II_TRILINOS_WITH_TPETRA
+# include <Teuchos_RCPDecl.hpp>
+#endif // DEAL_II_TRILINOS_WITH_TPETRA
+
DEAL_II_NAMESPACE_OPEN
/**
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
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())))
{}
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))
{}
Vector<Number>::Vector(const IndexSet ¶llel_partitioner,
const MPI_Comm communicator)
: Subscriptor()
- , vector(new VectorType(
+ , vector(Utilities::Trilinos::internal::make_rcp<VectorType>(
parallel_partitioner.make_tpetra_map_rcp(communicator, false)))
{}
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.);
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);
}
Tpetra::REPLACE);
}
else
- vector = Teuchos::rcp(new VectorType(V.trilinos_vector()));
+ vector = Utilities::Trilinos::internal::make_rcp<VectorType>(
+ V.trilinos_vector());
}
return *this;
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
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();
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
}
} // 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