Reference documentation for deal.II version 8.5.1
la_parallel_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 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__la_parallel_vector_h
17 #define dealii__la_parallel_vector_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/mpi.h>
21 #include <deal.II/base/numbers.h>
22 #include <deal.II/base/partitioner.h>
23 #include <deal.II/base/std_cxx11/shared_ptr.h>
24 #include <deal.II/base/thread_management.h>
25 #include <deal.II/lac/vector_space_vector.h>
26 #include <deal.II/lac/vector_type_traits.h>
27 
28 #include <iomanip>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace LinearAlgebra
33 {
34  namespace distributed
35  {
36  template <typename> class BlockVector;
37  }
38 }
39 
40 namespace LinearAlgebra
41 {
42  template <typename> class ReadWriteVector;
43 }
44 
45 #ifdef DEAL_II_WITH_PETSC
46 namespace PETScWrappers
47 {
48  namespace MPI
49  {
50  class Vector;
51  }
52 }
53 #endif
54 
55 #ifdef DEAL_II_WITH_TRILINOS
56 namespace TrilinosWrappers
57 {
58  namespace MPI
59  {
60  class Vector;
61  }
62 }
63 #endif
64 
65 namespace LinearAlgebra
66 {
67  namespace distributed
68  {
175  template <typename Number>
177  public Subscriptor
178  {
179  public:
180  typedef Number value_type;
181  typedef value_type *pointer;
182  typedef const value_type *const_pointer;
183  typedef value_type *iterator;
184  typedef const value_type *const_iterator;
185  typedef value_type &reference;
186  typedef const value_type &const_reference;
187  typedef types::global_dof_index size_type;
188  typedef typename numbers::NumberTraits<Number>::real_type real_type;
189 
203  static const bool supports_distributed_data DEAL_II_DEPRECATED = true;
204 
212  Vector ();
213 
217  Vector (const Vector<Number> &in_vector);
218 
223  Vector (const size_type size);
224 
241  Vector (const IndexSet &local_range,
242  const IndexSet &ghost_indices,
243  const MPI_Comm communicator);
244 
248  Vector (const IndexSet &local_range,
249  const MPI_Comm communicator);
250 
257  Vector (const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &partitioner);
258 
262  virtual ~Vector ();
263 
268  void reinit (const size_type size,
269  const bool omit_zeroing_entries = false);
270 
281  template <typename Number2>
282  void reinit(const Vector<Number2> &in_vector,
283  const bool omit_zeroing_entries = false);
284 
301  void reinit (const IndexSet &local_range,
302  const IndexSet &ghost_indices,
303  const MPI_Comm communicator);
304 
308  void reinit (const IndexSet &local_range,
309  const MPI_Comm communicator);
310 
317  void reinit (const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &partitioner);
318 
334  void swap (Vector<Number> &v);
335 
348  operator = (const Vector<Number> &in_vector);
349 
361  template <typename Number2>
363  operator = (const Vector<Number2> &in_vector);
364 
365 #ifdef DEAL_II_WITH_PETSC
366 
377  operator = (const PETScWrappers::MPI::Vector &petsc_vec) DEAL_II_DEPRECATED;
378 #endif
379 
380 #ifdef DEAL_II_WITH_TRILINOS
381 
393  operator = (const TrilinosWrappers::MPI::Vector &trilinos_vec) DEAL_II_DEPRECATED;
394 #endif
395 
396 
423  void compress (::VectorOperation::values operation);
424 
446  void update_ghost_values () const;
447 
462  void compress_start (const unsigned int communication_channel = 0,
464 
480  void compress_finish (::VectorOperation::values operation);
481 
496  void update_ghost_values_start (const unsigned int communication_channel = 0) const;
497 
498 
506  void update_ghost_values_finish () const;
507 
516  void zero_out_ghosts ();
517 
529  bool has_ghost_elements() const;
531 
536 
540  virtual Vector<Number> &operator*= (const Number factor);
541 
545  virtual Vector<Number> &operator/= (const Number factor);
546 
551 
556 
565  virtual void import(const LinearAlgebra::ReadWriteVector<Number> &V,
566  VectorOperation::values operation,
567  std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
568  std_cxx11::shared_ptr<const CommunicationPatternBase> ());
569 
573  virtual Number operator* (const VectorSpaceVector<Number> &V) const;
574 
578  virtual void add(const Number a);
579 
583  virtual void add(const Number a, const VectorSpaceVector<Number> &V);
584 
588  virtual void add(const Number a, const VectorSpaceVector<Number> &V,
589  const Number b, const VectorSpaceVector<Number> &W);
590 
595  virtual void add (const std::vector<size_type> &indices,
596  const std::vector<Number> &values);
597 
602  virtual void sadd(const Number s, const Number a,
603  const VectorSpaceVector<Number> &V);
604 
610  virtual void scale(const VectorSpaceVector<Number> &scaling_factors);
611 
615  virtual void equ(const Number a, const VectorSpaceVector<Number> &V);
616 
621  virtual real_type l1_norm() const;
622 
627  virtual real_type l2_norm() const;
628 
633  virtual real_type linfty_norm() const;
634 
651  virtual Number add_and_dot(const Number a,
652  const VectorSpaceVector<Number> &V,
653  const VectorSpaceVector<Number> &W);
654 
659  virtual size_type size() const;
660 
672  virtual ::IndexSet locally_owned_elements() const;
673 
677  virtual void print(std::ostream &out,
678  const unsigned int precision=3,
679  const bool scientific=true,
680  const bool across=true) const;
681 
685  virtual std::size_t memory_consumption() const;
687 
692 
698  Vector<Number> &operator = (const Number s);
699 
704  template <typename OtherNumber>
705  void add (const std::vector<size_type> &indices,
706  const ::Vector<OtherNumber> &values);
707 
712  template <typename OtherNumber>
713  void add (const size_type n_elements,
714  const size_type *indices,
715  const OtherNumber *values);
716 
721  void sadd (const Number s,
722  const Vector<Number> &V);
723 
729  void sadd (const Number s,
730  const Number a,
731  const Vector<Number> &V,
732  const Number b,
733  const Vector<Number> &W) DEAL_II_DEPRECATED;
734 
740  void equ (const Number a, const Vector<Number> &u,
741  const Number b, const Vector<Number> &v) DEAL_II_DEPRECATED;
742 
744 
745 
750 
755  size_type local_size() const;
756 
764  std::pair<size_type, size_type> local_range () const DEAL_II_DEPRECATED;
765 
772  bool in_local_range (const size_type global_index) const DEAL_II_DEPRECATED;
773 
779  size_type n_ghost_entries () const DEAL_II_DEPRECATED;
780 
788  const IndexSet &ghost_elements() const DEAL_II_DEPRECATED;
789 
797  bool is_ghost_entry (const types::global_dof_index global_index) const DEAL_II_DEPRECATED;
798 
806  iterator begin ();
807 
812  const_iterator begin () const;
813 
818  iterator end ();
819 
824  const_iterator end () const;
825 
835  Number operator () (const size_type global_index) const;
836 
846  Number &operator () (const size_type global_index);
847 
855  Number operator [] (const size_type global_index) const;
863  Number &operator [] (const size_type global_index);
864 
873  Number local_element (const size_type local_index) const;
874 
883  Number &local_element (const size_type local_index);
884 
891  template <typename OtherNumber>
892  void extract_subvector_to (const std::vector<size_type> &indices,
893  std::vector<OtherNumber> &values) const;
894 
899  template <typename ForwardIterator, typename OutputIterator>
900  void extract_subvector_to (ForwardIterator indices_begin,
901  const ForwardIterator indices_end,
902  OutputIterator values_begin) const;
908  bool all_zero () const;
909 
913  virtual Number mean_value () const;
914 
919  real_type lp_norm (const real_type p) const;
921 
926 
931  const MPI_Comm &get_mpi_communicator () const;
932 
939  const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &
940  get_partitioner () const;
941 
951  bool
952  partitioners_are_compatible (const Utilities::MPI::Partitioner &part) const;
953 
967  bool
968  partitioners_are_globally_compatible (const Utilities::MPI::Partitioner &part) const;
970 
977 
982  Number, Number, unsigned int,
983  << "Called compress(VectorOperation::insert), but"
984  << " the element received from a remote processor, value "
985  << std::setprecision(16) << arg1
986  << ", does not match with the value "
987  << std::setprecision(16) << arg2
988  << " on the owner processor " << arg3);
989 
994  size_type, size_type, size_type, size_type,
995  << "You tried to access element " << arg1
996  << " of a distributed vector, but this element is not "
997  << "stored on the current processor. Note: The range of "
998  << "locally owned elements is " << arg2 << " to "
999  << arg3 << ", and there are " << arg4 << " ghost elements "
1000  << "that this vector can access.");
1001 
1002  private:
1006  bool all_zero_local () const;
1007 
1011  template <typename Number2>
1012  Number inner_product_local (const Vector<Number2> &V) const;
1013 
1017  real_type norm_sqr_local () const;
1018 
1022  Number mean_value_local () const;
1023 
1027  real_type l1_norm_local () const;
1028 
1032  real_type lp_norm_local (const real_type p) const;
1033 
1037  real_type linfty_norm_local () const;
1038 
1043  Number add_and_dot_local (const Number a,
1044  const Vector<Number> &V,
1045  const Vector<Number> &W);
1046 
1052  std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> partitioner;
1053 
1057  size_type allocated_size;
1058 
1062  Number *val;
1063 
1068  mutable std_cxx11::shared_ptr< ::parallel::internal::TBBPartitioner> thread_loop_partitioner;
1069 
1075  mutable Number *import_data;
1076 
1084  mutable bool vector_is_ghosted;
1085 
1086 #ifdef DEAL_II_WITH_MPI
1087 
1095  std::vector<MPI_Request> compress_requests;
1096 
1101  mutable std::vector<MPI_Request> update_ghost_values_requests;
1102 #endif
1103 
1110 
1115  void clear_mpi_requests ();
1116 
1120  void resize_val (const size_type new_allocated_size);
1121 
1122  /*
1123  * Make all other vector types friends.
1124  */
1125  template <typename Number2> friend class Vector;
1126 
1130  template <typename Number2> friend class BlockVector;
1131  };
1135  /*----------------------- Inline functions ----------------------------------*/
1136 
1137 #ifndef DOXYGEN
1138 
1139 
1140  template <typename Number>
1141  inline
1142  bool
1144  {
1145  return vector_is_ghosted;
1146  }
1147 
1148 
1149 
1150  template <typename Number>
1151  inline
1152  typename Vector<Number>::size_type
1153  Vector<Number>::size () const
1154  {
1155  return partitioner->size();
1156  }
1157 
1158 
1159 
1160  template <typename Number>
1161  inline
1162  typename Vector<Number>::size_type
1164  {
1165  return partitioner->local_size();
1166  }
1167 
1168 
1169 
1170  template <typename Number>
1171  inline
1172  std::pair<typename Vector<Number>::size_type,
1173  typename Vector<Number>::size_type>
1175  {
1176  return partitioner->local_range();
1177  }
1178 
1179 
1180 
1181  template <typename Number>
1182  inline
1183  bool
1185  (const size_type global_index) const
1186  {
1187  return partitioner->in_local_range (global_index);
1188  }
1189 
1190 
1191 
1192  template <typename Number>
1193  inline
1194  IndexSet
1196  {
1197  IndexSet is (size());
1198 
1199  is.add_range (partitioner->local_range().first,
1200  partitioner->local_range().second);
1201 
1202  return is;
1203  }
1204 
1205 
1206 
1207  template <typename Number>
1208  inline
1209  typename Vector<Number>::size_type
1211  {
1212  return partitioner->n_ghost_indices();
1213  }
1214 
1215 
1216 
1217  template <typename Number>
1218  inline
1219  const IndexSet &
1221  {
1222  return partitioner->ghost_indices();
1223  }
1224 
1225 
1226 
1227  template <typename Number>
1228  inline
1229  bool
1230  Vector<Number>::is_ghost_entry (const size_type global_index) const
1231  {
1232  return partitioner->is_ghost_entry (global_index);
1233  }
1234 
1235 
1236 
1237  template <typename Number>
1238  inline
1239  typename Vector<Number>::iterator
1241  {
1242  return val;
1243  }
1244 
1245 
1246 
1247  template <typename Number>
1248  inline
1249  typename Vector<Number>::const_iterator
1250  Vector<Number>::begin () const
1251  {
1252  return val;
1253  }
1254 
1255 
1256 
1257  template <typename Number>
1258  inline
1259  typename Vector<Number>::iterator
1261  {
1262  return val+partitioner->local_size();
1263  }
1264 
1265 
1266 
1267  template <typename Number>
1268  inline
1269  typename Vector<Number>::const_iterator
1270  Vector<Number>::end () const
1271  {
1272  return val+partitioner->local_size();
1273  }
1274 
1275 
1276 
1277  template <typename Number>
1278  inline
1279  Number
1280  Vector<Number>::operator() (const size_type global_index) const
1281  {
1282  Assert (partitioner->in_local_range (global_index) ||
1283  partitioner->ghost_indices().is_element(global_index),
1284  ExcAccessToNonLocalElement(global_index, partitioner->local_range().first,
1285  partitioner->local_range().second,
1286  partitioner->ghost_indices().n_elements()));
1287  // do not allow reading a vector which is not in ghost mode
1288  Assert (partitioner->in_local_range (global_index) || vector_is_ghosted == true,
1289  ExcMessage("You tried to read a ghost element of this vector, "
1290  "but it has not imported its ghost values."));
1291  return val[partitioner->global_to_local(global_index)];
1292  }
1293 
1294 
1295 
1296  template <typename Number>
1297  inline
1298  Number &
1299  Vector<Number>::operator() (const size_type global_index)
1300  {
1301  Assert (partitioner->in_local_range (global_index) ||
1302  partitioner->ghost_indices().is_element(global_index),
1303  ExcAccessToNonLocalElement(global_index, partitioner->local_range().first,
1304  partitioner->local_range().second,
1305  partitioner->ghost_indices().n_elements()));
1306  // we would like to prevent reading ghosts from a vector that does not
1307  // have them imported, but this is not possible because we might be in a
1308  // part of the code where the vector has enabled ghosts but is non-const
1309  // (then, the compiler picks this method according to the C++ rule book
1310  // even if a human would pick the const method when this subsequent use
1311  // is just a read)
1312  return val[partitioner->global_to_local (global_index)];
1313  }
1314 
1315 
1316 
1317  template <typename Number>
1318  inline
1319  Number
1320  Vector<Number>::operator[] (const size_type global_index) const
1321  {
1322  return operator()(global_index);
1323  }
1324 
1325 
1326 
1327  template <typename Number>
1328  inline
1329  Number &
1330  Vector<Number>::operator[] (const size_type global_index)
1331  {
1332  return operator()(global_index);
1333  }
1334 
1335 
1336 
1337  template <typename Number>
1338  inline
1339  Number
1340  Vector<Number>::local_element (const size_type local_index) const
1341  {
1342  AssertIndexRange (local_index,
1343  partitioner->local_size()+
1344  partitioner->n_ghost_indices());
1345  // do not allow reading a vector which is not in ghost mode
1346  Assert (local_index < local_size() || vector_is_ghosted == true,
1347  ExcMessage("You tried to read a ghost element of this vector, "
1348  "but it has not imported its ghost values."));
1349  return val[local_index];
1350  }
1351 
1352 
1353 
1354  template <typename Number>
1355  inline
1356  Number &
1357  Vector<Number>::local_element (const size_type local_index)
1358  {
1359  AssertIndexRange (local_index,
1360  partitioner->local_size()+
1361  partitioner->n_ghost_indices());
1362  return val[local_index];
1363  }
1364 
1365 
1366 
1367  template <typename Number>
1368  template <typename OtherNumber>
1369  inline
1370  void Vector<Number>::extract_subvector_to (const std::vector<size_type> &indices,
1371  std::vector<OtherNumber> &values) const
1372  {
1373  for (size_type i = 0; i < indices.size(); ++i)
1374  values[i] = operator()(indices[i]);
1375  }
1376 
1377 
1378 
1379  template <typename Number>
1380  template <typename ForwardIterator, typename OutputIterator>
1381  inline
1382  void Vector<Number>::extract_subvector_to (ForwardIterator indices_begin,
1383  const ForwardIterator indices_end,
1384  OutputIterator values_begin) const
1385  {
1386  while (indices_begin != indices_end)
1387  {
1388  *values_begin = operator()(*indices_begin);
1389  indices_begin++;
1390  values_begin++;
1391  }
1392  }
1393 
1394 
1395 
1396  template <typename Number>
1397  template <typename OtherNumber>
1398  inline
1399  void
1400  Vector<Number>::add (const std::vector<size_type> &indices,
1401  const ::Vector<OtherNumber> &values)
1402  {
1403  AssertDimension (indices.size(), values.size());
1404  for (size_type i=0; i<indices.size(); ++i)
1405  {
1406  Assert (numbers::is_finite(values[i]),
1407  ExcMessage("The given value is not finite but either infinite or Not A Number (NaN)"));
1408  this->operator()(indices[i]) += values(i);
1409  }
1410  }
1411 
1412 
1413 
1414  template <typename Number>
1415  template <typename OtherNumber>
1416  inline
1417  void
1418  Vector<Number>::add (const size_type n_elements,
1419  const size_type *indices,
1420  const OtherNumber *values)
1421  {
1422  for (size_type i=0; i<n_elements; ++i, ++indices, ++values)
1423  {
1424  Assert (numbers::is_finite(*values),
1425  ExcMessage("The given value is not finite but either infinite or Not A Number (NaN)"));
1426  this->operator()(*indices) += *values;
1427  }
1428  }
1429 
1430 
1431 
1432  template <typename Number>
1433  inline
1434  const MPI_Comm &
1436  {
1437  return partitioner->get_mpi_communicator();
1438  }
1439 
1440 
1441 
1442  template <typename Number>
1443  inline
1444  const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &
1446  {
1447  return partitioner;
1448  }
1449 
1450 #endif
1451 
1452  }
1453 }
1454 
1455 
1464 template <typename Number>
1465 inline
1468 {
1469  u.swap (v);
1470 }
1471 
1472 
1478 template <typename Number>
1479 struct is_serial_vector< LinearAlgebra::distributed::Vector< Number > > : std_cxx11::false_type
1480 {
1481 };
1482 
1483 
1484 DEAL_II_NAMESPACE_CLOSE
1485 
1486 #endif
virtual Number add_and_dot(const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W)
void compress(::VectorOperation::values operation)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
virtual Number mean_value() const
std::vector< MPI_Request > compress_requests
virtual Vector< Number > & operator*=(const Number factor)
virtual void equ(const Number a, const VectorSpaceVector< Number > &V)
Number inner_product_local(const Vector< Number2 > &V) const
bool in_local_range(const size_type global_index) const 1
virtual size_type size() const
virtual void scale(const VectorSpaceVector< Number > &scaling_factors)
real_type lp_norm(const real_type p) const
void swap(LinearAlgebra::distributed::Vector< Number > &u, LinearAlgebra::distributed::Vector< Number > &v)
virtual std::size_t memory_consumption() const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
bool is_ghost_entry(const types::global_dof_index global_index) const 1
STL namespace.
std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > partitioner
virtual ::IndexSet locally_owned_elements() const
static ::ExceptionBase & ExcVectorTypeNotCompatible()
bool is_finite(const double x)
Definition: numbers.h:275
virtual void sadd(const Number s, const Number a, const VectorSpaceVector< Number > &V)
std::vector< MPI_Request > update_ghost_values_requests
virtual Vector< Number > & operator/=(const Number factor)
void swap(Vector< Number > &v)
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual void add(const Number a)
virtual real_type l2_norm() const
Definition: types.h:30
unsigned int global_dof_index
Definition: types.h:88
virtual Vector< Number > & operator+=(const VectorSpaceVector< Number > &V)
#define Assert(cond, exc)
Definition: exceptions.h:313
void resize_val(const size_type new_allocated_size)
const IndexSet & ghost_elements() const 1
virtual Number operator*(const VectorSpaceVector< Number > &V) const
Number operator[](const size_type global_index) const
#define DeclException0(Exception0)
Definition: exceptions.h:541
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
virtual real_type l1_norm() const
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< OtherNumber > &values) const
const MPI_Comm & get_mpi_communicator() const
Number add_and_dot_local(const Number a, const Vector< Number > &V, const Vector< Number > &W)
Vector< Number > & operator=(const Vector< Number > &in_vector)
void compress_start(const unsigned int communication_channel=0, ::VectorOperation::values operation=VectorOperation::add)
std_cxx11::shared_ptr< ::parallel::internal::TBBPartitioner > thread_loop_partitioner
std::pair< size_type, size_type > local_range() const 1
Definition: mpi.h:48
void compress_finish(::VectorOperation::values operation)
virtual real_type linfty_norm() const
Number operator()(const size_type global_index) const
bool partitioners_are_globally_compatible(const Utilities::MPI::Partitioner &part) const
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:600
static ::ExceptionBase & ExcNonMatchingElements(Number arg1, Number arg2, unsigned int arg3)
value_type operator()(const size_type i) const
bool partitioners_are_compatible(const Utilities::MPI::Partitioner &part) const
real_type lp_norm_local(const real_type p) const
Number local_element(const size_type local_index) const
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:588
const std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
size_type n_ghost_entries() const 1
virtual Vector< Number > & operator-=(const VectorSpaceVector< Number > &V)
void update_ghost_values_start(const unsigned int communication_channel=0) const