Reference documentation for deal.II version 9.6.0
\(\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
25
27# include <deal.II/lac/vector.h>
29
30# include <boost/serialization/split_member.hpp>
31# include <boost/serialization/utility.hpp>
32
33# include <petscvec.h>
34
35# include <utility>
36# include <vector>
37
39
40// forward declaration
41# ifndef DOXYGEN
42template <typename number>
43class Vector;
44
45namespace PETScWrappers
46{
47 class VectorBase;
48}
49# endif
50
57namespace PETScWrappers
58{
68 namespace internal
69 {
82 class VectorReference
83 {
84 public:
88 using size_type = types::global_dof_index;
89
90 private:
95 VectorReference(const VectorBase &vector, const size_type index);
96
97 public:
98 /*
99 * Copy constructor.
100 */
101 VectorReference(const VectorReference &vector) = default;
102
114 const VectorReference &
115 operator=(const VectorReference &r) const;
116
122 VectorReference &
123 operator=(const VectorReference &r);
124
128 const VectorReference &
129 operator=(const PetscScalar &s) const;
130
134 const VectorReference &
135 operator+=(const PetscScalar &s) const;
136
140 const VectorReference &
141 operator-=(const PetscScalar &s) const;
142
146 const VectorReference &
147 operator*=(const PetscScalar &s) const;
148
152 const VectorReference &
153 operator/=(const PetscScalar &s) const;
154
158 PetscReal
159 real() const;
160
167 PetscReal
168 imag() const;
169
174 operator PetscScalar() const;
179 ExcAccessToNonlocalElement,
180 int,
181 int,
182 int,
183 << "You tried to access element " << arg1
184 << " of a distributed vector, but only elements in range [" << arg2
185 << ',' << arg3 << "] are stored locally and can be accessed."
186 << "\n\n"
187 << "A common source for this kind of problem is that you "
188 << "are passing a 'fully distributed' vector into a function "
189 << "that needs read access to vector elements that correspond "
190 << "to degrees of freedom on ghost cells (or at least to "
191 << "'locally active' degrees of freedom that are not also "
192 << "'locally owned'). You need to pass a vector that has these "
193 << "elements as ghost entries.");
197 DeclException2(ExcWrongMode,
198 int,
199 int,
200 << "You tried to do a "
201 << (arg1 == 1 ? "'set'" : (arg1 == 2 ? "'add'" : "???"))
202 << " operation but the vector is currently in "
203 << (arg2 == 1 ? "'set'" : (arg2 == 2 ? "'add'" : "???"))
204 << " mode. You first have to call 'compress()'.");
205
206 private:
210 const VectorBase &vector;
211
215 const size_type index;
216
217 // Make the vector class a friend, so that it can create objects of the
218 // present type.
219 friend class ::PETScWrappers::VectorBase;
220 };
221 } // namespace internal
252 class VectorBase : public ReadVector<PetscScalar>, public Subscriptor
253 {
254 public:
260 using value_type = PetscScalar;
261 using real_type = PetscReal;
265
270 VectorBase();
271
276 VectorBase(const VectorBase &v);
277
282 explicit VectorBase(const Vec &v);
283
287 virtual ~VectorBase() override;
288
293 virtual void
294 clear();
295
306 void
307 compress(const VectorOperation::values operation);
308
312 VectorBase &
313 operator=(const VectorBase &);
314
328 VectorBase &
329 operator=(const PetscScalar s);
330
337 void
338 reinit(Vec v);
339
345 bool
346 operator==(const VectorBase &v) const;
347
353 bool
354 operator!=(const VectorBase &v) const;
355
360 size() const override;
361
371 locally_owned_size() const;
372
381 std::pair<size_type, size_type>
382 local_range() const;
383
388 bool
389 in_local_range(const size_type index) const;
390
406
413 bool
415
419 const IndexSet &
421
425 void
427
432 operator()(const size_type index);
433
437 PetscScalar
438 operator()(const size_type index) const;
439
446 operator[](const size_type index);
447
453 PetscScalar
454 operator[](const size_type index) const;
455
462 void
463 set(const std::vector<size_type> &indices,
464 const std::vector<PetscScalar> &values);
465
481 void
482 extract_subvector_to(const std::vector<size_type> &indices,
483 std::vector<PetscScalar> &values) const;
484
488 virtual void
491 ArrayView<PetscScalar> &elements) const override;
492
520 template <typename ForwardIterator, typename OutputIterator>
521 void
522 extract_subvector_to(const ForwardIterator indices_begin,
523 const ForwardIterator indices_end,
524 OutputIterator values_begin) const;
525
530 void
531 add(const std::vector<size_type> &indices,
532 const std::vector<PetscScalar> &values);
533
538 void
539 add(const std::vector<size_type> &indices,
540 const ::Vector<PetscScalar> &values);
541
547 void
548 add(const size_type n_elements,
549 const size_type *indices,
550 const PetscScalar *values);
551
558 PetscScalar
559 operator*(const VectorBase &vec) const;
560
565 norm_sqr() const;
566
570 PetscScalar
571 mean_value() const;
572
581 l1_norm() const;
582
588 l2_norm() const;
589
595 lp_norm(const real_type p) const;
596
602 linfty_norm() const;
603
623 PetscScalar
624 add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W);
625
631 bool
632 all_zero() const;
633
637 VectorBase &
638 operator*=(const PetscScalar factor);
639
643 VectorBase &
644 operator/=(const PetscScalar factor);
645
649 VectorBase &
650 operator+=(const VectorBase &V);
651
655 VectorBase &
656 operator-=(const VectorBase &V);
657
662 void
663 add(const PetscScalar s);
664
668 void
669 add(const PetscScalar a, const VectorBase &V);
670
674 void
675 add(const PetscScalar a,
676 const VectorBase &V,
677 const PetscScalar b,
678 const VectorBase &W);
679
683 void
684 sadd(const PetscScalar s, const VectorBase &V);
685
689 void
690 sadd(const PetscScalar s, const PetscScalar a, const VectorBase &V);
691
697 void
698 scale(const VectorBase &scaling_factors);
699
703 void
704 equ(const PetscScalar a, const VectorBase &V);
705
713 void
714 write_ascii(const PetscViewerFormat format = PETSC_VIEWER_DEFAULT);
715
723 void
724 print(std::ostream &out,
725 const unsigned int precision = 3,
726 const bool scientific = true,
727 const bool across = true) const;
728
736 template <class Archive>
737 void
738 save(Archive &ar, const unsigned int version) const;
739
745 template <class Archive>
746 void
747 load(Archive &ar, const unsigned int version);
748
749# ifdef DOXYGEN
755 template <class Archive>
756 void
757 serialize(Archive &archive, const unsigned int version);
758# else
759 // This macro defines the serialize() method that is compatible with
760 // the templated save() and load() method that have been implemented.
761 BOOST_SERIALIZATION_SPLIT_MEMBER()
762# endif
763
776 void
777 swap(VectorBase &v) noexcept;
778
786 operator const Vec &() const;
787
793 Vec &
794 petsc_vector();
795
799 std::size_t
800 memory_consumption() const;
801
807
808 protected:
814
821
828
835
836 // Make the reference class a friend.
838
844 void
845 do_set_add_operation(const size_type n_elements,
846 const size_type *indices,
847 const PetscScalar *values,
848 const bool add_values);
849
853 void
855 };
856
857
858
859 // ------------------- inline and template functions --------------
860
868 inline void
869 swap(VectorBase &u, VectorBase &v) noexcept
870 {
871 u.swap(v);
872 }
873
874# ifndef DOXYGEN
875 namespace internal
876 {
877 inline VectorReference::VectorReference(const VectorBase &vector,
878 const size_type index)
879 : vector(vector)
880 , index(index)
881 {}
882
883
884 inline const VectorReference &
885 VectorReference::operator=(const VectorReference &r) const
886 {
887 // as explained in the class
888 // documentation, this is not the copy
889 // operator. so simply pass on to the
890 // "correct" assignment operator
891 *this = static_cast<PetscScalar>(r);
892
893 return *this;
894 }
895
896
897
898 inline VectorReference &
899 VectorReference::operator=(const VectorReference &r)
900 {
901 // as explained in the class
902 // documentation, this is not the copy
903 // operator. so simply pass on to the
904 // "correct" assignment operator
905 *this = static_cast<PetscScalar>(r);
906
907 return *this;
908 }
909
910
911
912 inline const VectorReference &
913 VectorReference::operator=(const PetscScalar &value) const
914 {
915 Assert((vector.last_action == VectorOperation::insert) ||
916 (vector.last_action == VectorOperation::unknown),
917 ExcWrongMode(VectorOperation::insert, vector.last_action));
918
919 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
920
921 const PetscInt petsc_i = index;
922
923 const PetscErrorCode ierr =
924 VecSetValues(vector, 1, &petsc_i, &value, INSERT_VALUES);
925 AssertThrow(ierr == 0, ExcPETScError(ierr));
926
927 vector.last_action = VectorOperation::insert;
928
929 return *this;
930 }
931
932
933
934 inline const VectorReference &
935 VectorReference::operator+=(const PetscScalar &value) const
936 {
937 Assert((vector.last_action == VectorOperation::add) ||
938 (vector.last_action == VectorOperation::unknown),
939 ExcWrongMode(VectorOperation::add, vector.last_action));
940
941 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
942
943 vector.last_action = VectorOperation::add;
944
945 // we have to do above actions in any
946 // case to be consistent with the MPI
947 // communication model (see the
948 // comments in the documentation of
949 // PETScWrappers::MPI::Vector), but we
950 // can save some work if the addend is
951 // zero
952 if (value == PetscScalar())
953 return *this;
954
955 // use the PETSc function to add something
956 const PetscInt petsc_i = index;
957 const PetscErrorCode ierr =
958 VecSetValues(vector, 1, &petsc_i, &value, ADD_VALUES);
959 AssertThrow(ierr == 0, ExcPETScError(ierr));
960
961
962 return *this;
963 }
964
965
966
967 inline const VectorReference &
968 VectorReference::operator-=(const PetscScalar &value) const
969 {
970 Assert((vector.last_action == VectorOperation::add) ||
971 (vector.last_action == VectorOperation::unknown),
972 ExcWrongMode(VectorOperation::add, vector.last_action));
973
974 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
975
976 vector.last_action = VectorOperation::add;
977
978 // we have to do above actions in any
979 // case to be consistent with the MPI
980 // communication model (see the
981 // comments in the documentation of
982 // PETScWrappers::MPI::Vector), but we
983 // can save some work if the addend is
984 // zero
985 if (value == PetscScalar())
986 return *this;
987
988 // use the PETSc function to
989 // add something
990 const PetscInt petsc_i = index;
991 const PetscScalar subtractand = -value;
992 const PetscErrorCode ierr =
993 VecSetValues(vector, 1, &petsc_i, &subtractand, ADD_VALUES);
994 AssertThrow(ierr == 0, ExcPETScError(ierr));
995
996 return *this;
997 }
998
999
1000
1001 inline const VectorReference &
1002 VectorReference::operator*=(const PetscScalar &value) const
1003 {
1004 Assert((vector.last_action == VectorOperation::insert) ||
1005 (vector.last_action == VectorOperation::unknown),
1006 ExcWrongMode(VectorOperation::insert, vector.last_action));
1007
1008 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
1009
1010 vector.last_action = VectorOperation::insert;
1011
1012 // we have to do above actions in any
1013 // case to be consistent with the MPI
1014 // communication model (see the
1015 // comments in the documentation of
1016 // PETScWrappers::MPI::Vector), but we
1017 // can save some work if the factor is
1018 // one
1019 if (value == 1.)
1020 return *this;
1021
1022 const PetscInt petsc_i = index;
1023 const PetscScalar new_value = static_cast<PetscScalar>(*this) * value;
1024
1025 const PetscErrorCode ierr =
1026 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1027 AssertThrow(ierr == 0, ExcPETScError(ierr));
1028
1029 return *this;
1030 }
1031
1032
1033
1034 inline const VectorReference &
1035 VectorReference::operator/=(const PetscScalar &value) const
1036 {
1037 Assert((vector.last_action == VectorOperation::insert) ||
1038 (vector.last_action == VectorOperation::unknown),
1039 ExcWrongMode(VectorOperation::insert, vector.last_action));
1040
1041 Assert(!vector.has_ghost_elements(), ExcGhostsPresent());
1042
1043 vector.last_action = VectorOperation::insert;
1044
1045 // we have to do above actions in any
1046 // case to be consistent with the MPI
1047 // communication model (see the
1048 // comments in the documentation of
1049 // PETScWrappers::MPI::Vector), but we
1050 // can save some work if the factor is
1051 // one
1052 if (value == 1.)
1053 return *this;
1054
1055 const PetscInt petsc_i = index;
1056 const PetscScalar new_value = static_cast<PetscScalar>(*this) / value;
1057
1058 const PetscErrorCode ierr =
1059 VecSetValues(vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1060 AssertThrow(ierr == 0, ExcPETScError(ierr));
1061
1062 return *this;
1063 }
1064
1065
1066
1067 inline PetscReal
1068 VectorReference::real() const
1069 {
1070# ifndef PETSC_USE_COMPLEX
1071 return static_cast<PetscScalar>(*this);
1072# else
1073 return PetscRealPart(static_cast<PetscScalar>(*this));
1074# endif
1075 }
1076
1077
1078
1079 inline PetscReal
1080 VectorReference::imag() const
1081 {
1082# ifndef PETSC_USE_COMPLEX
1083 return PetscReal(0);
1084# else
1085 return PetscImaginaryPart(static_cast<PetscScalar>(*this));
1086# endif
1087 }
1088
1089 } // namespace internal
1090
1091 inline bool
1092 VectorBase::in_local_range(const size_type index) const
1093 {
1094 PetscInt begin, end;
1095 const PetscErrorCode ierr =
1096 VecGetOwnershipRange(static_cast<const Vec &>(vector), &begin, &end);
1097 AssertThrow(ierr == 0, ExcPETScError(ierr));
1098
1099 return ((index >= static_cast<size_type>(begin)) &&
1100 (index < static_cast<size_type>(end)));
1101 }
1102
1103
1104 inline IndexSet
1105 VectorBase::locally_owned_elements() const
1106 {
1107 IndexSet is(size());
1108
1109 // PETSc only allows for contiguous local ranges, so this is simple
1110 const std::pair<size_type, size_type> x = local_range();
1111 is.add_range(x.first, x.second);
1112 return is;
1113 }
1114
1115
1116
1117 inline bool
1118 VectorBase::has_ghost_elements() const
1119 {
1120 return ghosted;
1121 }
1122
1123
1124 inline const IndexSet &
1125 VectorBase::ghost_elements() const
1126 {
1127 return ghost_indices;
1128 }
1129
1130
1131 inline void
1132 VectorBase::update_ghost_values() const
1133 {
1134 if (ghosted)
1135 {
1136 PetscErrorCode ierr;
1137
1138 ierr = VecGhostUpdateBegin(vector, INSERT_VALUES, SCATTER_FORWARD);
1139 AssertThrow(ierr == 0, ExcPETScError(ierr));
1140 ierr = VecGhostUpdateEnd(vector, INSERT_VALUES, SCATTER_FORWARD);
1141 AssertThrow(ierr == 0, ExcPETScError(ierr));
1142 }
1143 }
1144
1145
1146
1147 inline internal::VectorReference
1148 VectorBase::operator()(const size_type index)
1149 {
1150 return internal::VectorReference(*this, index);
1151 }
1152
1153
1154
1155 inline PetscScalar
1156 VectorBase::operator()(const size_type index) const
1157 {
1158 return static_cast<PetscScalar>(internal::VectorReference(*this, index));
1159 }
1160
1161
1162
1163 inline internal::VectorReference
1164 VectorBase::operator[](const size_type index)
1165 {
1166 return operator()(index);
1167 }
1168
1169
1170
1171 inline PetscScalar
1172 VectorBase::operator[](const size_type index) const
1173 {
1174 return operator()(index);
1175 }
1176
1177 inline MPI_Comm
1178 VectorBase::get_mpi_communicator() const
1179 {
1180 return PetscObjectComm(reinterpret_cast<PetscObject>(vector));
1181 }
1182
1183 inline void
1184 VectorBase::extract_subvector_to(const std::vector<size_type> &indices,
1185 std::vector<PetscScalar> &values) const
1186 {
1187 Assert(indices.size() <= values.size(),
1188 ExcDimensionMismatch(indices.size(), values.size()));
1189 extract_subvector_to(indices.begin(), indices.end(), values.begin());
1190 }
1191
1192 inline void
1193 VectorBase::extract_subvector_to(
1195 ArrayView<PetscScalar> &elements) const
1196 {
1197 AssertDimension(indices.size(), elements.size());
1198 extract_subvector_to(indices.begin(), indices.end(), elements.begin());
1199 }
1200
1201
1202 template <typename ForwardIterator, typename OutputIterator>
1203 inline void
1204 VectorBase::extract_subvector_to(const ForwardIterator indices_begin,
1205 const ForwardIterator indices_end,
1206 OutputIterator values_begin) const
1207 {
1208 const PetscInt n_idx = static_cast<PetscInt>(indices_end - indices_begin);
1209 if (n_idx == 0)
1210 return;
1211
1212 // if we are dealing
1213 // with a parallel vector
1214 if (ghosted)
1215 {
1216 // there is the possibility
1217 // that the vector has
1218 // ghost elements. in that
1219 // case, we first need to
1220 // figure out which
1221 // elements we own locally,
1222 // then get a pointer to
1223 // the elements that are
1224 // stored here (both the
1225 // ones we own as well as
1226 // the ghost elements). in
1227 // this array, the locally
1228 // owned elements come
1229 // first followed by the
1230 // ghost elements whose
1231 // position we can get from
1232 // an index set
1233 PetscInt begin, end;
1234 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1235 AssertThrow(ierr == 0, ExcPETScError(ierr));
1236
1237 Vec locally_stored_elements = nullptr;
1238 ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1239 AssertThrow(ierr == 0, ExcPETScError(ierr));
1240
1241 PetscInt lsize;
1242 ierr = VecGetSize(locally_stored_elements, &lsize);
1243 AssertThrow(ierr == 0, ExcPETScError(ierr));
1244
1245 const PetscScalar *ptr;
1246 ierr = VecGetArrayRead(locally_stored_elements, &ptr);
1247 AssertThrow(ierr == 0, ExcPETScError(ierr));
1248
1249 for (PetscInt i = 0; i < n_idx; ++i)
1250 {
1251 const unsigned int index = *(indices_begin + i);
1252 if (index >= static_cast<unsigned int>(begin) &&
1253 index < static_cast<unsigned int>(end))
1254 {
1255 // local entry
1256 *(values_begin + i) = *(ptr + index - begin);
1257 }
1258 else
1259 {
1260 // ghost entry
1261 const unsigned int ghostidx =
1262 ghost_indices.index_within_set(index);
1263
1264 AssertIndexRange(ghostidx + end - begin, lsize);
1265 *(values_begin + i) = *(ptr + ghostidx + end - begin);
1266 }
1267 }
1268
1269 ierr = VecRestoreArrayRead(locally_stored_elements, &ptr);
1270 AssertThrow(ierr == 0, ExcPETScError(ierr));
1271
1272 ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1273 AssertThrow(ierr == 0, ExcPETScError(ierr));
1274 }
1275 // if the vector is local or the
1276 // caller, then simply access the
1277 // element we are interested in
1278 else
1279 {
1280 PetscInt begin, end;
1281 PetscErrorCode ierr = VecGetOwnershipRange(vector, &begin, &end);
1282 AssertThrow(ierr == 0, ExcPETScError(ierr));
1283
1284 const PetscScalar *ptr;
1285 ierr = VecGetArrayRead(vector, &ptr);
1286 AssertThrow(ierr == 0, ExcPETScError(ierr));
1287
1288 for (PetscInt i = 0; i < n_idx; ++i)
1289 {
1290 const unsigned int index = *(indices_begin + i);
1291
1292 Assert(index >= static_cast<unsigned int>(begin) &&
1293 index < static_cast<unsigned int>(end),
1294 ExcMessage("You are accessing elements of a vector without "
1295 "ghost elements that are not actually owned by "
1296 "this vector. A typical case where this may "
1297 "happen is if you are passing a non-ghosted "
1298 "(completely distributed) vector to a function "
1299 "that expects a vector that stores ghost "
1300 "elements for all locally relevant or locally "
1301 "active vector entries."));
1302
1303 *(values_begin + i) = *(ptr + index - begin);
1304 }
1305
1306 ierr = VecRestoreArrayRead(vector, &ptr);
1307 AssertThrow(ierr == 0, ExcPETScError(ierr));
1308 }
1309 }
1310
1311 template <class Archive>
1312 inline void
1313 VectorBase::save(Archive &ar, const unsigned int) const
1314 {
1315 // forward to serialization function in the base class.
1316 ar &static_cast<const Subscriptor &>(*this);
1317 ar &size();
1318 ar &local_range();
1319
1320 const PetscScalar *array = nullptr;
1321 int ierr = VecGetArrayRead(*this, &array);
1322 AssertThrow(ierr == 0, ExcPETScError(ierr));
1323
1324 boost::serialization::array_wrapper<const PetscScalar> wrapper(
1325 array, locally_owned_size());
1326 ar &wrapper;
1327
1328 ierr = VecRestoreArrayRead(*this, &array);
1329 AssertThrow(ierr == 0, ExcPETScError(ierr));
1330 }
1331
1332
1333
1334 template <class Archive>
1335 inline void
1336 VectorBase::load(Archive &ar, const unsigned int)
1337 {
1338 ar &static_cast<Subscriptor &>(*this);
1339
1340 size_type size = 0;
1341 std::pair<size_type, size_type> local_range;
1342
1343 ar &size;
1344 Assert(size == this->size(),
1345 ExcMessage("The serialized value of size (" + std::to_string(size) +
1346 ") does not match the current size (" +
1347 std::to_string(this->size()) + ")"));
1348 (void)size;
1349 ar &local_range;
1350 Assert(local_range == this->local_range(),
1351 ExcMessage("The serialized value of local_range (" +
1352 std::to_string(local_range.first) + ", " +
1353 std::to_string(local_range.second) +
1354 ") does not match the current local_range (" +
1355 std::to_string(this->local_range().first) + ", " +
1356 std::to_string(this->local_range().second) + ")"));
1357 (void)local_range;
1358
1359 PetscScalar *array = nullptr;
1360 int ierr = VecGetArray(petsc_vector(), &array);
1361 AssertThrow(ierr == 0, ExcPETScError(ierr));
1362
1363 boost::serialization::array_wrapper<PetscScalar> wrapper(
1364 array, locally_owned_size());
1365 ar &wrapper;
1366
1367 ierr = VecRestoreArray(petsc_vector(), &array);
1368 AssertThrow(ierr == 0, ExcPETScError(ierr));
1369 }
1370# endif // DOXYGEN
1371} // namespace PETScWrappers
1372
1374
1375#endif // DEAL_II_WITH_PETSC
1376
1377#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
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)
virtual void extract_subvector_to(const ArrayView< const types::global_dof_index > &indices, ArrayView< PetscScalar > &elements) const override
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:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
static ::ExceptionBase & ExcGhostsPresent()
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition exceptions.h:562
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int global_dof_index
Definition types.h:81