16 #ifndef dealii_petsc_vector_base_h
17 # define dealii_petsc_vector_base_h
22 # ifdef DEAL_II_WITH_PETSC
31 # include <petscvec.h>
40 template <
typename number>
94 VectorReference(
const VectorBase &vector,
const size_type index);
109 const VectorReference &
110 operator=(
const VectorReference &r)
const;
118 operator=(
const VectorReference &r);
123 const VectorReference &
124 operator=(
const PetscScalar &s)
const;
129 const VectorReference &
130 operator+=(
const PetscScalar &s)
const;
135 const VectorReference &
136 operator-=(
const PetscScalar &s)
const;
141 const VectorReference &
142 operator*=(
const PetscScalar &s)
const;
147 const VectorReference &
148 operator/=(
const PetscScalar &s)
const;
169 operator PetscScalar()
const;
177 <<
"You tried to access element " << arg1
178 <<
" of a distributed vector, but only elements " << arg2
179 <<
" through " << arg3
180 <<
" are stored locally and can be accessed.");
187 <<
"You tried to do a "
188 << (arg1 == 1 ?
"'set'" : (arg1 == 2 ?
"'add'" :
"???"))
189 <<
" operation but the vector is currently in "
190 << (arg2 == 1 ?
"'set'" : (arg2 == 2 ?
"'add'" :
"???"))
191 <<
" mode. You first have to call 'compress()'.");
197 const VectorBase &vector;
206 friend class ::PETScWrappers::VectorBase;
362 std::pair<size_type, size_type>
439 set(
const std::vector<size_type> & indices,
440 const std::vector<PetscScalar> &values);
459 std::vector<PetscScalar> & values)
const;
488 template <
typename ForwardIterator,
typename OutputIterator>
491 const ForwardIterator indices_end,
492 OutputIterator values_begin)
const;
499 add(
const std::vector<size_type> & indices,
500 const std::vector<PetscScalar> &values);
507 add(
const std::vector<size_type> & indices,
508 const ::Vector<PetscScalar> &values);
518 const PetscScalar *values);
650 add(
const PetscScalar s);
662 add(
const PetscScalar a,
701 write_ascii(
const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
711 print(std::ostream & out,
712 const unsigned int precision = 3,
713 const bool scientific =
true,
714 const bool across =
true)
const;
738 operator const Vec &()
const;
799 const PetscScalar *values,
800 const bool add_values);
824 inline VectorReference::VectorReference(
const VectorBase &vector,
831 inline const VectorReference &
832 VectorReference::operator=(
const VectorReference &r)
const
838 *
this =
static_cast<PetscScalar
>(r);
845 inline VectorReference &
846 VectorReference::operator=(
const VectorReference &r)
852 *
this =
static_cast<PetscScalar
>(r);
859 inline const VectorReference &
860 VectorReference::operator=(
const PetscScalar &
value)
const
868 const PetscInt petsc_i = index;
870 const PetscErrorCode ierr =
871 VecSetValues(vector, 1, &petsc_i, &
value, INSERT_VALUES);
881 inline const VectorReference &
882 VectorReference::operator+=(
const PetscScalar &
value)
const
899 if (
value == PetscScalar())
903 const PetscInt petsc_i = index;
904 const PetscErrorCode ierr =
905 VecSetValues(vector, 1, &petsc_i, &
value, ADD_VALUES);
914 inline const VectorReference &
915 VectorReference::operator-=(
const PetscScalar &
value)
const
932 if (
value == PetscScalar())
937 const PetscInt petsc_i = index;
938 const PetscScalar subtractand = -
value;
939 const PetscErrorCode ierr =
940 VecSetValues(vector, 1, &petsc_i, &subtractand, ADD_VALUES);
948 inline const VectorReference &
949 VectorReference::operator*=(
const PetscScalar &
value)
const
969 const PetscInt petsc_i = index;
970 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) *
value;
972 const PetscErrorCode ierr =
973 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
981 inline const VectorReference &
982 VectorReference::operator/=(
const PetscScalar &
value)
const
1002 const PetscInt petsc_i = index;
1003 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) /
value;
1005 const PetscErrorCode ierr =
1006 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1015 VectorReference::real()
const
1017 # ifndef PETSC_USE_COMPLEX
1018 return static_cast<PetscScalar
>(*this);
1020 return PetscRealPart(
static_cast<PetscScalar
>(*
this));
1027 VectorReference::imag()
const
1029 # ifndef PETSC_USE_COMPLEX
1030 return PetscReal(0);
1032 return PetscImaginaryPart(
static_cast<PetscScalar
>(*
this));
1039 VectorBase::in_local_range(
const size_type index)
const
1042 const PetscErrorCode ierr =
1043 VecGetOwnershipRange(
static_cast<const Vec &
>(vector), &
begin, &
end);
1052 VectorBase::locally_owned_elements()
const
1057 const std::pair<size_type, size_type> x = local_range();
1058 is.add_range(x.first, x.second);
1065 VectorBase::has_ghost_elements()
const
1073 VectorBase::update_ghost_values()
const
1078 inline internal::VectorReference
1079 VectorBase::operator()(
const size_type index)
1081 return internal::VectorReference(*
this, index);
1087 VectorBase::operator()(
const size_type index)
const
1089 return static_cast<PetscScalar
>(internal::VectorReference(*
this, index));
1094 inline internal::VectorReference VectorBase::operator[](
const size_type index)
1096 return operator()(index);
1101 inline PetscScalar VectorBase::operator[](
const size_type index)
const
1103 return operator()(index);
1107 VectorBase::get_mpi_communicator()
const
1110 PetscObjectGetComm(
reinterpret_cast<PetscObject
>(vector), &comm);
1115 VectorBase::extract_subvector_to(
const std::vector<size_type> &indices,
1116 std::vector<PetscScalar> & values)
const
1118 Assert(indices.size() <= values.size(),
1120 extract_subvector_to(indices.begin(), indices.end(), values.begin());
1123 template <
typename ForwardIterator,
typename OutputIterator>
1125 VectorBase::extract_subvector_to(
const ForwardIterator indices_begin,
1126 const ForwardIterator indices_end,
1127 OutputIterator values_begin)
const
1129 const PetscInt n_idx =
static_cast<PetscInt
>(indices_end - indices_begin);
1155 PetscErrorCode ierr = VecGetOwnershipRange(vector, &
begin, &
end);
1158 Vec locally_stored_elements =
nullptr;
1159 ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1163 ierr = VecGetSize(locally_stored_elements, &lsize);
1167 ierr = VecGetArray(locally_stored_elements, &ptr);
1170 for (PetscInt i = 0; i < n_idx; ++i)
1172 const unsigned int index = *(indices_begin + i);
1173 if (index >=
static_cast<unsigned int>(
begin) &&
1174 index < static_cast<unsigned int>(
end))
1177 *(values_begin + i) = *(ptr + index -
begin);
1182 const unsigned int ghostidx =
1183 ghost_indices.index_within_set(index);
1186 *(values_begin + i) = *(ptr + ghostidx +
end -
begin);
1190 ierr = VecRestoreArray(locally_stored_elements, &ptr);
1193 ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1202 PetscErrorCode ierr = VecGetOwnershipRange(vector, &
begin, &
end);
1206 ierr = VecGetArray(vector, &ptr);
1209 for (PetscInt i = 0; i < n_idx; ++i)
1211 const unsigned int index = *(indices_begin + i);
1213 Assert(index >=
static_cast<unsigned int>(
begin) &&
1214 index <
static_cast<unsigned int>(
end),
1217 *(values_begin + i) = *(ptr + index -
begin);
1220 ierr = VecRestoreArray(vector, &ptr);
1230 # endif // DEAL_II_WITH_PETSC