Reference documentation for deal.II version 8.4.2
petsc_vector_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2015 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 at
12 // the top level of the deal.II distribution.
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 
24 # include <deal.II/base/subscriptor.h>
25 # include <deal.II/lac/exceptions.h>
26 # include <deal.II/lac/vector.h>
27 
28 # include <vector>
29 # include <utility>
30 
31 # include <petscvec.h>
32 # include <deal.II/base/index_set.h>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 // forward declaration
37 template <typename number> class Vector;
38 
39 
47 namespace PETScWrappers
48 {
49  // forward declaration
50  class VectorBase;
51 
61  namespace internal
62  {
75  class VectorReference
76  {
77  public:
82 
83  private:
88  VectorReference (const VectorBase &vector,
89  const size_type index);
90 
91 
92  public:
93 
105  const VectorReference &operator = (const VectorReference &r) const;
106 
112  VectorReference &operator = (const VectorReference &r);
113 
117  const VectorReference &operator = (const PetscScalar &s) const;
118 
122  const VectorReference &operator += (const PetscScalar &s) const;
123 
127  const VectorReference &operator -= (const PetscScalar &s) const;
128 
132  const VectorReference &operator *= (const PetscScalar &s) const;
133 
137  const VectorReference &operator /= (const PetscScalar &s) const;
138 
142  PetscReal real () const;
143 
150  PetscReal imag () const;
151 
156  operator PetscScalar () const;
157 
161  DeclException1 (ExcPETScError,
162  int,
163  << "An error with error number " << arg1
164  << " occurred while calling a PETSc function");
168  DeclException3 (ExcAccessToNonlocalElement,
169  int, int, int,
170  << "You tried to access element " << arg1
171  << " of a distributed vector, but only elements "
172  << arg2 << " through " << arg3
173  << " are stored locally and can be accessed.");
177  DeclException2 (ExcWrongMode,
178  int, int,
179  << "You tried to do a "
180  << (arg1 == 1 ?
181  "'set'" :
182  (arg1 == 2 ?
183  "'add'" : "???"))
184  << " operation but the vector is currently in "
185  << (arg2 == 1 ?
186  "'set'" :
187  (arg2 == 2 ?
188  "'add'" : "???"))
189  << " mode. You first have to call 'compress()'.");
190 
191  private:
195  const VectorBase &vector;
196 
200  const size_type index;
201 
206  friend class ::PETScWrappers::VectorBase;
207  };
208  }
238  class VectorBase : public Subscriptor
239  {
240  public:
246  typedef PetscScalar value_type;
247  typedef PetscReal real_type;
248  typedef types::global_dof_index size_type;
249  typedef internal::VectorReference reference;
250  typedef const internal::VectorReference const_reference;
251 
256  VectorBase ();
257 
262  VectorBase (const VectorBase &v);
263 
269  explicit VectorBase (const Vec &v);
270 
274  virtual ~VectorBase ();
275 
280  virtual void clear ();
281 
292  void compress (const VectorOperation::values operation);
293 
307  VectorBase &operator = (const PetscScalar s);
308 
314  bool operator == (const VectorBase &v) const;
315 
321  bool operator != (const VectorBase &v) const;
322 
326  size_type size () const;
327 
336  size_type local_size () const;
337 
346  std::pair<size_type, size_type>
347  local_range () const;
348 
353  bool in_local_range (const size_type index) const;
354 
368  IndexSet locally_owned_elements () const;
369 
376  bool has_ghost_elements() const;
377 
381  reference
382  operator () (const size_type index);
383 
387  PetscScalar
388  operator () (const size_type index) const;
389 
395  reference
396  operator [] (const size_type index);
397 
403  PetscScalar
404  operator [] (const size_type index) const;
405 
412  void set (const std::vector<size_type> &indices,
413  const std::vector<PetscScalar> &values);
414 
421  void extract_subvector_to (const std::vector<size_type> &indices,
422  std::vector<PetscScalar> &values) const;
423 
428  template <typename ForwardIterator, typename OutputIterator>
429  void extract_subvector_to (const ForwardIterator indices_begin,
430  const ForwardIterator indices_end,
431  OutputIterator values_begin) const;
432 
437  void add (const std::vector<size_type> &indices,
438  const std::vector<PetscScalar> &values);
439 
444  void add (const std::vector<size_type> &indices,
445  const ::Vector<PetscScalar> &values);
446 
452  void add (const size_type n_elements,
453  const size_type *indices,
454  const PetscScalar *values);
455 
462  PetscScalar operator * (const VectorBase &vec) const;
463 
467  real_type norm_sqr () const;
468 
472  PetscScalar mean_value () const;
473 
477  real_type l1_norm () const;
478 
483  real_type l2_norm () const;
484 
489  real_type lp_norm (const real_type p) const;
490 
495  real_type linfty_norm () const;
496 
512  PetscScalar add_and_dot (const PetscScalar a,
513  const VectorBase &V,
514  const VectorBase &W);
515 
522  real_type normalize () const DEAL_II_DEPRECATED;
523 
527  real_type min () const;
528 
532  real_type max () const;
533 
540 
546  VectorBase &conjugate () DEAL_II_DEPRECATED;
547 
555  VectorBase &mult () DEAL_II_DEPRECATED;
556 
563  VectorBase &mult (const VectorBase &v) DEAL_II_DEPRECATED;
564 
571  VectorBase &mult (const VectorBase &u,
572  const VectorBase &v) DEAL_II_DEPRECATED;
573 
579  bool all_zero () const;
580 
586  bool is_non_negative () const;
587 
591  VectorBase &operator *= (const PetscScalar factor);
592 
596  VectorBase &operator /= (const PetscScalar factor);
597 
601  VectorBase &operator += (const VectorBase &V);
602 
606  VectorBase &operator -= (const VectorBase &V);
607 
612  void add (const PetscScalar s);
613 
619  void add (const VectorBase &V) DEAL_II_DEPRECATED;
620 
624  void add (const PetscScalar a, const VectorBase &V);
625 
629  void add (const PetscScalar a, const VectorBase &V,
630  const PetscScalar b, const VectorBase &W);
631 
635  void sadd (const PetscScalar s,
636  const VectorBase &V);
637 
641  void sadd (const PetscScalar s,
642  const PetscScalar a,
643  const VectorBase &V);
644 
650  void sadd (const PetscScalar s,
651  const PetscScalar a,
652  const VectorBase &V,
653  const PetscScalar b,
654  const VectorBase &W) DEAL_II_DEPRECATED;
655 
662  void sadd (const PetscScalar s,
663  const PetscScalar a,
664  const VectorBase &V,
665  const PetscScalar b,
666  const VectorBase &W,
667  const PetscScalar c,
668  const VectorBase &X) DEAL_II_DEPRECATED;
669 
675  void scale (const VectorBase &scaling_factors);
676 
680  void equ (const PetscScalar a, const VectorBase &V);
681 
687  void equ (const PetscScalar a, const VectorBase &V,
688  const PetscScalar b, const VectorBase &W) DEAL_II_DEPRECATED;
689 
700  void ratio (const VectorBase &a,
701  const VectorBase &b) DEAL_II_DEPRECATED;
702 
710  void write_ascii (const PetscViewerFormat format = PETSC_VIEWER_DEFAULT) ;
711 
719  void print (std::ostream &out,
720  const unsigned int precision = 3,
721  const bool scientific = true,
722  const bool across = true) const;
723 
736  void swap (VectorBase &v);
737 
745  operator const Vec &() const;
746 
750  std::size_t memory_consumption () const;
751 
756  virtual const MPI_Comm &get_mpi_communicator () const;
757 
758  protected:
763  Vec vector;
764 
770  bool ghosted;
771 
777  IndexSet ghost_indices;
778 
784  mutable VectorOperation::values last_action;
785 
789  friend class internal::VectorReference;
790 
796  bool attained_ownership;
797 
803  void do_set_add_operation (const size_type n_elements,
804  const size_type *indices,
805  const PetscScalar *values,
806  const bool add_values);
807 
808 
809  };
810 
811 
812 
813 // ------------------- inline and template functions --------------
814 
823  inline
824  void swap (VectorBase &u, VectorBase &v)
825  {
826  u.swap (v);
827  }
828 
829 #ifndef DOXYGEN
830  namespace internal
831  {
832  inline
833  VectorReference::VectorReference (const VectorBase &vector,
834  const size_type index)
835  :
836  vector (vector),
837  index (index)
838  {}
839 
840 
841  inline
842  const VectorReference &
843  VectorReference::operator = (const VectorReference &r) const
844  {
845  // as explained in the class
846  // documentation, this is not the copy
847  // operator. so simply pass on to the
848  // "correct" assignment operator
849  *this = static_cast<PetscScalar> (r);
850 
851  return *this;
852  }
853 
854 
855 
856  inline
857  VectorReference &
858  VectorReference::operator = (const VectorReference &r)
859  {
860  // as explained in the class
861  // documentation, this is not the copy
862  // operator. so simply pass on to the
863  // "correct" assignment operator
864  *this = static_cast<PetscScalar> (r);
865 
866  return *this;
867  }
868 
869 
870 
871  inline
872  const VectorReference &
873  VectorReference::operator = (const PetscScalar &value) const
874  {
875  Assert ((vector.last_action == VectorOperation::insert)
876  ||
877  (vector.last_action == VectorOperation::unknown),
878  ExcWrongMode (VectorOperation::insert,
879  vector.last_action));
880 
881  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
882 
883  const PetscInt petsc_i = index;
884 
885  const int ierr
886  = VecSetValues (vector, 1, &petsc_i, &value, INSERT_VALUES);
887  AssertThrow (ierr == 0, ExcPETScError(ierr));
888 
889  vector.last_action = VectorOperation::insert;
890 
891  return *this;
892  }
893 
894 
895 
896  inline
897  const VectorReference &
898  VectorReference::operator += (const PetscScalar &value) const
899  {
900  Assert ((vector.last_action == VectorOperation::add)
901  ||
902  (vector.last_action == VectorOperation::unknown),
903  ExcWrongMode (VectorOperation::add,
904  vector.last_action));
905 
906  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
907 
908  vector.last_action = VectorOperation::add;
909 
910  // we have to do above actions in any
911  // case to be consistent with the MPI
912  // communication model (see the
913  // comments in the documentation of
914  // PETScWrappers::MPI::Vector), but we
915  // can save some work if the addend is
916  // zero
917  if (value == PetscScalar())
918  return *this;
919 
920  // use the PETSc function to add something
921  const PetscInt petsc_i = index;
922  const int ierr
923  = VecSetValues (vector, 1, &petsc_i, &value, ADD_VALUES);
924  AssertThrow (ierr == 0, ExcPETScError(ierr));
925 
926 
927  return *this;
928  }
929 
930 
931 
932  inline
933  const VectorReference &
934  VectorReference::operator -= (const PetscScalar &value) const
935  {
936  Assert ((vector.last_action == VectorOperation::add)
937  ||
938  (vector.last_action == VectorOperation::unknown),
939  ExcWrongMode (VectorOperation::add,
940  vector.last_action));
941 
942  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
943 
944  vector.last_action = VectorOperation::add;
945 
946  // we have to do above actions in any
947  // case to be consistent with the MPI
948  // communication model (see the
949  // comments in the documentation of
950  // PETScWrappers::MPI::Vector), but we
951  // can save some work if the addend is
952  // zero
953  if (value == PetscScalar())
954  return *this;
955 
956  // use the PETSc function to
957  // add something
958  const PetscInt petsc_i = index;
959  const PetscScalar subtractand = -value;
960  const int ierr
961  = VecSetValues (vector, 1, &petsc_i, &subtractand, ADD_VALUES);
962  AssertThrow (ierr == 0, ExcPETScError(ierr));
963 
964  return *this;
965  }
966 
967 
968 
969  inline
970  const VectorReference &
971  VectorReference::operator *= (const PetscScalar &value) const
972  {
973  Assert ((vector.last_action == VectorOperation::insert)
974  ||
975  (vector.last_action == VectorOperation::unknown),
976  ExcWrongMode (VectorOperation::insert,
977  vector.last_action));
978 
979  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
980 
981  vector.last_action = VectorOperation::insert;
982 
983  // we have to do above actions in any
984  // case to be consistent with the MPI
985  // communication model (see the
986  // comments in the documentation of
987  // PETScWrappers::MPI::Vector), but we
988  // can save some work if the factor is
989  // one
990  if (value == 1.)
991  return *this;
992 
993  const PetscInt petsc_i = index;
994  const PetscScalar new_value
995  = static_cast<PetscScalar>(*this) * value;
996 
997  const int ierr
998  = VecSetValues (vector, 1, &petsc_i, &new_value, INSERT_VALUES);
999  AssertThrow (ierr == 0, ExcPETScError(ierr));
1000 
1001  return *this;
1002  }
1003 
1004 
1005 
1006  inline
1007  const VectorReference &
1008  VectorReference::operator /= (const PetscScalar &value) const
1009  {
1010  Assert ((vector.last_action == VectorOperation::insert)
1011  ||
1012  (vector.last_action == VectorOperation::unknown),
1013  ExcWrongMode (VectorOperation::insert,
1014  vector.last_action));
1015 
1016  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
1017 
1018  vector.last_action = VectorOperation::insert;
1019 
1020  // we have to do above actions in any
1021  // case to be consistent with the MPI
1022  // communication model (see the
1023  // comments in the documentation of
1024  // PETScWrappers::MPI::Vector), but we
1025  // can save some work if the factor is
1026  // one
1027  if (value == 1.)
1028  return *this;
1029 
1030  const PetscInt petsc_i = index;
1031  const PetscScalar new_value
1032  = static_cast<PetscScalar>(*this) / value;
1033 
1034  const int ierr
1035  = VecSetValues (vector, 1, &petsc_i, &new_value, INSERT_VALUES);
1036  AssertThrow (ierr == 0, ExcPETScError(ierr));
1037 
1038  return *this;
1039  }
1040 
1041 
1042 
1043  inline
1044  PetscReal
1045  VectorReference::real () const
1046  {
1047 #ifndef PETSC_USE_COMPLEX
1048  return static_cast<PetscScalar>(*this);
1049 #else
1050  return PetscRealPart (static_cast<PetscScalar>(*this));
1051 #endif
1052  }
1053 
1054 
1055 
1056  inline
1057  PetscReal
1058  VectorReference::imag () const
1059  {
1060 #ifndef PETSC_USE_COMPLEX
1061  return PetscReal (0);
1062 #else
1063  return PetscImaginaryPart (static_cast<PetscScalar>(*this));
1064 #endif
1065  }
1066 
1067  } // namespace internal
1068 
1069  inline
1070  bool
1071  VectorBase::in_local_range (const size_type index) const
1072  {
1073  PetscInt begin, end;
1074  const int ierr = VecGetOwnershipRange (static_cast<const Vec &>(vector),
1075  &begin, &end);
1076  AssertThrow (ierr == 0, ExcPETScError(ierr));
1077 
1078  return ((index >= static_cast<size_type>(begin)) &&
1079  (index < static_cast<size_type>(end)));
1080  }
1081 
1082 
1083  inline
1084  IndexSet
1086  {
1087  IndexSet is (size());
1088 
1089  // PETSc only allows for contiguous local ranges, so this is simple
1090  const std::pair<size_type, size_type> x = local_range();
1091  is.add_range (x.first, x.second);
1092  return is;
1093  }
1094 
1095 
1096 
1097  inline
1098  bool
1100  {
1101  return ghosted;
1102  }
1103 
1104 
1105 
1106  inline
1107  internal::VectorReference
1108  VectorBase::operator () (const size_type index)
1109  {
1110  return internal::VectorReference (*this, index);
1111  }
1112 
1113 
1114 
1115  inline
1116  PetscScalar
1117  VectorBase::operator () (const size_type index) const
1118  {
1119  return static_cast<PetscScalar>(internal::VectorReference (*this, index));
1120  }
1121 
1122 
1123 
1124  inline
1125  internal::VectorReference
1126  VectorBase::operator [] (const size_type index)
1127  {
1128  return operator()(index);
1129  }
1130 
1131 
1132 
1133  inline
1134  PetscScalar
1135  VectorBase::operator [] (const size_type index) const
1136  {
1137  return operator()(index);
1138  }
1139 
1140  inline
1141  const MPI_Comm &
1143  {
1144  static MPI_Comm comm;
1145  PetscObjectGetComm((PetscObject)vector, &comm);
1146  return comm;
1147  }
1148 
1149  inline
1150  void VectorBase::extract_subvector_to (const std::vector<size_type> &indices,
1151  std::vector<PetscScalar> &values) const
1152  {
1153  extract_subvector_to(&(indices[0]), &(indices[0]) + indices.size(), &(values[0]));
1154  }
1155 
1156  template <typename ForwardIterator, typename OutputIterator>
1157  inline
1158  void VectorBase::extract_subvector_to (const ForwardIterator indices_begin,
1159  const ForwardIterator indices_end,
1160  OutputIterator values_begin) const
1161  {
1162  const PetscInt n_idx = static_cast<PetscInt>(indices_end - indices_begin);
1163  if (n_idx == 0)
1164  return;
1165 
1166  // if we are dealing
1167  // with a parallel vector
1168  if (ghosted )
1169  {
1170 
1171  int ierr;
1172 
1173  // there is the possibility
1174  // that the vector has
1175  // ghost elements. in that
1176  // case, we first need to
1177  // figure out which
1178  // elements we own locally,
1179  // then get a pointer to
1180  // the elements that are
1181  // stored here (both the
1182  // ones we own as well as
1183  // the ghost elements). in
1184  // this array, the locally
1185  // owned elements come
1186  // first followed by the
1187  // ghost elements whose
1188  // position we can get from
1189  // an index set
1190  PetscInt begin, end;
1191  ierr = VecGetOwnershipRange (vector, &begin, &end);
1192  AssertThrow (ierr == 0, ExcPETScError(ierr));
1193 
1194  Vec locally_stored_elements = PETSC_NULL;
1195  ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1196  AssertThrow (ierr == 0, ExcPETScError(ierr));
1197 
1198  PetscInt lsize;
1199  ierr = VecGetSize(locally_stored_elements, &lsize);
1200  AssertThrow (ierr == 0, ExcPETScError(ierr));
1201 
1202  PetscScalar *ptr;
1203  ierr = VecGetArray(locally_stored_elements, &ptr);
1204  AssertThrow (ierr == 0, ExcPETScError(ierr));
1205 
1206  for (PetscInt i=0; i<n_idx; ++i)
1207  {
1208  const unsigned int index = *(indices_begin+i);
1209  if ( index>=static_cast<unsigned int>(begin)
1210  && index<static_cast<unsigned int>(end) )
1211  {
1212  //local entry
1213  *(values_begin+i) = *(ptr+index-begin);
1214  }
1215  else
1216  {
1217  //ghost entry
1218  const unsigned int ghostidx
1219  = ghost_indices.index_within_set(index);
1220 
1221  Assert(ghostidx+end-begin<(unsigned int)lsize, ExcInternalError());
1222  *(values_begin+i) = *(ptr+ghostidx+end-begin);
1223  }
1224  }
1225 
1226  ierr = VecRestoreArray(locally_stored_elements, &ptr);
1227  AssertThrow (ierr == 0, ExcPETScError(ierr));
1228 
1229  ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1230  AssertThrow (ierr == 0, ExcPETScError(ierr));
1231 
1232  }
1233  // if the vector is local or the
1234  // caller, then simply access the
1235  // element we are interested in
1236  else
1237  {
1238  int ierr;
1239 
1240  PetscInt begin, end;
1241  ierr = VecGetOwnershipRange (vector, &begin, &end);
1242  AssertThrow (ierr == 0, ExcPETScError(ierr));
1243 
1244  PetscScalar *ptr;
1245  ierr = VecGetArray(vector, &ptr);
1246  AssertThrow (ierr == 0, ExcPETScError(ierr));
1247 
1248  for (PetscInt i=0; i<n_idx; ++i)
1249  {
1250  const unsigned int index = *(indices_begin+i);
1251 
1252  Assert(index>=static_cast<unsigned int>(begin)
1253  && index<static_cast<unsigned int>(end), ExcInternalError());
1254 
1255  *(values_begin+i) = *(ptr+index-begin);
1256  }
1257 
1258  ierr = VecRestoreArray(vector, &ptr);
1259  AssertThrow (ierr == 0, ExcPETScError(ierr));
1260 
1261  }
1262  }
1263 
1264 #endif // DOXYGEN
1265 }
1266 
1267 DEAL_II_NAMESPACE_CLOSE
1268 
1269 #endif // DEAL_II_WITH_PETSC
1270 
1271 /*---------------------------- petsc_vector_base.h ---------------------------*/
1272 
1273 #endif
1274 /*---------------------------- petsc_vector_base.h ---------------------------*/
reference operator[](const size_type index)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:552
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
bool in_local_range(const size_type index) const
BlockDynamicSparsityPattern BlockCompressedSparsityPattern DEAL_II_DEPRECATED
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:542
unsigned int global_dof_index
Definition: types.h:88
VectorOperation::values last_action
#define Assert(cond, exc)
Definition: exceptions.h:294
void add_range(const size_type begin, const size_type end)
Definition: index_set.cc:86
bool has_ghost_elements() const
reference operator()(const size_type index)
IndexSet locally_owned_elements() const
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:562
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< PetscScalar > &values) const
virtual const MPI_Comm & get_mpi_communicator() const