16#ifndef dealii_petsc_vector_base_h
17# define dealii_petsc_vector_base_h
22# ifdef DEAL_II_WITH_PETSC
40template <
typename number>
93 VectorReference(
const VectorBase &vector,
const size_type index);
99 VectorReference(
const VectorReference &vector) =
default;
112 const VectorReference &
113 operator=(
const VectorReference &r)
const;
121 operator=(
const VectorReference &r);
126 const VectorReference &
127 operator=(
const PetscScalar &s)
const;
132 const VectorReference &
133 operator+=(
const PetscScalar &s)
const;
138 const VectorReference &
139 operator-=(
const PetscScalar &s)
const;
144 const VectorReference &
145 operator*=(
const PetscScalar &s)
const;
150 const VectorReference &
151 operator/=(
const PetscScalar &s)
const;
172 operator PetscScalar()
const;
177 ExcAccessToNonlocalElement,
181 <<
"You tried to access element " << arg1
182 <<
" of a distributed vector, but only elements in range [" << arg2
183 <<
"," << arg3 <<
"] are stored locally and can be accessed."
185 <<
"A common source for this kind of problem is that you "
186 <<
"are passing a 'fully distributed' vector into a function "
187 <<
"that needs read access to vector elements that correspond "
188 <<
"to degrees of freedom on ghost cells (or at least to "
189 <<
"'locally active' degrees of freedom that are not also "
190 <<
"'locally owned'). You need to pass a vector that has these "
191 <<
"elements as ghost entries.");
198 <<
"You tried to do a "
199 << (arg1 == 1 ?
"'set'" : (arg1 == 2 ?
"'add'" :
"???"))
200 <<
" operation but the vector is currently in "
201 << (arg2 == 1 ?
"'set'" : (arg2 == 2 ?
"'add'" :
"???"))
202 <<
" mode. You first have to call 'compress()'.");
208 const VectorBase &vector;
217 friend class ::PETScWrappers::VectorBase;
386 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;
512 template <
typename ForwardIterator,
typename OutputIterator>
515 const ForwardIterator indices_end,
516 OutputIterator values_begin)
const;
523 add(
const std::vector<size_type> & indices,
524 const std::vector<PetscScalar> &
values);
531 add(
const std::vector<size_type> & indices,
532 const ::Vector<PetscScalar> &
values);
542 const PetscScalar *
values);
688 add(
const PetscScalar s);
700 add(
const PetscScalar a,
739 write_ascii(
const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
749 print(std::ostream & out,
750 const unsigned int precision = 3,
751 const bool scientific =
true,
752 const bool across =
true)
const;
776 operator const Vec &()
const;
837 const PetscScalar *
values,
838 const bool add_values);
861 inline VectorReference::VectorReference(
const VectorBase &vector,
868 inline const VectorReference &
869 VectorReference::operator=(
const VectorReference &r)
const
875 *
this =
static_cast<PetscScalar
>(r);
882 inline VectorReference &
883 VectorReference::operator=(
const VectorReference &r)
889 *
this =
static_cast<PetscScalar
>(r);
896 inline const VectorReference &
897 VectorReference::operator=(
const PetscScalar &value)
const
905 const PetscInt petsc_i = index;
907 const PetscErrorCode ierr =
908 VecSetValues(vector, 1, &petsc_i, &value, INSERT_VALUES);
918 inline const VectorReference &
919 VectorReference::operator+=(
const PetscScalar &value)
const
936 if (value == PetscScalar())
940 const PetscInt petsc_i = index;
941 const PetscErrorCode ierr =
942 VecSetValues(vector, 1, &petsc_i, &value, ADD_VALUES);
951 inline const VectorReference &
952 VectorReference::operator-=(
const PetscScalar &value)
const
969 if (value == PetscScalar())
974 const PetscInt petsc_i = index;
975 const PetscScalar subtractand = -
value;
976 const PetscErrorCode ierr =
977 VecSetValues(vector, 1, &petsc_i, &subtractand, ADD_VALUES);
985 inline const VectorReference &
986 VectorReference::operator*=(
const PetscScalar &value)
const
1006 const PetscInt petsc_i = index;
1007 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) *
value;
1009 const PetscErrorCode ierr =
1010 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1018 inline const VectorReference &
1019 VectorReference::operator/=(
const PetscScalar &value)
const
1039 const PetscInt petsc_i = index;
1040 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) /
value;
1042 const PetscErrorCode ierr =
1043 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1052 VectorReference::real()
const
1054# ifndef PETSC_USE_COMPLEX
1055 return static_cast<PetscScalar
>(*this);
1057 return PetscRealPart(
static_cast<PetscScalar
>(*
this));
1064 VectorReference::imag()
const
1066# ifndef PETSC_USE_COMPLEX
1067 return PetscReal(0);
1069 return PetscImaginaryPart(
static_cast<PetscScalar
>(*
this));
1076 VectorBase::in_local_range(
const size_type index)
const
1079 const PetscErrorCode ierr =
1080 VecGetOwnershipRange(
static_cast<const Vec &
>(vector), &
begin, &
end);
1089 VectorBase::locally_owned_elements()
const
1094 const std::pair<size_type, size_type> x = local_range();
1095 is.add_range(x.first, x.second);
1102 VectorBase::has_ghost_elements()
const
1110 VectorBase::update_ghost_values()
const
1115 inline internal::VectorReference
1116 VectorBase::operator()(
const size_type index)
1118 return internal::VectorReference(*
this, index);
1124 VectorBase::operator()(
const size_type index)
const
1126 return static_cast<PetscScalar
>(internal::VectorReference(*
this, index));
1131 inline internal::VectorReference VectorBase::operator[](
const size_type index)
1133 return operator()(index);
1138 inline PetscScalar VectorBase::operator[](
const size_type index)
const
1140 return operator()(index);
1144 VectorBase::get_mpi_communicator()
const
1147 PetscObjectGetComm(
reinterpret_cast<PetscObject
>(vector), &
comm);
1152 VectorBase::extract_subvector_to(
const std::vector<size_type> &indices,
1153 std::vector<PetscScalar> &
values)
const
1157 extract_subvector_to(indices.begin(), indices.end(),
values.begin());
1160 template <
typename ForwardIterator,
typename OutputIterator>
1162 VectorBase::extract_subvector_to(
const ForwardIterator indices_begin,
1163 const ForwardIterator indices_end,
1164 OutputIterator values_begin)
const
1166 const PetscInt n_idx =
static_cast<PetscInt
>(indices_end - indices_begin);
1192 PetscErrorCode ierr = VecGetOwnershipRange(vector, &
begin, &
end);
1195 Vec locally_stored_elements =
nullptr;
1196 ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1200 ierr = VecGetSize(locally_stored_elements, &lsize);
1204 ierr = VecGetArray(locally_stored_elements, &ptr);
1207 for (PetscInt i = 0; i < n_idx; ++i)
1209 const unsigned int index = *(indices_begin + i);
1210 if (index >=
static_cast<unsigned int>(
begin) &&
1211 index < static_cast<unsigned int>(
end))
1214 *(values_begin + i) = *(ptr + index -
begin);
1219 const unsigned int ghostidx =
1220 ghost_indices.index_within_set(index);
1223 *(values_begin + i) = *(ptr + ghostidx +
end -
begin);
1227 ierr = VecRestoreArray(locally_stored_elements, &ptr);
1230 ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1239 PetscErrorCode ierr = VecGetOwnershipRange(vector, &
begin, &
end);
1243 ierr = VecGetArray(vector, &ptr);
1246 for (PetscInt i = 0; i < n_idx; ++i)
1248 const unsigned int index = *(indices_begin + i);
1250 Assert(index >=
static_cast<unsigned int>(
begin) &&
1251 index <
static_cast<unsigned int>(
end),
1254 *(values_begin + i) = *(ptr + index -
begin);
1257 ierr = VecRestoreArray(vector, &ptr);
real_type lp_norm(const real_type p) const
real_type l1_norm() const
VectorBase & operator+=(const VectorBase &V)
PetscScalar mean_value() const
size_type local_size() const
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)
bool operator!=(const VectorBase &v) const
IndexSet locally_owned_elements() const
bool in_local_range(const size_type index) const
std::size_t memory_consumption() const
internal::VectorReference reference
VectorBase & operator/=(const PetscScalar factor)
bool is_non_negative() const
friend class internal::VectorReference
PetscScalar operator[](const size_type index) const
virtual 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
real_type linfty_norm() const
reference operator()(const size_type index)
const internal::VectorReference const_reference
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)
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)
real_type l2_norm() const
VectorBase & operator=(const VectorBase &)=delete
void swap(VectorBase &u, VectorBase &v)
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_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DEPRECATED_EARLY
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcGhostsPresent()
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)
types::global_dof_index size_type
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int global_dof_index
*braid_SplitCommworld & comm