deal.II version GIT relicensing-2287-g6548a49e0a 2024-12-20 18:30:00+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_petsc_vector_base_h
16#define dealii_petsc_vector_base_h
17
18
19#include <deal.II/base/config.h>
20
21#ifdef DEAL_II_WITH_PETSC
22
24
26# include <deal.II/lac/vector.h>
28
29# include <boost/serialization/split_member.hpp>
30# include <boost/serialization/utility.hpp>
31
32# include <petscvec.h>
33
34# include <utility>
35# include <vector>
36
38
39// forward declaration
40# ifndef DOXYGEN
41template <typename number>
42class Vector;
43
44namespace PETScWrappers
45{
46 class VectorBase;
47}
48# endif
49
56namespace PETScWrappers
57{
67 namespace internal
68 {
81 class VectorReference
82 {
83 public:
88
89 private:
94 VectorReference(const VectorBase &vector, const size_type index);
95
96 public:
97 /*
98 * Copy constructor.
99 */
100 VectorReference(const VectorReference &vector) = default;
101
113 const VectorReference &
114 operator=(const VectorReference &r) const;
115
121 VectorReference &
122 operator=(const VectorReference &r);
123
127 const VectorReference &
128 operator=(const PetscScalar &s) const;
129
133 const VectorReference &
134 operator+=(const PetscScalar &s) const;
135
139 const VectorReference &
140 operator-=(const PetscScalar &s) const;
141
145 const VectorReference &
146 operator*=(const PetscScalar &s) const;
147
151 const VectorReference &
152 operator/=(const PetscScalar &s) const;
153
157 PetscReal
158 real() const;
159
166 PetscReal
167 imag() const;
168
173 operator PetscScalar() const;
178 ExcAccessToNonlocalElement,
179 int,
180 int,
181 int,
182 << "You tried to access element " << arg1
183 << " of a distributed vector, but only elements in range [" << arg2
184 << ',' << arg3 << "] are stored locally and can be accessed."
185 << "\n\n"
186 << "A common source for this kind of problem is that you "
187 << "are passing a 'fully distributed' vector into a function "
188 << "that needs read access to vector elements that correspond "
189 << "to degrees of freedom on ghost cells (or at least to "
190 << "'locally active' degrees of freedom that are not also "
191 << "'locally owned'). You need to pass a vector that has these "
192 << "elements as ghost entries.");
196 DeclException2(ExcWrongMode,
197 int,
198 int,
199 << "You tried to do a "
200 << (arg1 == 1 ? "'set'" : (arg1 == 2 ? "'add'" : "???"))
201 << " operation but the vector is currently in "
202 << (arg2 == 1 ? "'set'" : (arg2 == 2 ? "'add'" : "???"))
203 << " mode. You first have to call 'compress()'.");
204
205 private:
209 const VectorBase &vector;
210
214 const size_type index;
215
216 // Make the vector class a friend, so that it can create objects of the
217 // present type.
218 friend class ::PETScWrappers::VectorBase;
219 };
220 } // namespace internal
251 class VectorBase : public ReadVector<PetscScalar>
252 {
253 public:
259 using value_type = PetscScalar;
260 using real_type = PetscReal;
264
269 VectorBase();
270
275 VectorBase(const VectorBase &v);
276
281 explicit VectorBase(const Vec &v);
282
286 virtual ~VectorBase() override;
287
292 virtual void
293 clear();
294
305 void
307
311 VectorBase &
312 operator=(const VectorBase &);
313
327 VectorBase &
328 operator=(const PetscScalar s);
329
336 void
337 reinit(Vec v);
338
344 bool
345 operator==(const VectorBase &v) const;
346
352 bool
353 operator!=(const VectorBase &v) const;
354
359 size() const override;
360
370 locally_owned_size() const;
371
380 std::pair<size_type, size_type>
381 local_range() const;
382
387 bool
388 in_local_range(const size_type index) const;
389
405
412 bool
414
418 const IndexSet &
420
424 void
426
431 operator()(const size_type index);
432
436 PetscScalar
437 operator()(const size_type index) const;
438
445 operator[](const size_type index);
446
452 PetscScalar
453 operator[](const size_type index) const;
454
461 void
462 set(const std::vector<size_type> &indices,
463 const std::vector<PetscScalar> &values);
464
480 void
481 extract_subvector_to(const std::vector<size_type> &indices,
482 std::vector<PetscScalar> &values) const;
483
487 virtual void
490 const ArrayView<PetscScalar> &elements) const override;
491
519 template <typename ForwardIterator, typename OutputIterator>
520 void
524
529 void
530 add(const std::vector<size_type> &indices,
531 const std::vector<PetscScalar> &values);
532
537 void
538 add(const std::vector<size_type> &indices,
539 const ::Vector<PetscScalar> &values);
540
546 void
547 add(const size_type n_elements,
548 const size_type *indices,
549 const PetscScalar *values);
550
557 PetscScalar
558 operator*(const VectorBase &vec) const;
559
564 norm_sqr() const;
565
569 PetscScalar
570 mean_value() const;
571
580 l1_norm() const;
581
587 l2_norm() const;
588
594 lp_norm(const real_type p) const;
595
601 linfty_norm() const;
602
622 PetscScalar
623 add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W);
624
630 bool
631 all_zero() const;
632
636 VectorBase &
637 operator*=(const PetscScalar factor);
638
642 VectorBase &
643 operator/=(const PetscScalar factor);
644
648 VectorBase &
649 operator+=(const VectorBase &V);
650
654 VectorBase &
655 operator-=(const VectorBase &V);
656
661 void
662 add(const PetscScalar s);
663
667 void
668 add(const PetscScalar a, const VectorBase &V);
669
673 void
674 add(const PetscScalar a,
675 const VectorBase &V,
676 const PetscScalar b,
677 const VectorBase &W);
678
682 void
683 sadd(const PetscScalar s, const VectorBase &V);
684
688 void
689 sadd(const PetscScalar s, const PetscScalar a, const VectorBase &V);
690
696 void
698
702 void
703 equ(const PetscScalar a, const VectorBase &V);
704
712 void
714
722 void
723 print(std::ostream &out,
724 const unsigned int precision = 3,
725 const bool scientific = true,
726 const bool across = true) const;
727
735 template <class Archive>
736 void
737 save(Archive &ar, const unsigned int version) const;
738
744 template <class Archive>
745 void
746 load(Archive &ar, const unsigned int version);
747
748# ifdef DOXYGEN
754 template <class Archive>
755 void
756 serialize(Archive &archive, const unsigned int version);
757# else
758 // This macro defines the serialize() method that is compatible with
759 // the templated save() and load() method that have been implemented.
761# endif
762
775 void
776 swap(VectorBase &v) noexcept;
777
785 operator const Vec &() const;
786
792 Vec &
793 petsc_vector();
794
798 std::size_t
799 memory_consumption() const;
800
806
807 protected:
813
820
827
834
835 // Make the reference class a friend.
837
843 void
844 do_set_add_operation(const size_type n_elements,
845 const size_type *indices,
846 const PetscScalar *values,
847 const bool add_values);
848
852 void
854 };
855
856
857
858 // ------------------- inline and template functions --------------
859
867 inline void
868 swap(VectorBase &u, VectorBase &v) noexcept
869 {
870 u.swap(v);
871 }
872
873# ifndef DOXYGEN
874 namespace internal
875 {
876 inline VectorReference::VectorReference(const VectorBase &vector,
877 const size_type index)
878 : vector(vector)
879 , index(index)
880 {}
881
882
883 inline const VectorReference &
884 VectorReference::operator=(const VectorReference &r) const
885 {
886 // as explained in the class
887 // documentation, this is not the copy
888 // operator. so simply pass on to the
889 // "correct" assignment operator
890 *this = static_cast<PetscScalar>(r);
891
892 return *this;
893 }
894
895
896
897 inline VectorReference &
898 VectorReference::operator=(const VectorReference &r)
899 {
900 // as explained in the class
901 // documentation, this is not the copy
902 // operator. so simply pass on to the
903 // "correct" assignment operator
904 *this = static_cast<PetscScalar>(r);
905
906 return *this;
907 }
908
909
910
911 inline const VectorReference &
912 VectorReference::operator=(const PetscScalar &value) const
913 {
914 Assert((vector.last_action == VectorOperation::insert) ||
915 (vector.last_action == VectorOperation::unknown),
916 ExcWrongMode(VectorOperation::insert, vector.last_action));
917
918 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
919
920 const PetscInt petsc_i = index;
921
922 const PetscErrorCode ierr =
923 VecSetValues(vector, 1, &petsc_i, &value, INSERT_VALUES);
924 AssertThrow(ierr == 0, ExcPETScError(ierr));
925
926 vector.last_action = VectorOperation::insert;
927
928 return *this;
929 }
930
931
932
933 inline const VectorReference &
934 VectorReference::operator+=(const PetscScalar &value) const
935 {
936 Assert((vector.last_action == VectorOperation::add) ||
937 (vector.last_action == VectorOperation::unknown),
938 ExcWrongMode(VectorOperation::add, vector.last_action));
939
940 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
941
942 vector.last_action = VectorOperation::add;
943
944 // we have to do above actions in any
945 // case to be consistent with the MPI
946 // communication model (see the
947 // comments in the documentation of
948 // PETScWrappers::MPI::Vector), but we
949 // can save some work if the addend is
950 // zero
951 if (value == PetscScalar())
952 return *this;
953
954 // use the PETSc function to add something
955 const PetscInt petsc_i = index;
956 const PetscErrorCode ierr =
957 VecSetValues(vector, 1, &petsc_i, &value, ADD_VALUES);
958 AssertThrow(ierr == 0, ExcPETScError(ierr));
959
960
961 return *this;
962 }
963
964
965
966 inline const VectorReference &
967 VectorReference::operator-=(const PetscScalar &value) const
968 {
969 Assert((vector.last_action == VectorOperation::add) ||
970 (vector.last_action == VectorOperation::unknown),
971 ExcWrongMode(VectorOperation::add, vector.last_action));
972
973 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
974
975 vector.last_action = VectorOperation::add;
976
977 // we have to do above actions in any
978 // case to be consistent with the MPI
979 // communication model (see the
980 // comments in the documentation of
981 // PETScWrappers::MPI::Vector), but we
982 // can save some work if the addend is
983 // zero
984 if (value == PetscScalar())
985 return *this;
986
987 // use the PETSc function to
988 // add something
989 const PetscInt petsc_i = index;
990 const PetscScalar subtractand = -value;
991 const PetscErrorCode ierr =
992 VecSetValues(vector, 1, &petsc_i, &subtractand, ADD_VALUES);
993 AssertThrow(ierr == 0, ExcPETScError(ierr));
994
995 return *this;
996 }
997
998
999
1000 inline const VectorReference &
1001 VectorReference::operator*=(const PetscScalar &value) const
1002 {
1003 Assert((vector.last_action == VectorOperation::insert) ||
1004 (vector.last_action == VectorOperation::unknown),
1005 ExcWrongMode(VectorOperation::insert, vector.last_action));
1006
1007 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
1008
1009 vector.last_action = VectorOperation::insert;
1010
1011 // we have to do above actions in any
1012 // case to be consistent with the MPI
1013 // communication model (see the
1014 // comments in the documentation of
1015 // PETScWrappers::MPI::Vector), but we
1016 // can save some work if the factor is
1017 // one
1018 if (value == 1.)
1019 return *this;
1020
1021 const PetscInt petsc_i = index;
1022 const PetscScalar new_value = static_cast<PetscScalar>(*this) * value;
1023
1024 const PetscErrorCode ierr =
1025 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1026 AssertThrow(ierr == 0, ExcPETScError(ierr));
1027
1028 return *this;
1029 }
1030
1031
1032
1033 inline const VectorReference &
1034 VectorReference::operator/=(const PetscScalar &value) const
1035 {
1036 Assert((vector.last_action == VectorOperation::insert) ||
1037 (vector.last_action == VectorOperation::unknown),
1038 ExcWrongMode(VectorOperation::insert, vector.last_action));
1039
1040 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
1041
1042 vector.last_action = VectorOperation::insert;
1043
1044 // we have to do above actions in any
1045 // case to be consistent with the MPI
1046 // communication model (see the
1047 // comments in the documentation of
1048 // PETScWrappers::MPI::Vector), but we
1049 // can save some work if the factor is
1050 // one
1051 if (value == 1.)
1052 return *this;
1053
1054 const PetscInt petsc_i = index;
1055 const PetscScalar new_value = static_cast<PetscScalar>(*this) / value;
1056
1057 const PetscErrorCode ierr =
1058 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1059 AssertThrow(ierr == 0, ExcPETScError(ierr));
1060
1061 return *this;
1062 }
1063
1064
1065
1066 inline PetscReal
1067 VectorReference::real() const
1068 {
1069# ifndef PETSC_USE_COMPLEX
1070 return static_cast<PetscScalar>(*this);
1071# else
1072 return PetscRealPart(static_cast<PetscScalar>(*this));
1073# endif
1074 }
1075
1076
1077
1078 inline PetscReal
1079 VectorReference::imag() const
1080 {
1081# ifndef PETSC_USE_COMPLEX
1082 return PetscReal(0);
1083# else
1084 return PetscImaginaryPart(static_cast<PetscScalar>(*this));
1085# endif
1086 }
1087
1088 } // namespace internal
1089
1090 inline bool
1091 VectorBase::in_local_range(const size_type index) const
1092 {
1093 PetscInt begin, end;
1094 const PetscErrorCode ierr =
1095 VecGetOwnershipRange(static_cast<const Vec &>(vector), &begin, &end);
1096 AssertThrow(ierr == 0, ExcPETScError(ierr));
1097
1098 return ((index >= static_cast<size_type>(begin)) &&
1099 (index < static_cast<size_type>(end)));
1100 }
1101
1102
1103 inline IndexSet
1104 VectorBase::locally_owned_elements() const
1105 {
1106 IndexSet is(size());
1107
1108 // PETSc only allows for contiguous local ranges, so this is simple
1109 const std::pair<size_type, size_type> x = local_range();
1110 is.add_range(x.first, x.second);
1111 return is;
1112 }
1113
1114
1115
1116 inline bool
1117 VectorBase::has_ghost_elements() const
1118 {
1119 return ghosted;
1120 }
1121
1122
1123 inline const IndexSet &
1124 VectorBase::ghost_elements() const
1125 {
1126 return ghost_indices;
1127 }
1128
1129
1130 inline void
1131 VectorBase::update_ghost_values() const
1132 {
1133 if (ghosted)
1134 {
1135 PetscErrorCode ierr;
1136
1137 ierr = VecGhostUpdateBegin(vector, INSERT_VALUES, SCATTER_FORWARD);
1138 AssertThrow(ierr == 0, ExcPETScError(ierr));
1139 ierr = VecGhostUpdateEnd(vector, INSERT_VALUES, SCATTER_FORWARD);
1140 AssertThrow(ierr == 0, ExcPETScError(ierr));
1141 }
1142 }
1143
1144
1145
1146 inline internal::VectorReference
1147 VectorBase::operator()(const size_type index)
1148 {
1149 return internal::VectorReference(*this, index);
1150 }
1151
1152
1153
1154 inline PetscScalar
1155 VectorBase::operator()(const size_type index) const
1156 {
1157 return static_cast<PetscScalar>(internal::VectorReference(*this, index));
1158 }
1159
1160
1161
1162 inline internal::VectorReference
1163 VectorBase::operator[](const size_type index)
1164 {
1165 return operator()(index);
1166 }
1167
1168
1169
1170 inline PetscScalar
1171 VectorBase::operator[](const size_type index) const
1172 {
1173 return operator()(index);
1174 }
1175
1176 inline MPI_Comm
1177 VectorBase::get_mpi_communicator() const
1178 {
1179 return PetscObjectComm(reinterpret_cast<PetscObject>(vector));
1180 }
1181
1182 inline void
1183 VectorBase::extract_subvector_to(const std::vector<size_type> &indices,
1184 std::vector<PetscScalar> &values) const
1185 {
1186 Assert(indices.size() <= values.size(),
1187 ExcDimensionMismatch(indices.size(), values.size()));
1188 extract_subvector_to(indices.begin(), indices.end(), values.begin());
1189 }
1190
1191 inline void
1192 VectorBase::extract_subvector_to(
1194 const ArrayView<PetscScalar> &elements) const
1195 {
1196 AssertDimension(indices.size(), elements.size());
1197 extract_subvector_to(indices.begin(), indices.end(), elements.begin());
1198 }
1199
1200
1201 template <typename ForwardIterator, typename OutputIterator>
1202 inline void
1203 VectorBase::extract_subvector_to(const ForwardIterator indices_begin,
1204 const ForwardIterator indices_end,
1205 OutputIterator values_begin) const
1206 {
1207 if (indices_begin == indices_end)
1208 return;
1209
1210 // if we are dealing
1211 // with a parallel vector
1212 if (ghosted)
1213 {
1214 // there is the possibility
1215 // that the vector has
1216 // ghost elements. in that
1217 // case, we first need to
1218 // figure out which
1219 // elements we own locally,
1220 // then get a pointer to
1221 // the elements that are
1222 // stored here (both the
1223 // ones we own as well as
1224 // the ghost elements). in
1225 // this array, the locally
1226 // owned elements come
1227 // first followed by the
1228 // ghost elements whose
1229 // position we can get from
1230 // an index set
1231 PetscInt begin, end;
1232 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1233 AssertThrow(ierr == 0, ExcPETScError(ierr));
1234
1235 Vec locally_stored_elements = nullptr;
1236 ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1237 AssertThrow(ierr == 0, ExcPETScError(ierr));
1238
1239 PetscInt lsize;
1240 ierr = VecGetSize(locally_stored_elements, &lsize);
1241 AssertThrow(ierr == 0, ExcPETScError(ierr));
1242
1243 const PetscScalar *ptr;
1244 ierr = VecGetArrayRead(locally_stored_elements, &ptr);
1245 AssertThrow(ierr == 0, ExcPETScError(ierr));
1246
1247 auto input = indices_begin;
1248 auto output = values_begin;
1249 while (input != indices_end)
1250 {
1251 const auto index = static_cast<PetscInt>(*input);
1252 AssertIntegerConversion(index, *input);
1253 if (index >= begin && index < end)
1254 {
1255 // local entry
1256 *output = *(ptr + index - begin);
1257 }
1258 else
1259 {
1260 // ghost entry
1261 const auto ghost_index = ghost_indices.index_within_set(*input);
1262
1263 AssertIndexRange(ghost_index + end - begin, lsize);
1264 *output = *(ptr + ghost_index + end - begin);
1265 }
1266
1267 ++input;
1268 ++output;
1269 }
1270
1271 ierr = VecRestoreArrayRead(locally_stored_elements, &ptr);
1272 AssertThrow(ierr == 0, ExcPETScError(ierr));
1273
1274 ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1275 AssertThrow(ierr == 0, ExcPETScError(ierr));
1276 }
1277 // if the vector is local or the
1278 // caller, then simply access the
1279 // element we are interested in
1280 else
1281 {
1282 PetscInt begin, end;
1283 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1284 AssertThrow(ierr == 0, ExcPETScError(ierr));
1285
1286 const PetscScalar *ptr;
1287 ierr = VecGetArrayRead(vector, &ptr);
1288 AssertThrow(ierr == 0, ExcPETScError(ierr));
1289
1290 auto input = indices_begin;
1291 auto output = values_begin;
1292 while (input != indices_end)
1293 {
1294 const auto index = static_cast<PetscInt>(*input);
1295 AssertIntegerConversion(index, *input);
1296
1297 Assert(index >= begin && index < end,
1298 ExcMessage("You are accessing elements of a vector without "
1299 "ghost elements that are not actually owned by "
1300 "this vector. A typical case where this may "
1301 "happen is if you are passing a non-ghosted "
1302 "(completely distributed) vector to a function "
1303 "that expects a vector that stores ghost "
1304 "elements for all locally relevant or locally "
1305 "active vector entries."));
1306
1307 *output = *(ptr + index - begin);
1308
1309 ++input;
1310 ++output;
1311 }
1312
1313 ierr = VecRestoreArrayRead(vector, &ptr);
1314 AssertThrow(ierr == 0, ExcPETScError(ierr));
1315 }
1316 }
1317
1318 template <class Archive>
1319 inline void
1320 VectorBase::save(Archive &ar, const unsigned int) const
1321 {
1322 // forward to serialization function in the base class.
1323 ar &static_cast<const EnableObserverPointer &>(*this);
1324 ar &size();
1325 ar &local_range();
1326
1327 const PetscScalar *array = nullptr;
1328 int ierr = VecGetArrayRead(*this, &array);
1329 AssertThrow(ierr == 0, ExcPETScError(ierr));
1330
1331 boost::serialization::array_wrapper<const PetscScalar> wrapper(
1332 array, locally_owned_size());
1333 ar &wrapper;
1334
1335 ierr = VecRestoreArrayRead(*this, &array);
1336 AssertThrow(ierr == 0, ExcPETScError(ierr));
1337 }
1338
1339
1340
1341 template <class Archive>
1342 inline void
1343 VectorBase::load(Archive &ar, const unsigned int)
1344 {
1345 ar &static_cast<EnableObserverPointer &>(*this);
1346
1347 size_type size = 0;
1348 std::pair<size_type, size_type> local_range;
1349
1350 ar &size;
1351 Assert(size == this->size(),
1352 ExcMessage("The serialized value of size (" + std::to_string(size) +
1353 ") does not match the current size (" +
1354 std::to_string(this->size()) + ")"));
1355 (void)size;
1356 ar &local_range;
1357 Assert(local_range == this->local_range(),
1358 ExcMessage("The serialized value of local_range (" +
1359 std::to_string(local_range.first) + ", " +
1360 std::to_string(local_range.second) +
1361 ") does not match the current local_range (" +
1362 std::to_string(this->local_range().first) + ", " +
1363 std::to_string(this->local_range().second) + ")"));
1364 (void)local_range;
1365
1366 PetscScalar *array = nullptr;
1367 int ierr = VecGetArray(petsc_vector(), &array);
1368 AssertThrow(ierr == 0, ExcPETScError(ierr));
1369
1370 boost::serialization::array_wrapper<PetscScalar> wrapper(
1371 array, locally_owned_size());
1372 ar &wrapper;
1373
1374 ierr = VecRestoreArray(petsc_vector(), &array);
1375 AssertThrow(ierr == 0, ExcPETScError(ierr));
1376 }
1377# endif // DOXYGEN
1378} // namespace PETScWrappers
1379
1381
1382#endif // DEAL_II_WITH_PETSC
1383
1384#endif
iterator begin() const
Definition array_view.h:702
iterator end() const
Definition array_view.h:711
std::size_t size() const
Definition array_view.h:684
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)
void save(Archive &ar, const unsigned int version) const
bool operator!=(const VectorBase &v) const
IndexSet locally_owned_elements() const
bool in_local_range(const size_type index) const
void load(Archive &ar, const unsigned int version)
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
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)
const IndexSet & ghost_elements() const
void swap(VectorBase &v) noexcept
reference operator()(const size_type index)
void swap(VectorBase &u, VectorBase &v) noexcept
virtual void extract_subvector_to(const ArrayView< const types::global_dof_index > &indices, const ArrayView< PetscScalar > &elements) const override
const internal::VectorReference const_reference
void serialize(Archive &archive, const unsigned int version)
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 &)
size_type size() const override
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_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
static ::ExceptionBase & ExcGhostsPresent()
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:534
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition exceptions.h:557
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
#define AssertIntegerConversion(index1, index2)
std::pair< types::global_dof_index, types::global_dof_index > local_range
Definition mpi.cc:815
std::size_t size
Definition mpi.cc:734
types::global_dof_index locally_owned_size
Definition mpi.cc:822
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int global_dof_index
Definition types.h:81