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;
319 operator=(
const PetscScalar s);
362 std::pair<size_type, size_type>
370 in_local_range(
const size_type index)
const;
386 locally_owned_elements()
const;
395 has_ghost_elements()
const;
404 update_ghost_values()
const;
430 PetscScalar operator[](
const size_type index)
const;
439 set(
const std::vector<size_type> & indices,
440 const std::vector<PetscScalar> &values);
458 extract_subvector_to(
const std::vector<size_type> &indices,
459 std::vector<PetscScalar> & values)
const;
488 template <
typename ForwardIterator,
typename OutputIterator>
490 extract_subvector_to(
const ForwardIterator indices_begin,
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);
625 operator*=(
const PetscScalar factor);
631 operator/=(
const PetscScalar factor);
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;
750 virtual const MPI_Comm &
751 get_mpi_communicator()
const;
782 friend class internal::VectorReference;
797 do_set_add_operation(
const size_type n_elements,
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));
1042 const PetscErrorCode ierr =
1043 VecGetOwnershipRange(static_cast<const Vec &>(vector), &begin, &end);
1046 return ((index >= static_cast<size_type>(begin)) &&
1047 (index < static_cast<size_type>(end)));
1057 const std::pair<size_type, size_type> x = local_range();
1078 inline internal::VectorReference
1081 return internal::VectorReference(*
this, index);
1089 return static_cast<PetscScalar
>(internal::VectorReference(*
this, index));
1096 return operator()(index);
1103 return operator()(index);
1106 inline const MPI_Comm &
1109 static MPI_Comm comm;
1110 PetscObjectGetComm(reinterpret_cast<PetscObject>(vector), &comm);
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>
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 reference operator[](const size_type index)
#define DeclException2(Exception2, type1, type2, outsequence)
types::global_dof_index size_type
void update_ghost_values() const
const internal::VectorReference const_reference
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
std::vector< value_type > l2_norm(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
#define AssertIndexRange(index, range)
void swap(VectorBase &u, VectorBase &v)
__global__ void add_and_dot(Number *res, Number *v1, const Number *v2, const Number *v3, const Number a, const size_type N)
#define AssertThrow(cond, exc)
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
bool in_local_range(const size_type index) const
bool is_non_negative(const T &t)
std::string compress(const std::string &input)
VectorOperation::values last_action
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Number linfty_norm(const Tensor< 2, dim, Number > &t)
#define DEAL_II_NAMESPACE_CLOSE
VectorType::value_type * end(VectorType &V)
__global__ void equ(Number *val, const Number a, const Number *V_val, const size_type N)
__global__ void sadd(const Number s, Number *val, const Number a, const Number *V_val, const size_type N)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
static ::ExceptionBase & ExcGhostsPresent()
void add_range(const size_type begin, const size_type end)
unsigned int global_dof_index
bool has_ghost_elements() const
#define DEAL_II_NAMESPACE_OPEN
VectorType::value_type * begin(VectorType &V)
reference operator()(const size_type index)
T min(const T &t, const MPI_Comm &mpi_communicator)
Number l1_norm(const Tensor< 2, dim, Number > &t)
IndexSet locally_owned_elements() const
internal::VectorReference reference
#define DeclException3(Exception3, type1, type2, type3, outsequence)
std::enable_if< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typename ProductType< std::complex< T >, std::complex< U > >::type >::type operator*(const std::complex< T > &left, const std::complex< U > &right)
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< PetscScalar > &values) const
virtual const MPI_Comm & get_mpi_communicator() const
T max(const T &t, const MPI_Comm &mpi_communicator)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()