Skip to content

Commit

Permalink
Merge pull request #16270 from kinnewig/make_teuchos_rcp
Browse files Browse the repository at this point in the history
Introduce make_rcp(...)
  • Loading branch information
bangerth committed Nov 19, 2023
2 parents 6d5f499 + 2a2db68 commit 4f8abe4
Show file tree
Hide file tree
Showing 6 changed files with 93 additions and 29 deletions.
4 changes: 4 additions & 0 deletions doc/news/changes/minor/20231115SebastianKinnewig
Original file line number Diff line number Diff line change
@@ -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)
1 change: 1 addition & 0 deletions include/deal.II/base/index_set.h
Original file line number Diff line number Diff line change
Expand Up @@ -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>

Expand Down
33 changes: 33 additions & 0 deletions include/deal.II/base/trilinos_utilities.h
Original file line number Diff line number Diff line change
Expand Up @@ -32,6 +32,10 @@
# endif
#endif

#ifdef DEAL_II_TRILINOS_WITH_TPETRA
# include <Teuchos_RCPDecl.hpp>
#endif // DEAL_II_TRILINOS_WITH_TPETRA

DEAL_II_NAMESPACE_OPEN

/**
Expand Down Expand Up @@ -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
Expand Down
28 changes: 17 additions & 11 deletions include/deal.II/lac/trilinos_tpetra_vector.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -47,16 +47,21 @@ 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())))
{}



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))
{}


Expand All @@ -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)))
{}

Expand All @@ -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.);
Expand All @@ -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);
}


Expand Down Expand Up @@ -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;
Expand Down Expand Up @@ -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
Expand Down
39 changes: 21 additions & 18 deletions source/base/index_set.cc
Original file line number Diff line number Diff line change
Expand Up @@ -995,35 +995,38 @@ 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();
std::vector<types::signed_global_dof_index> int_indices(indices.size());
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
Expand Down
17 changes: 17 additions & 0 deletions source/base/trilinos_utilities.cc
Original file line number Diff line number Diff line change
Expand Up @@ -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

0 comments on commit 4f8abe4

Please sign in to comment.