16 #ifndef dealii_petsc_vector_base_h 17 # define dealii_petsc_vector_base_h 20 # include <deal.II/base/config.h> 22 # ifdef DEAL_II_WITH_PETSC 24 # include <deal.II/base/index_set.h> 25 # include <deal.II/base/subscriptor.h> 27 # include <deal.II/lac/exceptions.h> 28 # include <deal.II/lac/vector.h> 29 # include <deal.II/lac/vector_operation.h> 31 # include <petscvec.h> 36 DEAL_II_NAMESPACE_OPEN
39 template <
typename number>
91 VectorReference(
const VectorBase &vector,
const size_type index);
106 const VectorReference &
107 operator=(
const VectorReference &r)
const;
115 operator=(
const VectorReference &r);
120 const VectorReference &
121 operator=(
const PetscScalar &s)
const;
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;
166 operator PetscScalar()
const;
174 <<
"You tried to access element " << arg1
175 <<
" of a distributed vector, but only elements " << arg2
176 <<
" through " << arg3
177 <<
" are stored locally and can be accessed.");
184 <<
"You tried to do a " 185 << (arg1 == 1 ?
"'set'" : (arg1 == 2 ?
"'add'" :
"???"))
186 <<
" operation but the vector is currently in " 187 << (arg2 == 1 ?
"'set'" : (arg2 == 2 ?
"'add'" :
"???"))
188 <<
" mode. You first have to call 'compress()'.");
194 const VectorBase &vector;
199 const size_type index;
205 friend class ::PETScWrappers::VectorBase;
248 using real_type = PetscReal;
361 std::pair<size_type, size_type>
429 PetscScalar
operator[](
const size_type index)
const;
438 set(
const std::vector<size_type> & indices,
439 const std::vector<PetscScalar> &values);
458 std::vector<PetscScalar> & values)
const;
487 template <
typename ForwardIterator,
typename OutputIterator>
490 const ForwardIterator indices_end,
491 OutputIterator values_begin)
const;
498 add(
const std::vector<size_type> & indices,
499 const std::vector<PetscScalar> &values);
506 add(
const std::vector<size_type> & indices,
507 const ::Vector<PetscScalar> &values);
515 add(
const size_type n_elements,
516 const size_type * indices,
517 const PetscScalar *values);
561 lp_norm(
const real_type p)
const;
649 add(
const PetscScalar s);
661 add(
const PetscScalar a,
676 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;
751 operator const Vec &()
const;
763 virtual const MPI_Comm &
813 const size_type * indices,
814 const PetscScalar *values,
815 const bool add_values);
839 inline VectorReference::VectorReference(
const VectorBase &vector,
840 const size_type index)
846 inline const VectorReference &
847 VectorReference::operator=(
const VectorReference &r)
const 853 *
this =
static_cast<PetscScalar
>(r);
860 inline VectorReference &
861 VectorReference::operator=(
const VectorReference &r)
867 *
this =
static_cast<PetscScalar
>(r);
874 inline const VectorReference &
875 VectorReference::operator=(
const PetscScalar &value)
const 883 const PetscInt petsc_i = index;
885 const PetscErrorCode ierr =
886 VecSetValues(vector, 1, &petsc_i, &value, INSERT_VALUES);
896 inline const VectorReference &
897 VectorReference::operator+=(
const PetscScalar &value)
const 914 if (value == PetscScalar())
918 const PetscInt petsc_i = index;
919 const PetscErrorCode ierr =
920 VecSetValues(vector, 1, &petsc_i, &value, ADD_VALUES);
929 inline const VectorReference &
930 VectorReference::operator-=(
const PetscScalar &value)
const 947 if (value == PetscScalar())
952 const PetscInt petsc_i = index;
953 const PetscScalar subtractand = -value;
954 const PetscErrorCode ierr =
955 VecSetValues(vector, 1, &petsc_i, &subtractand, ADD_VALUES);
963 inline const VectorReference &
964 VectorReference::operator*=(
const PetscScalar &value)
const 984 const PetscInt petsc_i = index;
985 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) * value;
987 const PetscErrorCode ierr =
988 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
996 inline const VectorReference &
997 VectorReference::operator/=(
const PetscScalar &value)
const 1017 const PetscInt petsc_i = index;
1018 const PetscScalar new_value =
static_cast<PetscScalar
>(*this) / value;
1020 const PetscErrorCode ierr =
1021 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1030 VectorReference::real()
const 1032 # ifndef PETSC_USE_COMPLEX 1033 return static_cast<PetscScalar
>(*this);
1035 return PetscRealPart(static_cast<PetscScalar>(*
this));
1042 VectorReference::imag()
const 1044 # ifndef PETSC_USE_COMPLEX 1045 return PetscReal(0);
1047 return PetscImaginaryPart(static_cast<PetscScalar>(*
this));
1054 VectorBase::in_local_range(
const size_type index)
const 1056 PetscInt begin, end;
1057 const PetscErrorCode ierr =
1058 VecGetOwnershipRange(static_cast<const Vec &>(vector), &begin, &end);
1061 return ((index >= static_cast<size_type>(begin)) &&
1062 (index < static_cast<size_type>(end)));
1067 VectorBase::locally_owned_elements()
const 1072 const std::pair<size_type, size_type> x = local_range();
1073 is.add_range(x.first, x.second);
1080 VectorBase::has_ghost_elements()
const 1088 VectorBase::update_ghost_values()
const 1093 inline internal::VectorReference
1094 VectorBase::operator()(
const size_type index)
1096 return internal::VectorReference(*
this, index);
1102 VectorBase::operator()(
const size_type index)
const 1104 return static_cast<PetscScalar
>(internal::VectorReference(*
this, index));
1109 inline internal::VectorReference VectorBase::operator[](
const size_type index)
1111 return operator()(index);
1116 inline PetscScalar VectorBase::operator[](
const size_type index)
const 1118 return operator()(index);
1121 inline const MPI_Comm &
1122 VectorBase::get_mpi_communicator()
const 1124 static MPI_Comm comm;
1125 PetscObjectGetComm(reinterpret_cast<PetscObject>(vector), &comm);
1130 VectorBase::extract_subvector_to(
const std::vector<size_type> &indices,
1131 std::vector<PetscScalar> & values)
const 1133 Assert(indices.size() <= values.size(),
1135 extract_subvector_to(indices.begin(), indices.end(), values.begin());
1138 template <
typename ForwardIterator,
typename OutputIterator>
1140 VectorBase::extract_subvector_to(
const ForwardIterator indices_begin,
1141 const ForwardIterator indices_end,
1142 OutputIterator values_begin)
const 1144 const PetscInt n_idx =
static_cast<PetscInt
>(indices_end - indices_begin);
1169 PetscInt begin, end;
1170 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1173 Vec locally_stored_elements =
nullptr;
1174 ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1178 ierr = VecGetSize(locally_stored_elements, &lsize);
1182 ierr = VecGetArray(locally_stored_elements, &ptr);
1185 for (PetscInt i = 0; i < n_idx; ++i)
1187 const unsigned int index = *(indices_begin + i);
1188 if (index >= static_cast<unsigned int>(begin) &&
1189 index < static_cast<unsigned int>(end))
1192 *(values_begin + i) = *(ptr + index - begin);
1197 const unsigned int ghostidx =
1198 ghost_indices.index_within_set(index);
1201 *(values_begin + i) = *(ptr + ghostidx + end - begin);
1205 ierr = VecRestoreArray(locally_stored_elements, &ptr);
1208 ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1216 PetscInt begin, end;
1217 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1221 ierr = VecGetArray(vector, &ptr);
1224 for (PetscInt i = 0; i < n_idx; ++i)
1226 const unsigned int index = *(indices_begin + i);
1228 Assert(index >= static_cast<unsigned int>(begin) &&
1229 index < static_cast<unsigned int>(end),
1232 *(values_begin + i) = *(ptr + index - begin);
1235 ierr = VecRestoreArray(vector, &ptr);
1243 DEAL_II_NAMESPACE_CLOSE
1245 # endif // DEAL_II_WITH_PETSC reference operator[](const size_type index)
#define DeclException2(Exception2, type1, type2, outsequence)
VectorBase & operator-=(const VectorBase &V)
void sadd(const PetscScalar s, const VectorBase &V)
real_type l2_norm() const
void update_ghost_values() const
void ratio(const VectorBase &a, const VectorBase &b)
#define AssertIndexRange(index, range)
void swap(VectorBase &u, VectorBase &v)
bool is_non_negative() const
VectorBase & operator/=(const PetscScalar factor)
#define AssertThrow(cond, exc)
real_type norm_sqr() const
bool in_local_range(const size_type index) const
void scale(const VectorBase &scaling_factors)
real_type linfty_norm() const
VectorBase & operator=(const VectorBase &)=delete
LinearAlgebra::distributed::Vector< Number > Vector
std::pair< size_type, size_type > local_range() const
real_type lp_norm(const real_type p) const
size_type local_size() const
void compress(const VectorOperation::values operation)
void do_set_add_operation(const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
VectorOperation::values last_action
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
std::size_t memory_consumption() const
void equ(const PetscScalar a, const VectorBase &V)
friend class internal::VectorReference
VectorBase & operator*=(const PetscScalar factor)
real_type l1_norm() const
static ::ExceptionBase & ExcGhostsPresent()
unsigned int global_dof_index
VectorBase & operator+=(const VectorBase &V)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
PetscScalar operator*(const VectorBase &vec) const
bool has_ghost_elements() const
PetscScalar mean_value() const
reference operator()(const size_type index)
bool operator==(const VectorBase &v) const
IndexSet locally_owned_elements() const
void add(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< PetscScalar > &values) const
virtual const MPI_Comm & get_mpi_communicator() const
virtual ~VectorBase() override
static ::ExceptionBase & ExcInternalError()
bool operator!=(const VectorBase &v) const