Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
petsc_vector_base.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2004 - 2022 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_petsc_vector_base_h
17# define dealii_petsc_vector_base_h
18
19
20# include <deal.II/base/config.h>
21
22# ifdef DEAL_II_WITH_PETSC
23
26
28# include <deal.II/lac/vector.h>
30
31# include <petscvec.h>
32
33# include <utility>
34# include <vector>
35
37
38// forward declaration
39# ifndef DOXYGEN
40template <typename number>
41class Vector;
42
43namespace PETScWrappers
44{
45 class VectorBase;
46}
47# endif
48
55namespace PETScWrappers
56{
66 namespace internal
67 {
80 class VectorReference
81 {
82 public:
87
88 private:
93 VectorReference(const VectorBase &vector, const size_type index);
94
95 public:
96 /*
97 * Copy constructor.
98 */
99 VectorReference(const VectorReference &vector) = default;
100
112 const VectorReference &
113 operator=(const VectorReference &r) const;
114
120 VectorReference &
121 operator=(const VectorReference &r);
122
126 const VectorReference &
127 operator=(const PetscScalar &s) const;
128
132 const VectorReference &
133 operator+=(const PetscScalar &s) const;
134
138 const VectorReference &
139 operator-=(const PetscScalar &s) const;
140
144 const VectorReference &
145 operator*=(const PetscScalar &s) const;
146
150 const VectorReference &
151 operator/=(const PetscScalar &s) const;
152
156 PetscReal
157 real() const;
158
165 PetscReal
166 imag() const;
167
172 operator PetscScalar() const;
177 ExcAccessToNonlocalElement,
178 int,
179 int,
180 int,
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."
184 << "\n\n"
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.");
195 DeclException2(ExcWrongMode,
196 int,
197 int,
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()'.");
203
204 private:
208 const VectorBase &vector;
209
213 const size_type index;
214
215 // Make the vector class a friend, so that it can create objects of the
216 // present type.
217 friend class ::PETScWrappers::VectorBase;
218 };
219 } // namespace internal
250 class VectorBase : public Subscriptor
251 {
252 public:
258 using value_type = PetscScalar;
259 using real_type = PetscReal;
263
268 VectorBase();
269
274 VectorBase(const VectorBase &v);
275
281 explicit VectorBase(const Vec &v);
282
287 VectorBase &
288 operator=(const VectorBase &) = delete;
289
293 virtual ~VectorBase() override;
294
299 virtual void
300 clear();
301
312 void
313 compress(const VectorOperation::values operation);
314
328 VectorBase &
329 operator=(const PetscScalar s);
330
336 bool
337 operator==(const VectorBase &v) const;
338
344 bool
345 operator!=(const VectorBase &v) const;
346
351 size() const;
352
365 local_size() const;
366
376 locally_owned_size() const;
377
386 std::pair<size_type, size_type>
387 local_range() const;
388
393 bool
394 in_local_range(const size_type index) const;
395
411
418 bool
420
427 void
429
434 operator()(const size_type index);
435
439 PetscScalar
440 operator()(const size_type index) const;
441
448 operator[](const size_type index);
449
455 PetscScalar
456 operator[](const size_type index) const;
457
464 void
465 set(const std::vector<size_type> & indices,
466 const std::vector<PetscScalar> &values);
467
483 void
484 extract_subvector_to(const std::vector<size_type> &indices,
485 std::vector<PetscScalar> & values) const;
486
514 template <typename ForwardIterator, typename OutputIterator>
515 void
516 extract_subvector_to(const ForwardIterator indices_begin,
517 const ForwardIterator indices_end,
518 OutputIterator values_begin) const;
519
524 void
525 add(const std::vector<size_type> & indices,
526 const std::vector<PetscScalar> &values);
527
532 void
533 add(const std::vector<size_type> & indices,
534 const ::Vector<PetscScalar> &values);
535
541 void
542 add(const size_type n_elements,
543 const size_type * indices,
544 const PetscScalar *values);
545
552 PetscScalar
553 operator*(const VectorBase &vec) const;
554
559 norm_sqr() const;
560
564 PetscScalar
565 mean_value() const;
566
575 l1_norm() const;
576
582 l2_norm() const;
583
589 lp_norm(const real_type p) const;
590
596 linfty_norm() const;
597
617 PetscScalar
618 add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W);
619
629 min() const;
630
640 max() const;
641
647 bool
648 all_zero() const;
649
659 bool
660 is_non_negative() const;
661
665 VectorBase &
666 operator*=(const PetscScalar factor);
667
671 VectorBase &
672 operator/=(const PetscScalar factor);
673
677 VectorBase &
678 operator+=(const VectorBase &V);
679
683 VectorBase &
684 operator-=(const VectorBase &V);
685
690 void
691 add(const PetscScalar s);
692
696 void
697 add(const PetscScalar a, const VectorBase &V);
698
702 void
703 add(const PetscScalar a,
704 const VectorBase &V,
705 const PetscScalar b,
706 const VectorBase &W);
707
711 void
712 sadd(const PetscScalar s, const VectorBase &V);
713
717 void
718 sadd(const PetscScalar s, const PetscScalar a, const VectorBase &V);
719
725 void
726 scale(const VectorBase &scaling_factors);
727
731 void
732 equ(const PetscScalar a, const VectorBase &V);
733
741 void
742 write_ascii(const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
743
751 void
752 print(std::ostream & out,
753 const unsigned int precision = 3,
754 const bool scientific = true,
755 const bool across = true) const;
756
769 void
770 swap(VectorBase &v);
771
779 operator const Vec &() const;
780
784 std::size_t
785 memory_consumption() const;
786
791 virtual const MPI_Comm &
793
794 protected:
800
807
814
821
822 // Make the reference class a friend.
824
831
837 void
838 do_set_add_operation(const size_type n_elements,
839 const size_type * indices,
840 const PetscScalar *values,
841 const bool add_values);
842 };
843
844
845
846 // ------------------- inline and template functions --------------
847
855 inline void
857 {
858 u.swap(v);
859 }
860
861# ifndef DOXYGEN
862 namespace internal
863 {
864 inline VectorReference::VectorReference(const VectorBase &vector,
865 const size_type index)
866 : vector(vector)
867 , index(index)
868 {}
869
870
871 inline const VectorReference &
872 VectorReference::operator=(const VectorReference &r) const
873 {
874 // as explained in the class
875 // documentation, this is not the copy
876 // operator. so simply pass on to the
877 // "correct" assignment operator
878 *this = static_cast<PetscScalar>(r);
879
880 return *this;
881 }
882
883
884
885 inline VectorReference &
886 VectorReference::operator=(const VectorReference &r)
887 {
888 // as explained in the class
889 // documentation, this is not the copy
890 // operator. so simply pass on to the
891 // "correct" assignment operator
892 *this = static_cast<PetscScalar>(r);
893
894 return *this;
895 }
896
897
898
899 inline const VectorReference &
900 VectorReference::operator=(const PetscScalar &value) const
901 {
902 Assert((vector.last_action == VectorOperation::insert) ||
903 (vector.last_action == VectorOperation::unknown),
904 ExcWrongMode(VectorOperation::insert, vector.last_action));
905
906 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
907
908 const PetscInt petsc_i = index;
909
910 const PetscErrorCode ierr =
911 VecSetValues(vector, 1, &petsc_i, &value, INSERT_VALUES);
912 AssertThrow(ierr == 0, ExcPETScError(ierr));
913
914 vector.last_action = VectorOperation::insert;
915
916 return *this;
917 }
918
919
920
921 inline const VectorReference &
922 VectorReference::operator+=(const PetscScalar &value) const
923 {
924 Assert((vector.last_action == VectorOperation::add) ||
925 (vector.last_action == VectorOperation::unknown),
926 ExcWrongMode(VectorOperation::add, vector.last_action));
927
928 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
929
930 vector.last_action = VectorOperation::add;
931
932 // we have to do above actions in any
933 // case to be consistent with the MPI
934 // communication model (see the
935 // comments in the documentation of
936 // PETScWrappers::MPI::Vector), but we
937 // can save some work if the addend is
938 // zero
939 if (value == PetscScalar())
940 return *this;
941
942 // use the PETSc function to add something
943 const PetscInt petsc_i = index;
944 const PetscErrorCode ierr =
945 VecSetValues(vector, 1, &petsc_i, &value, ADD_VALUES);
946 AssertThrow(ierr == 0, ExcPETScError(ierr));
947
948
949 return *this;
950 }
951
952
953
954 inline const VectorReference &
955 VectorReference::operator-=(const PetscScalar &value) const
956 {
957 Assert((vector.last_action == VectorOperation::add) ||
958 (vector.last_action == VectorOperation::unknown),
959 ExcWrongMode(VectorOperation::add, vector.last_action));
960
961 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
962
963 vector.last_action = VectorOperation::add;
964
965 // we have to do above actions in any
966 // case to be consistent with the MPI
967 // communication model (see the
968 // comments in the documentation of
969 // PETScWrappers::MPI::Vector), but we
970 // can save some work if the addend is
971 // zero
972 if (value == PetscScalar())
973 return *this;
974
975 // use the PETSc function to
976 // add something
977 const PetscInt petsc_i = index;
978 const PetscScalar subtractand = -value;
979 const PetscErrorCode ierr =
980 VecSetValues(vector, 1, &petsc_i, &subtractand, ADD_VALUES);
981 AssertThrow(ierr == 0, ExcPETScError(ierr));
982
983 return *this;
984 }
985
986
987
988 inline const VectorReference &
989 VectorReference::operator*=(const PetscScalar &value) const
990 {
991 Assert((vector.last_action == VectorOperation::insert) ||
992 (vector.last_action == VectorOperation::unknown),
993 ExcWrongMode(VectorOperation::insert, vector.last_action));
994
995 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
996
997 vector.last_action = VectorOperation::insert;
998
999 // we have to do above actions in any
1000 // case to be consistent with the MPI
1001 // communication model (see the
1002 // comments in the documentation of
1003 // PETScWrappers::MPI::Vector), but we
1004 // can save some work if the factor is
1005 // one
1006 if (value == 1.)
1007 return *this;
1008
1009 const PetscInt petsc_i = index;
1010 const PetscScalar new_value = static_cast<PetscScalar>(*this) * value;
1011
1012 const PetscErrorCode ierr =
1013 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1014 AssertThrow(ierr == 0, ExcPETScError(ierr));
1015
1016 return *this;
1017 }
1018
1019
1020
1021 inline const VectorReference &
1022 VectorReference::operator/=(const PetscScalar &value) const
1023 {
1024 Assert((vector.last_action == VectorOperation::insert) ||
1025 (vector.last_action == VectorOperation::unknown),
1026 ExcWrongMode(VectorOperation::insert, vector.last_action));
1027
1028 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
1029
1030 vector.last_action = VectorOperation::insert;
1031
1032 // we have to do above actions in any
1033 // case to be consistent with the MPI
1034 // communication model (see the
1035 // comments in the documentation of
1036 // PETScWrappers::MPI::Vector), but we
1037 // can save some work if the factor is
1038 // one
1039 if (value == 1.)
1040 return *this;
1041
1042 const PetscInt petsc_i = index;
1043 const PetscScalar new_value = static_cast<PetscScalar>(*this) / value;
1044
1045 const PetscErrorCode ierr =
1046 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1047 AssertThrow(ierr == 0, ExcPETScError(ierr));
1048
1049 return *this;
1050 }
1051
1052
1053
1054 inline PetscReal
1055 VectorReference::real() const
1056 {
1057# ifndef PETSC_USE_COMPLEX
1058 return static_cast<PetscScalar>(*this);
1059# else
1060 return PetscRealPart(static_cast<PetscScalar>(*this));
1061# endif
1062 }
1063
1064
1065
1066 inline PetscReal
1067 VectorReference::imag() const
1068 {
1069# ifndef PETSC_USE_COMPLEX
1070 return PetscReal(0);
1071# else
1072 return PetscImaginaryPart(static_cast<PetscScalar>(*this));
1073# endif
1074 }
1075
1076 } // namespace internal
1077
1078 inline bool
1079 VectorBase::in_local_range(const size_type index) const
1080 {
1081 PetscInt begin, end;
1082 const PetscErrorCode ierr =
1083 VecGetOwnershipRange(static_cast<const Vec &>(vector), &begin, &end);
1084 AssertThrow(ierr == 0, ExcPETScError(ierr));
1085
1086 return ((index >= static_cast<size_type>(begin)) &&
1087 (index < static_cast<size_type>(end)));
1088 }
1089
1090
1091 inline IndexSet
1092 VectorBase::locally_owned_elements() const
1093 {
1094 IndexSet is(size());
1095
1096 // PETSc only allows for contiguous local ranges, so this is simple
1097 const std::pair<size_type, size_type> x = local_range();
1098 is.add_range(x.first, x.second);
1099 return is;
1100 }
1101
1102
1103
1104 inline bool
1105 VectorBase::has_ghost_elements() const
1106 {
1107 return ghosted;
1108 }
1109
1110
1111
1112 inline void
1113 VectorBase::update_ghost_values() const
1114 {}
1115
1116
1117
1118 inline internal::VectorReference
1119 VectorBase::operator()(const size_type index)
1120 {
1121 return internal::VectorReference(*this, index);
1122 }
1123
1124
1125
1126 inline PetscScalar
1127 VectorBase::operator()(const size_type index) const
1128 {
1129 return static_cast<PetscScalar>(internal::VectorReference(*this, index));
1130 }
1131
1132
1133
1134 inline internal::VectorReference
1135 VectorBase::operator[](const size_type index)
1136 {
1137 return operator()(index);
1138 }
1139
1140
1141
1142 inline PetscScalar
1143 VectorBase::operator[](const size_type index) const
1144 {
1145 return operator()(index);
1146 }
1147
1148 inline const MPI_Comm &
1149 VectorBase::get_mpi_communicator() const
1150 {
1151 static MPI_Comm comm;
1152 PetscObjectGetComm(reinterpret_cast<PetscObject>(vector), &comm);
1153 return comm;
1154 }
1155
1156 inline void
1157 VectorBase::extract_subvector_to(const std::vector<size_type> &indices,
1158 std::vector<PetscScalar> & values) const
1159 {
1160 Assert(indices.size() <= values.size(),
1161 ExcDimensionMismatch(indices.size(), values.size()));
1162 extract_subvector_to(indices.begin(), indices.end(), values.begin());
1163 }
1164
1165 template <typename ForwardIterator, typename OutputIterator>
1166 inline void
1167 VectorBase::extract_subvector_to(const ForwardIterator indices_begin,
1168 const ForwardIterator indices_end,
1169 OutputIterator values_begin) const
1170 {
1171 const PetscInt n_idx = static_cast<PetscInt>(indices_end - indices_begin);
1172 if (n_idx == 0)
1173 return;
1174
1175 // if we are dealing
1176 // with a parallel vector
1177 if (ghosted)
1178 {
1179 // there is the possibility
1180 // that the vector has
1181 // ghost elements. in that
1182 // case, we first need to
1183 // figure out which
1184 // elements we own locally,
1185 // then get a pointer to
1186 // the elements that are
1187 // stored here (both the
1188 // ones we own as well as
1189 // the ghost elements). in
1190 // this array, the locally
1191 // owned elements come
1192 // first followed by the
1193 // ghost elements whose
1194 // position we can get from
1195 // an index set
1196 PetscInt begin, end;
1197 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1198 AssertThrow(ierr == 0, ExcPETScError(ierr));
1199
1200 Vec locally_stored_elements = nullptr;
1201 ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1202 AssertThrow(ierr == 0, ExcPETScError(ierr));
1203
1204 PetscInt lsize;
1205 ierr = VecGetSize(locally_stored_elements, &lsize);
1206 AssertThrow(ierr == 0, ExcPETScError(ierr));
1207
1208 PetscScalar *ptr;
1209 ierr = VecGetArray(locally_stored_elements, &ptr);
1210 AssertThrow(ierr == 0, ExcPETScError(ierr));
1211
1212 for (PetscInt i = 0; i < n_idx; ++i)
1213 {
1214 const unsigned int index = *(indices_begin + i);
1215 if (index >= static_cast<unsigned int>(begin) &&
1216 index < static_cast<unsigned int>(end))
1217 {
1218 // local entry
1219 *(values_begin + i) = *(ptr + index - begin);
1220 }
1221 else
1222 {
1223 // ghost entry
1224 const unsigned int ghostidx =
1225 ghost_indices.index_within_set(index);
1226
1227 AssertIndexRange(ghostidx + end - begin, lsize);
1228 *(values_begin + i) = *(ptr + ghostidx + end - begin);
1229 }
1230 }
1231
1232 ierr = VecRestoreArray(locally_stored_elements, &ptr);
1233 AssertThrow(ierr == 0, ExcPETScError(ierr));
1234
1235 ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1236 AssertThrow(ierr == 0, ExcPETScError(ierr));
1237 }
1238 // if the vector is local or the
1239 // caller, then simply access the
1240 // element we are interested in
1241 else
1242 {
1243 PetscInt begin, end;
1244 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1245 AssertThrow(ierr == 0, ExcPETScError(ierr));
1246
1247 PetscScalar *ptr;
1248 ierr = VecGetArray(vector, &ptr);
1249 AssertThrow(ierr == 0, ExcPETScError(ierr));
1250
1251 for (PetscInt i = 0; i < n_idx; ++i)
1252 {
1253 const unsigned int index = *(indices_begin + i);
1254
1255 Assert(index >= static_cast<unsigned int>(begin) &&
1256 index < static_cast<unsigned int>(end),
1257 ExcMessage("You are accessing elements of a vector without "
1258 "ghost elements that are not actually owned by "
1259 "this vector. A typical case where this may "
1260 "happen is if you are passing a non-ghosted "
1261 "(completely distributed) vector to a function "
1262 "that expects a vector that stores ghost "
1263 "elements for all locally relevant or locally "
1264 "active vector entries."));
1265
1266 *(values_begin + i) = *(ptr + index - begin);
1267 }
1268
1269 ierr = VecRestoreArray(vector, &ptr);
1270 AssertThrow(ierr == 0, ExcPETScError(ierr));
1271 }
1272 }
1273
1274# endif // DOXYGEN
1275} // namespace PETScWrappers
1276
1278
1279# endif // DEAL_II_WITH_PETSC
1280
1281#endif
1282/*---------------------------- petsc_vector_base.h --------------------------*/
real_type lp_norm(const real_type p) const
VectorBase & operator+=(const VectorBase &V)
PetscScalar mean_value() 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)
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)
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
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)
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
Definition: vector.h:109
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & ExcGhostsPresent()
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:532
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:555
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1583
types::global_dof_index size_type
Definition: cuda_kernels.h:45
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int global_dof_index
Definition: types.h:76
const MPI_Comm & comm