Skip to content

Commit

Permalink
TpetraWrappers::Vector: Don't use rcp pointers internally
Browse files Browse the repository at this point in the history
  • Loading branch information
masterleinad committed Jan 26, 2024
1 parent ec7c27e commit 23ecc6c
Show file tree
Hide file tree
Showing 4 changed files with 143 additions and 198 deletions.
2 changes: 1 addition & 1 deletion include/deal.II/lac/affine_constraints.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -2671,7 +2671,7 @@ namespace internal
parallel_partitioner.add_indices(needed_elements);

const MPI_Comm mpi_comm = Utilities::Trilinos::teuchos_comm_to_mpi_comm(
vec.trilinos_rcp()->getMap()->getComm());
vec.trilinos_vector().getMap()->getComm());

output.reinit(locally_owned_elements, needed_elements, mpi_comm);

Expand Down
30 changes: 16 additions & 14 deletions include/deal.II/lac/trilinos_tpetra_sparse_matrix.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -590,11 +590,11 @@ namespace LinearAlgebra
{
Assert(&src != &dst, ExcSourceEqualsDestination());
Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
Assert(src.trilinos_rcp()->getMap()->isSameAs(*matrix->getDomainMap()),
Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
ExcColMapMissmatch());
Assert(dst.trilinos_rcp()->getMap()->isSameAs(*matrix->getRangeMap()),
Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
ExcDomainMapMissmatch());
matrix->apply(*src.trilinos_rcp(), *dst.trilinos_rcp());
matrix->apply(src.trilinos_vector(), dst.trilinos_vector());
}


Expand All @@ -606,11 +606,13 @@ namespace LinearAlgebra
{
Assert(&src != &dst, ExcSourceEqualsDestination());
Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
Assert(dst.trilinos_rcp()->getMap()->isSameAs(*matrix->getDomainMap()),
Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
ExcColMapMissmatch());
Assert(src.trilinos_rcp()->getMap()->isSameAs(*matrix->getRangeMap()),
Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
ExcDomainMapMissmatch());
matrix->apply(*src.trilinos_rcp(), *dst.trilinos_rcp(), Teuchos::TRANS);
matrix->apply(src.trilinos_vector(),
dst.trilinos_vector(),
Teuchos::TRANS);
}


