Reference documentation for deal.II version 8.5.1
petsc_vector_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 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;
160  DeclException3 (ExcAccessToNonlocalElement,
161  int, int, int,
162  << "You tried to access element " << arg1
163  << " of a distributed vector, but only elements "
164  << arg2 << " through " << arg3
165  << " are stored locally and can be accessed.");
169  DeclException2 (ExcWrongMode,
170  int, int,
171  << "You tried to do a "
172  << (arg1 == 1 ?
173  "'set'" :
174  (arg1 == 2 ?
175  "'add'" : "???"))
176  << " operation but the vector is currently in "
177  << (arg2 == 1 ?
178  "'set'" :
179  (arg2 == 2 ?
180  "'add'" : "???"))
181  << " mode. You first have to call 'compress()'.");
182 
183  private:
187  const VectorBase &vector;
188 
192  const size_type index;
193 
198  friend class ::PETScWrappers::VectorBase;
199  };
200  }
230  class VectorBase : public Subscriptor
231  {
232  public:
238  typedef PetscScalar value_type;
239  typedef PetscReal real_type;
240  typedef types::global_dof_index size_type;
241  typedef internal::VectorReference reference;
242  typedef const internal::VectorReference const_reference;
243 
248  VectorBase ();
249 
254  VectorBase (const VectorBase &v);
255 
261  explicit VectorBase (const Vec &v);
262 
266  virtual ~VectorBase ();
267 
272  virtual void clear ();
273 
284  void compress (const VectorOperation::values operation);
285 
299  VectorBase &operator = (const PetscScalar s);
300 
306  bool operator == (const VectorBase &v) const;
307 
313  bool operator != (const VectorBase &v) const;
314 
318  size_type size () const;
319 
328  size_type local_size () const;
329 
338  std::pair<size_type, size_type>
339  local_range () const;
340 
345  bool in_local_range (const size_type index) const;
346 
361 
368  bool has_ghost_elements() const;
369 
373  reference
374  operator () (const size_type index);
375 
379  PetscScalar
380  operator () (const size_type index) const;
381 
387  reference
388  operator [] (const size_type index);
389 
395  PetscScalar
396  operator [] (const size_type index) const;
397 
404  void set (const std::vector<size_type> &indices,
405  const std::vector<PetscScalar> &values);
406 
413  void extract_subvector_to (const std::vector<size_type> &indices,
414  std::vector<PetscScalar> &values) const;
415 
420  template <typename ForwardIterator, typename OutputIterator>
421  void extract_subvector_to (const ForwardIterator indices_begin,
422  const ForwardIterator indices_end,
423  OutputIterator values_begin) const;
424 
429  void add (const std::vector<size_type> &indices,
430  const std::vector<PetscScalar> &values);
431 
436  void add (const std::vector<size_type> &indices,
437  const ::Vector<PetscScalar> &values);
438 
444  void add (const size_type n_elements,
445  const size_type *indices,
446  const PetscScalar *values);
447 
454  PetscScalar operator * (const VectorBase &vec) const;
455 
459  real_type norm_sqr () const;
460 
464  PetscScalar mean_value () const;
465 
473  real_type l1_norm () const;
474 
479  real_type l2_norm () const;
480 
485  real_type lp_norm (const real_type p) const;
486 
491  real_type linfty_norm () const;
492 
508  PetscScalar add_and_dot (const PetscScalar a,
509  const VectorBase &V,
510  const VectorBase &W);
511 
518  real_type normalize () const DEAL_II_DEPRECATED;
519 
523  real_type min () const;
524 
528  real_type max () const;
529 
535  VectorBase &abs () DEAL_II_DEPRECATED;
536 
542  VectorBase &conjugate () DEAL_II_DEPRECATED;
543 
551  VectorBase &mult () DEAL_II_DEPRECATED;
552 
559  VectorBase &mult (const VectorBase &v) DEAL_II_DEPRECATED;
560 
567  VectorBase &mult (const VectorBase &u,
568  const VectorBase &v) DEAL_II_DEPRECATED;
569 
575  bool all_zero () const;
576 
582  bool is_non_negative () const;
583 
587  VectorBase &operator *= (const PetscScalar factor);
588 
592  VectorBase &operator /= (const PetscScalar factor);
593 
597  VectorBase &operator += (const VectorBase &V);
598 
602  VectorBase &operator -= (const VectorBase &V);
603 
608  void add (const PetscScalar s);
609 
615  void add (const VectorBase &V) DEAL_II_DEPRECATED;
616 
620  void add (const PetscScalar a, const VectorBase &V);
621 
625  void add (const PetscScalar a, const VectorBase &V,
626  const PetscScalar b, const VectorBase &W);
627 
631  void sadd (const PetscScalar s,
632  const VectorBase &V);
633 
637  void sadd (const PetscScalar s,
638  const PetscScalar a,
639  const VectorBase &V);
640 
646  void sadd (const PetscScalar s,
647  const PetscScalar a,
648  const VectorBase &V,
649  const PetscScalar b,
650  const VectorBase &W) DEAL_II_DEPRECATED;
651 
658  void sadd (const PetscScalar s,
659  const PetscScalar a,
660  const VectorBase &V,
661  const PetscScalar b,
662  const VectorBase &W,
663  const PetscScalar c,
664  const VectorBase &X) DEAL_II_DEPRECATED;
665 
671  void scale (const VectorBase &scaling_factors);
672 
676  void equ (const PetscScalar a, const VectorBase &V);
677 
683  void equ (const PetscScalar a, const VectorBase &V,
684  const PetscScalar b, const VectorBase &W) DEAL_II_DEPRECATED;
685 
696  void ratio (const VectorBase &a,
697  const VectorBase &b) DEAL_II_DEPRECATED;
698 
706  void write_ascii (const PetscViewerFormat format = PETSC_VIEWER_DEFAULT) ;
707 
715  void print (std::ostream &out,
716  const unsigned int precision = 3,
717  const bool scientific = true,
718  const bool across = true) const;
719 
732  void swap (VectorBase &v);
733 
741  operator const Vec &() const;
742 
746  std::size_t memory_consumption () const;
747 
752  virtual const MPI_Comm &get_mpi_communicator () const;
753 
754  protected:
759  Vec vector;
760 
766  bool ghosted;
767 
774 
780  mutable VectorOperation::values last_action;
781 
785  friend class internal::VectorReference;
786 
793 
799  void do_set_add_operation (const size_type n_elements,
800  const size_type *indices,
801  const PetscScalar *values,
802  const bool add_values);
803 
804  private:
810  VectorBase &operator=(const VectorBase &);
811  };
812 
813 
814 
815 // ------------------- inline and template functions --------------
816 
825  inline
826  void swap (VectorBase &u, VectorBase &v)
827  {
828  u.swap (v);
829  }
830 
831 #ifndef DOXYGEN
832  namespace internal
833  {
834  inline
835  VectorReference::VectorReference (const VectorBase &vector,
836  const size_type index)
837  :
838  vector (vector),
839  index (index)
840  {}
841 
842 
843  inline
844  const VectorReference &
845  VectorReference::operator = (const VectorReference &r) const
846  {
847  // as explained in the class
848  // documentation, this is not the copy
849  // operator. so simply pass on to the
850  // "correct" assignment operator
851  *this = static_cast<PetscScalar> (r);
852 
853  return *this;
854  }
855 
856 
857 
858  inline
859  VectorReference &
860  VectorReference::operator = (const VectorReference &r)
861  {
862  // as explained in the class
863  // documentation, this is not the copy
864  // operator. so simply pass on to the
865  // "correct" assignment operator
866  *this = static_cast<PetscScalar> (r);
867 
868  return *this;
869  }
870 
871 
872 
873  inline
874  const VectorReference &
875  VectorReference::operator = (const PetscScalar &value) const
876  {
877  Assert ((vector.last_action == VectorOperation::insert)
878  ||
879  (vector.last_action == VectorOperation::unknown),
880  ExcWrongMode (VectorOperation::insert,
881  vector.last_action));
882 
883  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
884 
885  const PetscInt petsc_i = index;
886 
887  const PetscErrorCode ierr = VecSetValues (vector, 1, &petsc_i, &value,
888  INSERT_VALUES);
889  AssertThrow (ierr == 0, ExcPETScError(ierr));
890 
891  vector.last_action = VectorOperation::insert;
892 
893  return *this;
894  }
895 
896 
897 
898  inline
899  const VectorReference &
900  VectorReference::operator += (const PetscScalar &value) const
901  {
902  Assert ((vector.last_action == VectorOperation::add)
903  ||
904  (vector.last_action == VectorOperation::unknown),
905  ExcWrongMode (VectorOperation::add,
906  vector.last_action));
907 
908  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
909 
910  vector.last_action = VectorOperation::add;
911 
912  // we have to do above actions in any
913  // case to be consistent with the MPI
914  // communication model (see the
915  // comments in the documentation of
916  // PETScWrappers::MPI::Vector), but we
917  // can save some work if the addend is
918  // zero
919  if (value == PetscScalar())
920  return *this;
921 
922  // use the PETSc function to add something
923  const PetscInt petsc_i = index;
924  const PetscErrorCode ierr = VecSetValues (vector, 1, &petsc_i, &value,
925  ADD_VALUES);
926  AssertThrow (ierr == 0, ExcPETScError(ierr));
927 
928 
929  return *this;
930  }
931 
932 
933 
934  inline
935  const VectorReference &
936  VectorReference::operator -= (const PetscScalar &value) const
937  {
938  Assert ((vector.last_action == VectorOperation::add)
939  ||
940  (vector.last_action == VectorOperation::unknown),
941  ExcWrongMode (VectorOperation::add,
942  vector.last_action));
943 
944  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
945 
946  vector.last_action = VectorOperation::add;
947 
948  // we have to do above actions in any
949  // case to be consistent with the MPI
950  // communication model (see the
951  // comments in the documentation of
952  // PETScWrappers::MPI::Vector), but we
953  // can save some work if the addend is
954  // zero
955  if (value == PetscScalar())
956  return *this;
957 
958  // use the PETSc function to
959  // add something
960  const PetscInt petsc_i = index;
961  const PetscScalar subtractand = -value;
962  const PetscErrorCode ierr = VecSetValues (vector, 1, &petsc_i, &subtractand,
963  ADD_VALUES);
964  AssertThrow (ierr == 0, ExcPETScError(ierr));
965 
966  return *this;
967  }
968 
969 
970 
971  inline
972  const VectorReference &
973  VectorReference::operator *= (const PetscScalar &value) const
974  {
975  Assert ((vector.last_action == VectorOperation::insert)
976  ||
977  (vector.last_action == VectorOperation::unknown),
978  ExcWrongMode (VectorOperation::insert,
979  vector.last_action));
980 
981  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
982 
983  vector.last_action = VectorOperation::insert;
984 
985  // we have to do above actions in any
986  // case to be consistent with the MPI
987  // communication model (see the
988  // comments in the documentation of
989  // PETScWrappers::MPI::Vector), but we
990  // can save some work if the factor is
991  // one
992  if (value == 1.)
993  return *this;
994 
995  const PetscInt petsc_i = index;
996  const PetscScalar new_value
997  = static_cast<PetscScalar>(*this) * value;
998 
999  const PetscErrorCode ierr = VecSetValues (vector, 1, &petsc_i, &new_value,
1000  INSERT_VALUES);
1001  AssertThrow (ierr == 0, ExcPETScError(ierr));
1002 
1003  return *this;
1004  }
1005 
1006 
1007 
1008  inline
1009  const VectorReference &
1010  VectorReference::operator /= (const PetscScalar &value) const
1011  {
1012  Assert ((vector.last_action == VectorOperation::insert)
1013  ||
1014  (vector.last_action == VectorOperation::unknown),
1015  ExcWrongMode (VectorOperation::insert,
1016  vector.last_action));
1017 
1018  Assert (!vector.has_ghost_elements(), ExcGhostsPresent());
1019 
1020  vector.last_action = VectorOperation::insert;
1021 
1022  // we have to do above actions in any
1023  // case to be consistent with the MPI
1024  // communication model (see the
1025  // comments in the documentation of
1026  // PETScWrappers::MPI::Vector), but we
1027  // can save some work if the factor is
1028  // one
1029  if (value == 1.)
1030  return *this;
1031 
1032  const PetscInt petsc_i = index;
1033  const PetscScalar new_value
1034  = static_cast<PetscScalar>(*this) / value;
1035 
1036  const PetscErrorCode ierr = VecSetValues (vector, 1, &petsc_i, &new_value,
1037  INSERT_VALUES);
1038  AssertThrow (ierr == 0, ExcPETScError(ierr));
1039 
1040  return *this;
1041  }
1042 
1043 
1044 
1045  inline
1046  PetscReal
1047  VectorReference::real () const
1048  {
1049 #ifndef PETSC_USE_COMPLEX
1050  return static_cast<PetscScalar>(*this);
1051 #else
1052  return PetscRealPart (static_cast<PetscScalar>(*this));
1053 #endif
1054  }
1055 
1056 
1057 
1058  inline
1059  PetscReal
1060  VectorReference::imag () const
1061  {
1062 #ifndef PETSC_USE_COMPLEX
1063  return PetscReal (0);
1064 #else
1065  return PetscImaginaryPart (static_cast<PetscScalar>(*this));
1066 #endif
1067  }
1068 
1069  } // namespace internal
1070 
1071  inline
1072  bool
1073  VectorBase::in_local_range (const size_type index) const
1074  {
1075  PetscInt begin, end;
1076  const PetscErrorCode ierr = VecGetOwnershipRange (static_cast<const Vec &>(vector),
1077  &begin, &end);
1078  AssertThrow (ierr == 0, ExcPETScError(ierr));
1079 
1080  return ((index >= static_cast<size_type>(begin)) &&
1081  (index < static_cast<size_type>(end)));
1082  }
1083 
1084 
1085  inline
1086  IndexSet
1087  VectorBase::locally_owned_elements() const
1088  {
1089  IndexSet is (size());
1090 
1091  // PETSc only allows for contiguous local ranges, so this is simple
1092  const std::pair<size_type, size_type> x = local_range();
1093  is.add_range (x.first, x.second);
1094  return is;
1095  }
1096 
1097 
1098 
1099  inline
1100  bool
1101  VectorBase::has_ghost_elements() const
1102  {
1103  return ghosted;
1104  }
1105 
1106 
1107 
1108  inline
1109  internal::VectorReference
1110  VectorBase::operator () (const size_type index)
1111  {
1112  return internal::VectorReference (*this, index);
1113  }
1114 
1115 
1116 
1117  inline
1118  PetscScalar
1119  VectorBase::operator () (const size_type index) const
1120  {
1121  return static_cast<PetscScalar>(internal::VectorReference (*this, index));
1122  }
1123 
1124 
1125 
1126  inline
1127  internal::VectorReference
1128  VectorBase::operator [] (const size_type index)
1129  {
1130  return operator()(index);
1131  }
1132 
1133 
1134 
1135  inline
1136  PetscScalar
1137  VectorBase::operator [] (const size_type index) const
1138  {
1139  return operator()(index);
1140  }
1141 
1142  inline
1143  const MPI_Comm &
1144  VectorBase::get_mpi_communicator () const
1145  {
1146  static MPI_Comm comm;
1147  PetscObjectGetComm((PetscObject)vector, &comm);
1148  return comm;
1149  }
1150 
1151  inline
1152  void VectorBase::extract_subvector_to (const std::vector<size_type> &indices,
1153  std::vector<PetscScalar> &values) const
1154  {
1155  extract_subvector_to(&(indices[0]), &(indices[0]) + indices.size(), &(values[0]));
1156  }
1157 
1158  template <typename ForwardIterator, typename OutputIterator>
1159  inline
1160  void VectorBase::extract_subvector_to (const ForwardIterator indices_begin,
1161  const ForwardIterator indices_end,
1162  OutputIterator values_begin) const
1163  {
1164  const PetscInt n_idx = static_cast<PetscInt>(indices_end - indices_begin);
1165  if (n_idx == 0)
1166  return;
1167 
1168  // if we are dealing
1169  // with a parallel vector
1170  if (ghosted )
1171  {
1172  // there is the possibility
1173  // that the vector has
1174  // ghost elements. in that
1175  // case, we first need to
1176  // figure out which
1177  // elements we own locally,
1178  // then get a pointer to
1179  // the elements that are
1180  // stored here (both the
1181  // ones we own as well as
1182  // the ghost elements). in
1183  // this array, the locally
1184  // owned elements come
1185  // first followed by the
1186  // ghost elements whose
1187  // position we can get from
1188  // an index set
1189  PetscInt begin, end;
1190  PetscErrorCode ierr = VecGetOwnershipRange (vector, &begin, &end);
1191  AssertThrow (ierr == 0, ExcPETScError(ierr));
1192 
1193  Vec locally_stored_elements = PETSC_NULL;
1194  ierr = VecGhostGetLocalForm(vector, &locally_stored_elements);
1195  AssertThrow (ierr == 0, ExcPETScError(ierr));
1196 
1197  PetscInt lsize;
1198  ierr = VecGetSize(locally_stored_elements, &lsize);
1199  AssertThrow (ierr == 0, ExcPETScError(ierr));
1200 
1201  PetscScalar *ptr;
1202  ierr = VecGetArray(locally_stored_elements, &ptr);
1203  AssertThrow (ierr == 0, ExcPETScError(ierr));
1204 
1205  for (PetscInt i=0; i<n_idx; ++i)
1206  {
1207  const unsigned int index = *(indices_begin+i);
1208  if ( index>=static_cast<unsigned int>(begin)
1209  && index<static_cast<unsigned int>(end) )
1210  {
1211  //local entry
1212  *(values_begin+i) = *(ptr+index-begin);
1213  }
1214  else
1215  {
1216  //ghost entry
1217  const unsigned int ghostidx
1218  = ghost_indices.index_within_set(index);
1219 
1220  Assert(ghostidx+end-begin<(unsigned int)lsize, ExcInternalError());
1221  *(values_begin+i) = *(ptr+ghostidx+end-begin);
1222  }
1223  }
1224 
1225  ierr = VecRestoreArray(locally_stored_elements, &ptr);
1226  AssertThrow (ierr == 0, ExcPETScError(ierr));
1227 
1228  ierr = VecGhostRestoreLocalForm(vector, &locally_stored_elements);
1229  AssertThrow (ierr == 0, ExcPETScError(ierr));
1230 
1231  }
1232  // if the vector is local or the
1233  // caller, then simply access the
1234  // element we are interested in
1235  else
1236  {
1237  PetscInt begin, end;
1238  PetscErrorCode ierr = VecGetOwnershipRange (vector, &begin, &end);
1239  AssertThrow (ierr == 0, ExcPETScError(ierr));
1240 
1241  PetscScalar *ptr;
1242  ierr = VecGetArray(vector, &ptr);
1243  AssertThrow (ierr == 0, ExcPETScError(ierr));
1244 
1245  for (PetscInt i=0; i<n_idx; ++i)
1246  {
1247  const unsigned int index = *(indices_begin+i);
1248 
1249  Assert(index>=static_cast<unsigned int>(begin)
1250  && index<static_cast<unsigned int>(end), ExcInternalError());
1251 
1252  *(values_begin+i) = *(ptr+index-begin);
1253  }
1254 
1255  ierr = VecRestoreArray(vector, &ptr);
1256  AssertThrow (ierr == 0, ExcPETScError(ierr));
1257 
1258  }
1259  }
1260 
1261 #endif // DOXYGEN
1262 }
1263 
1264 DEAL_II_NAMESPACE_CLOSE
1265 
1266 #endif // DEAL_II_WITH_PETSC
1267 
1268 /*---------------------------- petsc_vector_base.h ---------------------------*/
1269 
1270 #endif
1271 /*---------------------------- petsc_vector_base.h ---------------------------*/
reference operator[](const size_type index)
real_type normalize() const 1
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
void ratio(const VectorBase &a, const VectorBase &b) 1
void sadd(const PetscScalar s, const VectorBase &V)
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
VectorBase & operator=(const PetscScalar s)
bool in_local_range(const size_type index) const
void scale(const VectorBase &scaling_factors)
std::pair< size_type, size_type > local_range() const
real_type lp_norm(const real_type p) const
void compress(const VectorOperation::values operation)
unsigned int global_dof_index
Definition: types.h:88
void do_set_add_operation(const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
VectorOperation::values last_action
#define Assert(cond, exc)
Definition: exceptions.h:313
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
std::size_t memory_consumption() const
void equ(const PetscScalar a, const VectorBase &V)
friend class internal::VectorReference
static ::ExceptionBase & ExcGhostsPresent()
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
PetscScalar operator*(const VectorBase &vec) const
bool has_ghost_elements() const
PetscScalar mean_value() const
reference operator()(const size_type index)
bool operator==(const VectorBase &v) const
IndexSet locally_owned_elements() const
void add(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:588
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< PetscScalar > &values) const
virtual const MPI_Comm & get_mpi_communicator() const
static ::ExceptionBase & ExcInternalError()
bool operator!=(const VectorBase &v) const