15#ifndef dealii_petsc_vector_base_h
16#define dealii_petsc_vector_base_h
21#ifdef DEAL_II_WITH_PETSC
30# include <boost/serialization/split_member.hpp>
31# include <boost/serialization/utility.hpp>
42template <
typename number>
95 VectorReference(
const VectorBase &vector,
const size_type index);
101 VectorReference(
const VectorReference &vector) =
default;
114 const VectorReference &
115 operator=(
const VectorReference &r)
const;
123 operator=(
const VectorReference &r);
128 const VectorReference &
129 operator=(
const PetscScalar &s)
const;
134 const VectorReference &
135 operator+=(
const PetscScalar &s)
const;
140 const VectorReference &
141 operator-=(
const PetscScalar &s)
const;
146 const VectorReference &
147 operator*=(
const PetscScalar &s)
const;
152 const VectorReference &
153 operator/=(
const PetscScalar &s)
const;
174 operator PetscScalar()
const;
179 ExcAccessToNonlocalElement,
183 <<
"You tried to access element " << arg1
184 <<
" of a distributed vector, but only elements in range [" << arg2
185 <<
',' << arg3 <<
"] are stored locally and can be accessed."
187 <<
"A common source for this kind of problem is that you "
188 <<
"are passing a 'fully distributed' vector into a function "
189 <<
"that needs read access to vector elements that correspond "
190 <<
"to degrees of freedom on ghost cells (or at least to "
191 <<
"'locally active' degrees of freedom that are not also "
192 <<
"'locally owned'). You need to pass a vector that has these "
193 <<
"elements as ghost entries.");
200 <<
"You tried to do a "
201 << (arg1 == 1 ?
"'set'" : (arg1 == 2 ?
"'add'" :
"???"))
202 <<
" operation but the vector is currently in "
203 << (arg2 == 1 ?
"'set'" : (arg2 == 2 ?
"'add'" :
"???"))
204 <<
" mode. You first have to call 'compress()'.");
210 const VectorBase &vector;
215 const size_type
index;
219 friend class ::PETScWrappers::VectorBase;
360 size()
const override;
381 std::pair<size_type, size_type>
463 set(
const std::vector<size_type> &indices,
464 const std::vector<PetscScalar> &values);
483 std::vector<PetscScalar> &values)
const;
520 template <
typename ForwardIterator,
typename OutputIterator>
523 const ForwardIterator indices_end,
524 OutputIterator values_begin)
const;
531 add(
const std::vector<size_type> &indices,
532 const std::vector<PetscScalar> &values);
539 add(
const std::vector<size_type> &indices,
540 const ::Vector<PetscScalar> &values);
550 const PetscScalar *values);
663 add(
const PetscScalar s);
675 add(
const PetscScalar a,
690 sadd(
const PetscScalar s,
const PetscScalar a,
const VectorBase &V);
714 write_ascii(
const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
724 print(std::ostream &out,
725 const unsigned int precision = 3,
726 const bool scientific =
true,
727 const bool across =
true)
const;
736 template <
class Archive>
738 save(Archive &ar,
const unsigned int version)
const;
745 template <
class Archive>
747 load(Archive &ar,
const unsigned int version);
755 template <
class Archive>
757 serialize(Archive &archive,
const unsigned int version);
761 BOOST_SERIALIZATION_SPLIT_MEMBER()
786 operator const Vec &()
const;
847 const PetscScalar *values,
848 const bool add_values);
877 inline VectorReference::VectorReference(
const VectorBase &vector,
878 const size_type index)
884 inline const VectorReference &
885 VectorReference::operator=(
const VectorReference &r)
const
891 *
this =
static_cast<PetscScalar
>(r);
898 inline VectorReference &
899 VectorReference::operator=(
const VectorReference &r)
905 *
this =
static_cast<PetscScalar
>(r);
912 inline const VectorReference &
913 VectorReference::operator=(
const PetscScalar &value)
const
921 const PetscInt petsc_i =
index;
923 const PetscErrorCode ierr =
924 VecSetValues(vector, 1, &petsc_i, &value, INSERT_VALUES);
934 inline const VectorReference &
935 VectorReference::operator+=(
const PetscScalar &value)
const
952 if (value == PetscScalar())
956 const PetscInt petsc_i =
index;
957 const PetscErrorCode ierr =
958 VecSetValues(vector, 1, &petsc_i, &value, ADD_VALUES);
967 inline const VectorReference &
968 VectorReference::operator-=(
const PetscScalar &value)
const
985 if (value == PetscScalar())
990 const PetscInt petsc_i =
index;
991 const PetscScalar subtractand = -
value;
992 const PetscErrorCode ierr =
993 VecSetValues(vector, 1, &petsc_i, &subtractand, ADD_VALUES);
1001 inline const VectorReference &
1002 VectorReference::operator*=(
const PetscScalar &value)
const
1022 const PetscInt petsc_i =
index;
1023 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) *
value;
1025 const PetscErrorCode ierr =
1026 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1034 inline const VectorReference &
1035 VectorReference::operator/=(
const PetscScalar &value)
const
1055 const PetscInt petsc_i =
index;
1056 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) /
value;
1058 const PetscErrorCode ierr =
1059 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1068 VectorReference::real()
const
1070# ifndef PETSC_USE_COMPLEX
1071 return static_cast<PetscScalar
>(*this);
1073 return PetscRealPart(
static_cast<PetscScalar
>(*
this));
1080 VectorReference::imag()
const
1082# ifndef PETSC_USE_COMPLEX
1083 return PetscReal(0);
1085 return PetscImaginaryPart(
static_cast<PetscScalar
>(*
this));
1092 VectorBase::in_local_range(
const size_type index)
const
1095 const PetscErrorCode ierr =
1096 VecGetOwnershipRange(
static_cast<const Vec &
>(vector), &begin, &end);
1099 return ((index >=
static_cast<size_type
>(begin)) &&
1100 (index <
static_cast<size_type
>(end)));
1105 VectorBase::locally_owned_elements()
const
1110 const std::pair<size_type, size_type> x = local_range();
1111 is.add_range(x.first, x.second);
1118 VectorBase::has_ghost_elements()
const
1125 VectorBase::ghost_elements()
const
1127 return ghost_indices;
1132 VectorBase::update_ghost_values()
const
1136 PetscErrorCode ierr;
1138 ierr = VecGhostUpdateBegin(vector, INSERT_VALUES, SCATTER_FORWARD);
1140 ierr = VecGhostUpdateEnd(vector, INSERT_VALUES, SCATTER_FORWARD);
1147 inline internal::VectorReference
1148 VectorBase::operator()(
const size_type index)
1150 return internal::VectorReference(*
this, index);
1156 VectorBase::operator()(
const size_type index)
const
1158 return static_cast<PetscScalar
>(internal::VectorReference(*
this, index));
1163 inline internal::VectorReference
1164 VectorBase::operator[](
const size_type index)
1166 return operator()(index);
1172 VectorBase::operator[](
const size_type index)
const
1174 return operator()(index);
1178 VectorBase::get_mpi_communicator()
const
1180 return PetscObjectComm(
reinterpret_cast<PetscObject
>(vector));
1184 VectorBase::extract_subvector_to(
const std::vector<size_type> &indices,
1185 std::vector<PetscScalar> &values)
const
1187 Assert(indices.size() <= values.size(),
1189 extract_subvector_to(indices.begin(), indices.end(), values.begin());
1193 VectorBase::extract_subvector_to(
1198 extract_subvector_to(indices.
begin(), indices.
end(), elements.
begin());
1202 template <
typename ForwardIterator,
typename OutputIterator>
1204 VectorBase::extract_subvector_to(
const ForwardIterator indices_begin,
1205 const ForwardIterator indices_end,
1206 OutputIterator values_begin)
const
1208 const PetscInt n_idx =
static_cast<PetscInt
>(indices_end - indices_begin);
1234 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1237 Vec locally_stored_elements =
nullptr;
1238 ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1242 ierr = VecGetSize(locally_stored_elements, &lsize);
1245 const PetscScalar *ptr;
1246 ierr = VecGetArrayRead(locally_stored_elements, &ptr);
1249 for (PetscInt i = 0; i < n_idx; ++i)
1251 const unsigned int index = *(indices_begin + i);
1252 if (index >=
static_cast<unsigned int>(begin) &&
1253 index < static_cast<unsigned int>(end))
1256 *(values_begin + i) = *(ptr + index - begin);
1261 const unsigned int ghostidx =
1262 ghost_indices.index_within_set(index);
1265 *(values_begin + i) = *(ptr + ghostidx + end - begin);
1269 ierr = VecRestoreArrayRead(locally_stored_elements, &ptr);
1272 ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1281 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1284 const PetscScalar *ptr;
1285 ierr = VecGetArrayRead(vector, &ptr);
1288 for (PetscInt i = 0; i < n_idx; ++i)
1290 const unsigned int index = *(indices_begin + i);
1292 Assert(index >=
static_cast<unsigned int>(begin) &&
1293 index <
static_cast<unsigned int>(end),
1294 ExcMessage(
"You are accessing elements of a vector without "
1295 "ghost elements that are not actually owned by "
1296 "this vector. A typical case where this may "
1297 "happen is if you are passing a non-ghosted "
1298 "(completely distributed) vector to a function "
1299 "that expects a vector that stores ghost "
1300 "elements for all locally relevant or locally "
1301 "active vector entries."));
1303 *(values_begin + i) = *(ptr + index - begin);
1306 ierr = VecRestoreArrayRead(vector, &ptr);
1311 template <
class Archive>
1313 VectorBase::save(Archive &ar,
const unsigned int)
const
1320 const PetscScalar *array =
nullptr;
1321 int ierr = VecGetArrayRead(*
this, &array);
1324 boost::serialization::array_wrapper<const PetscScalar> wrapper(
1325 array, locally_owned_size());
1328 ierr = VecRestoreArrayRead(*
this, &array);
1334 template <
class Archive>
1336 VectorBase::load(Archive &ar,
const unsigned int)
1341 std::pair<size_type, size_type> local_range;
1344 Assert(size == this->size(),
1345 ExcMessage(
"The serialized value of size (" + std::to_string(size) +
1346 ") does not match the current size (" +
1347 std::to_string(this->size()) +
")"));
1350 Assert(local_range == this->local_range(),
1351 ExcMessage(
"The serialized value of local_range (" +
1352 std::to_string(local_range.first) +
", " +
1353 std::to_string(local_range.second) +
1354 ") does not match the current local_range (" +
1355 std::to_string(this->local_range().first) +
", " +
1356 std::to_string(this->local_range().second) +
")"));
1359 PetscScalar *array =
nullptr;
1360 int ierr = VecGetArray(petsc_vector(), &array);
1363 boost::serialization::array_wrapper<PetscScalar> wrapper(
1364 array, locally_owned_size());
1367 ierr = VecRestoreArray(petsc_vector(), &array);
real_type lp_norm(const real_type p) const
real_type l1_norm() const
VectorBase & operator+=(const VectorBase &V)
PetscScalar mean_value() const
void determine_ghost_indices()
VectorOperation::values last_action
PetscScalar operator*(const VectorBase &vec) const
bool operator==(const VectorBase &v) const
VectorBase & operator*=(const PetscScalar factor)
VectorBase & operator-=(const VectorBase &V)
void save(Archive &ar, const unsigned int version) const
bool operator!=(const VectorBase &v) const
IndexSet locally_owned_elements() const
bool in_local_range(const size_type index) const
void load(Archive &ar, const unsigned int version)
std::size_t memory_consumption() const
internal::VectorReference reference
VectorBase & operator/=(const PetscScalar factor)
friend class internal::VectorReference
PetscScalar operator[](const size_type index) const
MPI_Comm get_mpi_communicator() const
void scale(const VectorBase &scaling_factors)
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< PetscScalar > &values) const
void update_ghost_values() const
std::pair< size_type, size_type > local_range() const
void sadd(const PetscScalar s, const VectorBase &V)
real_type norm_sqr() const
void compress(const VectorOperation::values operation)
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
virtual ~VectorBase() override
const IndexSet & ghost_elements() const
real_type linfty_norm() const
void swap(VectorBase &v) noexcept
reference operator()(const size_type index)
void swap(VectorBase &u, VectorBase &v) noexcept
const internal::VectorReference const_reference
void serialize(Archive &archive, const unsigned int version)
void add(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
bool has_ghost_elements() const
reference operator[](const size_type index)
virtual void extract_subvector_to(const ArrayView< const types::global_dof_index > &indices, ArrayView< PetscScalar > &elements) const override
size_type locally_owned_size() const
void equ(const PetscScalar a, const VectorBase &V)
void do_set_add_operation(const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
VectorBase & operator=(const VectorBase &)
real_type l2_norm() const
size_type size() const override
PetscScalar operator()(const size_type index) const
void extract_subvector_to(const ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcGhostsPresent()
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int global_dof_index