15#ifndef dealii_trilinos_tpetra_vector_h
16#define dealii_trilinos_tpetra_vector_h
25#ifdef DEAL_II_TRILINOS_WITH_TPETRA
37# include <Teuchos_Comm.hpp>
38# include <Teuchos_OrdinalTraits.hpp>
39# include <Tpetra_Core.hpp>
40# include <Tpetra_Vector.hpp>
41# include <Tpetra_Version.hpp>
53template <
typename Number>
57# ifdef HAVE_TPETRA_INST_FLOAT
63# ifdef HAVE_TPETRA_INST_DOUBLE
69# ifdef DEAL_II_WITH_COMPLEX_VALUES
70# ifdef HAVE_TPETRA_INST_COMPLEX_FLOAT
76# ifdef HAVE_TPETRA_INST_COMPLEX_DOUBLE
87 template <
typename Number>
88 class ReadWriteVector;
105 namespace TpetraWrappers
144 template <
typename Number,
146 class VectorReference
154 const size_type index);
160 VectorReference(
const VectorReference &) =
default;
173 const VectorReference &
174 operator=(
const VectorReference &r)
const;
180 operator=(
const VectorReference &r);
185 const VectorReference &
186 operator=(
const Number &s)
const;
191 const VectorReference &
192 operator+=(
const Number &s)
const;
197 const VectorReference &
198 operator-=(
const Number &s)
const;
203 const VectorReference &
204 operator*=(
const Number &s)
const;
209 const VectorReference &
210 operator/=(
const Number &s)
const;
216 operator Number()
const;
223 <<
"An error with error number " << arg1
224 <<
" occurred while calling a Trilinos function");
232 ExcAccessToNonLocalElement,
237 <<
"You are trying to access element " << arg1
238 <<
" of a distributed vector, but this element is not stored "
239 <<
"on the current processor. Note: There are " << arg2
240 <<
" elements stored "
241 <<
"on the current processor from within the range [" << arg3 <<
','
242 << arg4 <<
"] but Trilinos vectors need not store contiguous "
243 <<
"ranges on each processor, and not every element in "
244 <<
"this range may in fact be stored locally."
246 <<
"A common source for this kind of problem is that you "
247 <<
"are passing a 'fully distributed' vector into a function "
248 <<
"that needs read access to vector elements that correspond "
249 <<
"to degrees of freedom on ghost cells (or at least to "
250 <<
"'locally active' degrees of freedom that are not also "
251 <<
"'locally owned'). You need to pass a vector that has these "
252 <<
"elements as ghost entries.");
263 const size_type index;
287 template <
typename Number,
typename MemorySpace = ::MemorySpace::Host>
297 using reference = internal::VectorReference<Number, MemorySpace>;
299 const internal::VectorReference<Number, MemorySpace>;
351 const bool vector_writable =
false);
368 const MPI_Comm communicator = MPI_COMM_WORLD,
369 const bool omit_zeroing_entries =
false);
387 const IndexSet &locally_relevant_or_ghost_entries,
388 const MPI_Comm communicator = MPI_COMM_WORLD,
389 const bool vector_writable =
false);
397 const bool omit_zeroing_entries =
false);
462 template <
typename OtherNumber>
485 const Teuchos::RCP<const Utilities::MPI::CommunicationPatternBase>
486 &communication_pattern);
496 const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
497 &communication_pattern);
516 std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
517 communication_pattern = {})
633 add(
const std::vector<size_type> &indices,
634 const std::vector<Number> &values);
645 const Number *values);
665 const Number *values);
834 std::pair<size_type, size_type>
918 Teuchos::RCP<const TpetraTypes::VectorType<Number, MemorySpace>>
925 Teuchos::RCP<TpetraTypes::VectorType<Number, MemorySpace>>
933 const unsigned int precision = 3,
934 const bool scientific =
true,
935 const bool across =
true)
const;
978 <<
"You are trying to access element " << arg1
979 <<
" of a distributed vector, but this element is not stored "
980 <<
"on the current processor. Note: There are " << arg2
981 <<
" elements stored "
982 <<
"on the current processor from within the range [" << arg3 <<
','
983 << arg4 <<
"] but Trilinos vectors need not store contiguous "
984 <<
"ranges on each processor, and not every element in "
985 <<
"this range may in fact be stored locally."
987 <<
"A common source for this kind of problem is that you "
988 <<
"are passing a 'fully distributed' vector into a function "
989 <<
"that needs read access to vector elements that correspond "
990 <<
"to degrees of freedom on ghost cells (or at least to "
991 <<
"'locally active' degrees of freedom that are not also "
992 <<
"'locally owned'). You need to pass a vector that has these "
993 <<
"elements as ghost entries.");
1001 "To compress a vector, a locally_relevant_dofs "
1002 "index set, and a locally_owned_dofs index set "
1003 "must be provided. These index sets must be "
1004 "provided either when the vector is initialized "
1005 "or when compress is called. See the documentation "
1006 "of compress() for more information.");
1015 <<
"An error with error number " << arg1
1016 <<
" occurred while calling a Trilinos function");
1047 Teuchos::RCP<TpetraTypes::VectorType<Number, MemorySpace>>
vector;
1054 Teuchos::RCP<TpetraTypes::VectorType<Number, MemorySpace>>
1066 Teuchos::RCP<const TpetraWrappers::CommunicationPattern<MemorySpace>>
1070 friend class internal::VectorReference<Number,
MemorySpace>;
1076 template <
typename Number,
typename MemorySpace>
1085 template <
typename Number,
typename MemorySpace>
1094 template <
typename Number,
typename MemorySpace>
1101 template <
typename Number,
typename MemorySpace>
1105 vector.swap(v.vector);
1109 template <
typename Number,
typename MemorySpace>
1112 const std::vector<Number> &values)
1118 add(indices.size(), indices.data(), values.data());
1123 template <
typename Number,
typename MemorySpace>
1127 const Number *values)
1133# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1134 auto vector_2d_local = vector->template getLocalView<Kokkos::HostSpace>(
1135 Tpetra::Access::ReadWrite);
1137 vector->template sync<Kokkos::HostSpace>();
1139 auto vector_2d_local = vector->template getLocalView<Kokkos::HostSpace>();
1148 auto vector_1d_local = Kokkos::subview(vector_2d_local, Kokkos::ALL(), 0);
1149 using ViewType1d =
decltype(vector_1d_local);
1150 std::optional<ViewType1d> vector_1d_nonlocal;
1152# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1155 vector->template modify<Kokkos::HostSpace>();
1158 for (
size_type i = 0; i < n_elements; ++i)
1166 vector->getMap()->getLocalElement(row);
1167 local_row != Teuchos::OrdinalTraits<int>::invalid())
1169 vector_1d_local(local_row) += values[i];
1178 nonlocal_vector->getMap()->getLocalElement(row);
1180# if DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0)
1181 Assert(nonlocal_row != Teuchos::OrdinalTraits<int>::invalid(),
1182 ExcAccessToNonLocalElement(
1184 vector->getMap()->getLocalNumElements(),
1185 vector->getMap()->getMinLocalIndex(),
1186 vector->getMap()->getMaxLocalIndex()));
1188 Assert(nonlocal_row != Teuchos::OrdinalTraits<int>::invalid(),
1189 ExcAccessToNonLocalElement(
1191 vector->getMap()->getNodeNumElements(),
1192 vector->getMap()->getMinLocalIndex(),
1193 vector->getMap()->getMaxLocalIndex()));
1201 if (!vector_1d_nonlocal)
1203# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1204 auto vector_2d_nonlocal =
1205 nonlocal_vector->template getLocalView<Kokkos::HostSpace>(
1206 Tpetra::Access::ReadWrite);
1208 auto vector_2d_nonlocal =
1209 nonlocal_vector->template getLocalView<Kokkos::HostSpace>();
1212 vector_1d_nonlocal =
1213 Kokkos::subview(vector_2d_nonlocal, Kokkos::ALL(), 0);
1215# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1217 nonlocal_vector->template modify<Kokkos::HostSpace>();
1220 (*vector_1d_nonlocal)(nonlocal_row) += values[i];
1225# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1226 vector->template sync<
1227 typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
1228 device_type::memory_space>();
1232 if (vector_1d_nonlocal)
1233 nonlocal_vector->template sync<
1234 typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
1235 device_type::memory_space>();
1241 template <
typename Number,
typename MemorySpace>
1245 const Number *values)
1251# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1252 auto vector_2d_local = vector->template getLocalView<Kokkos::HostSpace>(
1253 Tpetra::Access::ReadWrite);
1255 vector->template sync<Kokkos::HostSpace>();
1257 auto vector_2d_local = vector->template getLocalView<Kokkos::HostSpace>();
1266 auto vector_1d_local = Kokkos::subview(vector_2d_local, Kokkos::ALL(), 0);
1267 using ViewType1d =
decltype(vector_1d_local);
1268 std::optional<ViewType1d> vector_1d_nonlocal;
1270# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1273 vector->template modify<Kokkos::HostSpace>();
1276 for (
size_type i = 0; i < n_elements; ++i)
1284 vector->getMap()->getLocalElement(row);
1285 local_row != Teuchos::OrdinalTraits<int>::invalid())
1287 vector_1d_local(local_row) = values[i];
1296 nonlocal_vector->getMap()->getLocalElement(row);
1298# if DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0)
1299 Assert(nonlocal_row != Teuchos::OrdinalTraits<int>::invalid(),
1300 ExcAccessToNonLocalElement(
1302 vector->getMap()->getLocalNumElements(),
1303 vector->getMap()->getMinLocalIndex(),
1304 vector->getMap()->getMaxLocalIndex()));
1306 Assert(nonlocal_row != Teuchos::OrdinalTraits<int>::invalid(),
1307 ExcAccessToNonLocalElement(
1309 vector->getMap()->getNodeNumElements(),
1310 vector->getMap()->getMinLocalIndex(),
1311 vector->getMap()->getMaxLocalIndex()));
1319 if (!vector_1d_nonlocal)
1321# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1322 auto vector_2d_nonlocal =
1323 nonlocal_vector->template getLocalView<Kokkos::HostSpace>(
1324 Tpetra::Access::ReadWrite);
1326 auto vector_2d_nonlocal =
1327 nonlocal_vector->template getLocalView<Kokkos::HostSpace>();
1330 vector_1d_nonlocal =
1331 Kokkos::subview(vector_2d_nonlocal, Kokkos::ALL(), 0);
1333# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1335 nonlocal_vector->template modify<Kokkos::HostSpace>();
1338 (*vector_1d_nonlocal)(nonlocal_row) = values[i];
1343# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
1344 vector->template sync<
1345 typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
1346 device_type::memory_space>();
1350 if (vector_1d_nonlocal)
1351 nonlocal_vector->template sync<
1352 typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
1353 device_type::memory_space>();
1359 template <
typename Number,
typename MemorySpace>
1360 inline internal::VectorReference<Number, MemorySpace>
1363 return internal::VectorReference(*
this, index);
1366 template <
typename Number,
typename MemorySpace>
1367 inline internal::VectorReference<Number, MemorySpace>
1370 return operator()(index);
1373 template <
typename Number,
typename MemorySpace>
1377 return operator()(index);
1385 template <
typename Number,
typename MemorySpace>
1386 inline VectorReference<Number, MemorySpace>::VectorReference(
1388 const size_type index)
1395 template <
typename Number,
typename MemorySpace>
1396 inline const VectorReference<Number, MemorySpace> &
1397 VectorReference<Number, MemorySpace>::operator=(
1398 const VectorReference<Number, MemorySpace> &r)
const
1404 *
this =
static_cast<Number
>(r);
1411 template <
typename Number,
typename MemorySpace>
1412 inline VectorReference<Number, MemorySpace> &
1413 VectorReference<Number, MemorySpace>::operator=(
1414 const VectorReference<Number, MemorySpace> &r)
1417 *
this =
static_cast<Number
>(r);
1424 template <
typename Number,
typename MemorySpace>
1425 inline const VectorReference<Number, MemorySpace> &
1426 VectorReference<Number, MemorySpace>::operator=(
const Number &value)
const
1428 vector.set(1, &index, &value);
1434 template <
typename Number,
typename MemorySpace>
1435 inline const VectorReference<Number, MemorySpace> &
1436 VectorReference<Number, MemorySpace>::operator+=(
1437 const Number &value)
const
1439 vector.add(1, &index, &value);
1445 template <
typename Number,
typename MemorySpace>
1446 inline const VectorReference<Number, MemorySpace> &
1447 VectorReference<Number, MemorySpace>::operator-=(
1448 const Number &value)
const
1450 Number new_value = -
value;
1451 vector.add(1, &index, &new_value);
1457 template <
typename Number,
typename MemorySpace>
1458 inline const VectorReference<Number, MemorySpace> &
1459 VectorReference<Number, MemorySpace>::operator*=(
1460 const Number &value)
const
1462 Number new_value =
static_cast<Number
>(*this) *
value;
1463 vector.set(1, &index, &new_value);
1469 template <
typename Number,
typename MemorySpace>
1470 inline const VectorReference<Number, MemorySpace> &
1471 VectorReference<Number, MemorySpace>::operator/=(
1472 const Number &value)
const
1474 Number new_value =
static_cast<Number
>(*this) /
value;
1475 vector.set(1, &index, &new_value);
1491 namespace LinearOperatorImplementation
1500 template <
typename Number,
typename MemorySpace>
1505 template <
typename Matrix>
1508 const Matrix &matrix,
1510 bool omit_zeroing_entries)
1512 v.
reinit(matrix.locally_owned_range_indices(),
1513 matrix.get_mpi_communicator(),
1514 omit_zeroing_entries);
1517 template <
typename Matrix>
1520 const Matrix &matrix,
1522 bool omit_zeroing_entries)
1524 v.
reinit(matrix.locally_owned_domain_indices(),
1525 matrix.get_mpi_communicator(),
1526 omit_zeroing_entries);
1536template <
typename Number,
typename MemorySpace>
1538 LinearAlgebra::TpetraWrappers::Vector<Number, MemorySpace>> : std::false_type
void reinit(const Vector< Number, MemorySpace > &V, const bool omit_zeroing_entries=false)
void equ(const Number a, const Vector< Number, MemorySpace > &V)
Teuchos::RCP< const TpetraWrappers::CommunicationPattern< MemorySpace > > tpetra_comm_pattern
void reinit(const IndexSet &locally_owned_entries, const IndexSet &locally_relevant_or_ghost_entries, const MPI_Comm communicator=MPI_COMM_WORLD, const bool vector_writable=false)
void add(const Number a, const Vector< Number, MemorySpace > &V, const Number b, const Vector< Number, MemorySpace > &W)
bool is_non_negative() const
Number add_and_dot(const Number a, const Vector< Number, MemorySpace > &V, const Vector< Number, MemorySpace > &W)
TpetraTypes::VectorType< Number, MemorySpace > & trilinos_vector()
void reinit(const IndexSet ¶llel_partitioner, const MPI_Comm communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false)
MPI_Comm get_mpi_communicator() const
Teuchos::RCP< TpetraTypes::VectorType< Number, MemorySpace > > vector
Number operator[](const size_type index) const
bool operator==(const Vector< Number, MemorySpace > &v) const
Vector(const IndexSet ¶llel_partitioner, const MPI_Comm communicator)
void scale(const Vector< Number, MemorySpace > &scaling_factors)
void set(const size_type n_elements, const size_type *indices, const Number *values)
void import_elements(const ReadWriteVector< Number > &V, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern)
MPI_Comm mpi_comm() const
Vector(const Teuchos::RCP< TpetraTypes::VectorType< Number, MemorySpace > > V)
internal::VectorReference< Number, MemorySpace > reference
size_type locally_owned_size() const
reference operator()(const size_type index)
virtual void swap(Vector &v) noexcept
void compress(const VectorOperation::values operation)
bool is_compressed() const
bool operator!=(const Vector< Number, MemorySpace > &v) const
Number operator()(const size_type index) const
Vector & operator/=(const Number factor)
std::pair< size_type, size_type > local_range() const
virtual void extract_subvector_to(const ArrayView< const types::global_dof_index > &indices, ArrayView< Number > &elements) const override
::IndexSet source_stored_elements
void sadd(const Number s, const Number a, const Vector< Number, MemorySpace > &V)
reference operator[](const size_type index)
void add(const Number a, const Vector< Number, MemorySpace > &V)
Vector & operator*=(const Number factor)
void add(const size_type n_elements, const size_type *indices, const Number *values)
void import_elements(const ReadWriteVector< Number > &V, VectorOperation::values operation, const Teuchos::RCP< const Utilities::MPI::CommunicationPatternBase > &communication_pattern)
real_type l2_norm() const
void add(const std::vector< size_type > &indices, const std::vector< Number > &values)
typename numbers::NumberTraits< Number >::real_type real_type
Vector & operator-=(const Vector< Number, MemorySpace > &V)
Vector & operator=(const ::Vector< OtherNumber > &V)
::IndexSet locally_owned_elements() const
Vector & operator=(const Vector &V)
Vector & operator=(const Number s)
Vector & operator+=(const Vector< Number, MemorySpace > &V)
const TpetraTypes::VectorType< Number, MemorySpace > & trilinos_vector() const
real_type linfty_norm() const
real_type norm_sqr() const
virtual size_type size() const override
void import_elements(const ReadWriteVector< Number > &V, VectorOperation::values operation)
void create_tpetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm mpi_comm)
real_type l1_norm() const
Teuchos::RCP< const TpetraTypes::VectorType< Number, MemorySpace > > trilinos_rcp() const
const internal::VectorReference< Number, MemorySpace > const_reference
Number operator*(const Vector< Number, MemorySpace > &V) const
Number mean_value() const
std::size_t memory_consumption() const
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
bool has_ghost_elements() const
bool in_local_range(const size_type index) const
Vector(const IndexSet &locally_owned_entries, const IndexSet &ghost_entries, const MPI_Comm communicator, const bool vector_writable=false)
Teuchos::RCP< TpetraTypes::VectorType< Number, MemorySpace > > nonlocal_vector
Teuchos::RCP< TpetraTypes::VectorType< Number, MemorySpace > > trilinos_rcp()
types::global_dof_index size_type
static void reinit_range_vector(const Matrix &matrix, LinearAlgebra::TpetraWrappers::Vector< Number, MemorySpace > &v, bool omit_zeroing_entries)
static void reinit_domain_vector(const Matrix &matrix, LinearAlgebra::TpetraWrappers::Vector< Number, MemorySpace > &v, bool omit_zeroing_entries)
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcGhostsPresent()
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
static ::ExceptionBase & ExcVectorTypeNotCompatible()
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define Assert(cond, exc)
static ::ExceptionBase & ExcMissingIndexSet()
#define AssertDimension(dim1, dim2)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcDifferentParallelPartitioning()
static ::ExceptionBase & ExcTrilinosError(int arg1)
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
void swap(Vector< Number, MemorySpace > &u, Vector< Number, MemorySpace > &v) noexcept
unsigned int global_dof_index