Expand All @@ -623,12 +625,12 @@ namespace LinearAlgebra
{
Assert(&src != &dst, ExcSourceEqualsDestination());
Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
Assert(src.trilinos_rcp()->getMap()->isSameAs(*matrix->getDomainMap()),
Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
ExcColMapMissmatch());
Assert(dst.trilinos_rcp()->getMap()->isSameAs(*matrix->getRangeMap()),
Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
ExcDomainMapMissmatch());
matrix->apply(*src.trilinos_rcp(),
*dst.trilinos_rcp(),
matrix->apply(src.trilinos_vector(),
dst.trilinos_vector(),
Teuchos::NO_TRANS,
Teuchos::ScalarTraits<Number>::one(),
Teuchos::ScalarTraits<Number>::one());
Expand All @@ -644,12 +646,12 @@ namespace LinearAlgebra
{
Assert(&src != &dst, ExcSourceEqualsDestination());
Assert(matrix->isFillComplete(), ExcMatrixNotCompressed());
Assert(dst.trilinos_rcp()->getMap()->isSameAs(*matrix->getDomainMap()),
Assert(dst.trilinos_vector().getMap()->isSameAs(*matrix->getDomainMap()),
ExcColMapMissmatch());
Assert(src.trilinos_rcp()->getMap()->isSameAs(*matrix->getRangeMap()),
Assert(src.trilinos_vector().getMap()->isSameAs(*matrix->getRangeMap()),
ExcDomainMapMissmatch());
matrix->apply(*src.trilinos_rcp(),
*dst.trilinos_rcp(),
matrix->apply(src.trilinos_vector(),
dst.trilinos_vector(),
Teuchos::TRANS,
Teuchos::ScalarTraits<Number>::one(),
Teuchos::ScalarTraits<Number>::one());
Expand Down
82 changes: 31 additions & 51 deletions include/deal.II/lac/trilinos_tpetra_vector.h
Original file line number Diff line number Diff line change
Expand Up @@ -282,11 +282,6 @@ namespace LinearAlgebra
*/
Vector(const Vector &V);

/**
* Copy constructor from Teuchos::RCP<Tpetra::Vector>.
*/
Vector(const Teuchos::RCP<VectorType> V);

/**
* TODO: This is not used
* This constructor takes an IndexSet that defines how to distribute the
Expand Down Expand Up @@ -747,21 +742,6 @@ namespace LinearAlgebra
Tpetra::Vector<Number, int, types::signed_global_dof_index> &
trilinos_vector();

/**
* Return a const Teuchos::RCP to the underlying Trilinos
* Tpetra::Vector class.
*/
Teuchos::RCP<
const Tpetra::Vector<Number, int, types::signed_global_dof_index>>
trilinos_rcp() const;

/**
* Return a (modifiable) Teuchos::RCP to the underlying Trilinos
* Tpetra::Vector class.
*/
Teuchos::RCP<Tpetra::Vector<Number, int, types::signed_global_dof_index>>
trilinos_rcp();

/**
* Prints the vector to the output stream @p out.
*/
Expand Down Expand Up @@ -881,14 +861,14 @@ namespace LinearAlgebra
/**
* Teuchos::RCP to the actual Tpetra vector object.
*/
Teuchos::RCP<VectorType> vector;
VectorType vector;

/**
* A vector object in Trilinos to be used for collecting the non-local
* elements if the vector was constructed with an additional IndexSet
* describing ghost elements.
*/
Teuchos::RCP<VectorType> nonlocal_vector;
VectorType nonlocal_vector;

/**
* IndexSet of the elements of the last imported vector.
Expand Down Expand Up @@ -948,37 +928,37 @@ namespace LinearAlgebra
const Number *values)
{
# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
auto vector_2d_local = vector->template getLocalView<Kokkos::HostSpace>(
auto vector_2d_local = vector.template getLocalView<Kokkos::HostSpace>(
Tpetra::Access::ReadWrite);
auto vector_2d_nonlocal =
nonlocal_vector->template getLocalView<Kokkos::HostSpace>(
nonlocal_vector.template getLocalView<Kokkos::HostSpace>(
Tpetra::Access::ReadWrite);
# else
vector->template sync<Kokkos::HostSpace>();
auto vector_2d_local = vector->template getLocalView<Kokkos::HostSpace>();
vector.template sync<Kokkos::HostSpace>();
auto vector_2d_local = vector.template getLocalView<Kokkos::HostSpace>();
auto vector_2d_nonlocal =
nonlocal_vector->template getLocalView<Kokkos::HostSpace>();
nonlocal_vector.template getLocalView<Kokkos::HostSpace>();
# endif

auto vector_1d_local = Kokkos::subview(vector_2d_local, Kokkos::ALL(), 0);
auto vector_1d_nonlocal =
Kokkos::subview(vector_2d_nonlocal, Kokkos::ALL(), 0);
# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
vector->template modify<Kokkos::HostSpace>();
nonlocal_vector->template modify<Kokkos::HostSpace>();
vector.template modify<Kokkos::HostSpace>();
nonlocal_vector.template modify<Kokkos::HostSpace>();
# endif

for (size_type i = 0; i < n_elements; ++i)
{
const size_type row = indices[i];
TrilinosWrappers::types::int_type local_row =
vector->getMap()->getLocalElement(row);
vector.getMap()->getLocalElement(row);

// check if the index is in the non local set
bool nonlocal = false;
if (local_row == Teuchos::OrdinalTraits<int>::invalid())
{
local_row = nonlocal_vector->getMap()->getLocalElement(row);
local_row = nonlocal_vector.getMap()->getLocalElement(row);
nonlocal = true;
compressed = false;
}
Expand All @@ -987,16 +967,16 @@ namespace LinearAlgebra
Assert(
local_row != Teuchos::OrdinalTraits<int>::invalid(),
ExcAccessToNonLocalElement(row,
vector->getMap()->getLocalNumElements(),
vector->getMap()->getMinLocalIndex(),
vector->getMap()->getMaxLocalIndex()));
vector.getMap()->getLocalNumElements(),
vector.getMap()->getMinLocalIndex(),
vector.getMap()->getMaxLocalIndex()));
# else
Assert(
local_row != Teuchos::OrdinalTraits<int>::invalid(),
ExcAccessToNonLocalElement(row,
vector->getMap()->getNodeNumElements(),
vector->getMap()->getMinLocalIndex(),
vector->getMap()->getMaxLocalIndex()));
vector.getMap()->getNodeNumElements(),
vector.getMap()->getMinLocalIndex(),
vector.getMap()->getMaxLocalIndex()));

# endif

Expand All @@ -1010,10 +990,10 @@ namespace LinearAlgebra
}

# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
vector->template sync<
vector.template sync<
typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
device_type::memory_space>();
nonlocal_vector->template sync<
nonlocal_vector.template sync<
typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
device_type::memory_space>();
# endif
Expand All @@ -1032,45 +1012,45 @@ namespace LinearAlgebra
Assert(!has_ghost_elements(), ExcGhostsPresent());

# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
auto vector_2d = vector->template getLocalView<Kokkos::HostSpace>(
auto vector_2d = vector.template getLocalView<Kokkos::HostSpace>(
Tpetra::Access::ReadWrite);
# else
vector->template sync<Kokkos::HostSpace>();
auto vector_2d = vector->template getLocalView<Kokkos::HostSpace>();
vector.template sync<Kokkos::HostSpace>();
auto vector_2d = vector.template getLocalView<Kokkos::HostSpace>();
# endif
auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0);
# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
vector->template modify<Kokkos::HostSpace>();
vector.template modify<Kokkos::HostSpace>();
# endif

for (size_type i = 0; i < n_elements; ++i)
{
const size_type row = indices[i];
TrilinosWrappers::types::int_type local_row =
vector->getMap()->getLocalElement(row);
vector.getMap()->getLocalElement(row);

# if DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0)
Assert(
local_row != Teuchos::OrdinalTraits<int>::invalid(),
ExcAccessToNonLocalElement(row,
vector->getMap()->getLocalNumElements(),
vector->getMap()->getMinLocalIndex(),
vector->getMap()->getMaxLocalIndex()));
vector.getMap()->getLocalNumElements(),
vector.getMap()->getMinLocalIndex(),
vector.getMap()->getMaxLocalIndex()));
# else
Assert(
local_row != Teuchos::OrdinalTraits<int>::invalid(),
ExcAccessToNonLocalElement(row,
vector->getMap()->getNodeNumElements(),
vector->getMap()->getMinLocalIndex(),
vector->getMap()->getMaxLocalIndex()));
vector.getMap()->getNodeNumElements(),
vector.getMap()->getMinLocalIndex(),
vector.getMap()->getMaxLocalIndex()));
# endif

if (local_row != Teuchos::OrdinalTraits<int>::invalid())
vector_1d(local_row) = values[i];
}

# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
vector->template sync<
vector.template sync<
typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
device_type::memory_space>();
# endif
Expand Down

0 comments on commit 23ecc6c

Please sign in to comment.