Reference documentation for deal.II version 8.5.1
trilinos_vector_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 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__trilinos_vector_base_h
17 #define dealii__trilinos_vector_base_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
24 #include <deal.II/base/utilities.h>
25 # include <deal.II/base/std_cxx11/shared_ptr.h>
26 # include <deal.II/base/subscriptor.h>
27 # include <deal.II/lac/exceptions.h>
28 # include <deal.II/lac/vector.h>
29 # include <deal.II/base/mpi.h>
30 
31 # include <vector>
32 # include <utility>
33 # include <memory>
34 
35 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
36 # include "Epetra_ConfigDefs.h"
37 # ifdef DEAL_II_WITH_MPI // only if MPI is installed
38 # include "mpi.h"
39 # include "Epetra_MpiComm.h"
40 # else
41 # include "Epetra_SerialComm.h"
42 # endif
43 # include "Epetra_FEVector.h"
44 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
45 
46 DEAL_II_NAMESPACE_OPEN
47 
48 // forward declaration
49 template <typename number> class Vector;
50 
51 
62 namespace TrilinosWrappers
63 {
64  // forward declaration
65  class VectorBase;
66 
77  namespace internal
78  {
82  typedef ::types::global_dof_index size_type;
83 
93  class VectorReference
94  {
95  private:
100  VectorReference (VectorBase &vector,
101  const size_type index);
102 
103  public:
104 
116  const VectorReference &
117  operator = (const VectorReference &r) const;
118 
122  VectorReference &
123  operator = (const VectorReference &r);
124 
128  const VectorReference &
129  operator = (const TrilinosScalar &s) const;
130 
134  const VectorReference &
135  operator += (const TrilinosScalar &s) const;
136 
140  const VectorReference &
141  operator -= (const TrilinosScalar &s) const;
142 
146  const VectorReference &
147  operator *= (const TrilinosScalar &s) const;
148 
152  const VectorReference &
153  operator /= (const TrilinosScalar &s) const;
154 
159  operator TrilinosScalar () const;
160 
165  int,
166  << "An error with error number " << arg1
167  << " occurred while calling a Trilinos function");
168 
169  private:
173  VectorBase &vector;
174 
178  const size_type index;
179 
184  friend class ::TrilinosWrappers::VectorBase;
185  };
186  }
222  class VectorBase : public Subscriptor
223  {
224  public:
230  typedef TrilinosScalar value_type;
231  typedef TrilinosScalar real_type;
232  typedef ::types::global_dof_index size_type;
233  typedef value_type *iterator;
234  typedef const value_type *const_iterator;
235  typedef internal::VectorReference reference;
236  typedef const internal::VectorReference const_reference;
237 
242 
248  VectorBase ();
249 
254  VectorBase (const VectorBase &v);
255 
259  virtual ~VectorBase ();
260 
265  void clear ();
266 
272  void reinit (const VectorBase &v,
273  const bool omit_zeroing_entries = false);
274 
291  void compress (::VectorOperation::values operation);
292 
299  bool is_compressed () const DEAL_II_DEPRECATED;
300 
313  VectorBase &
314  operator = (const TrilinosScalar s);
315 
321  VectorBase &
322  operator = (const VectorBase &v);
323 
332  template <typename Number>
333  VectorBase &
334  operator = (const ::Vector<Number> &v);
335 
341  bool operator == (const VectorBase &v) const;
342 
348  bool operator != (const VectorBase &v) const;
349 
353  size_type size () const;
354 
366  size_type local_size () const;
367 
389  std::pair<size_type, size_type> local_range () const;
390 
398  bool in_local_range (const size_type index) const;
399 
414 
422  bool has_ghost_elements() const;
423 
428  TrilinosScalar operator * (const VectorBase &vec) const;
429 
433  real_type norm_sqr () const;
434 
438  TrilinosScalar mean_value () const;
439 
445  TrilinosScalar minimal_value () const DEAL_II_DEPRECATED;
446 
450  TrilinosScalar min () const;
451 
455  TrilinosScalar max () const;
456 
460  real_type l1_norm () const;
461 
466  real_type l2_norm () const;
467 
472  real_type lp_norm (const TrilinosScalar p) const;
473 
477  real_type linfty_norm () const;
478 
494  TrilinosScalar add_and_dot (const TrilinosScalar a,
495  const VectorBase &V,
496  const VectorBase &W);
497 
503  bool all_zero () const;
504 
510  bool is_non_negative () const;
512 
513 
518 
528  reference
529  operator () (const size_type index);
530 
540  TrilinosScalar
541  operator () (const size_type index) const;
542 
548  reference
549  operator [] (const size_type index);
550 
556  TrilinosScalar
557  operator [] (const size_type index) const;
558 
565  void extract_subvector_to (const std::vector<size_type> &indices,
566  std::vector<TrilinosScalar> &values) const;
567 
572  template <typename ForwardIterator, typename OutputIterator>
573  void extract_subvector_to (ForwardIterator indices_begin,
574  const ForwardIterator indices_end,
575  OutputIterator values_begin) const;
576 
587  TrilinosScalar el (const size_type index) const DEAL_II_DEPRECATED;
588 
599  iterator begin ();
600 
605  const_iterator begin () const;
606 
611  iterator end ();
612 
617  const_iterator end () const;
618 
620 
621 
626 
633  void set (const std::vector<size_type> &indices,
634  const std::vector<TrilinosScalar> &values);
635 
640  void set (const std::vector<size_type> &indices,
641  const ::Vector<TrilinosScalar> &values);
642 
648  void set (const size_type n_elements,
649  const size_type *indices,
650  const TrilinosScalar *values);
651 
656  void add (const std::vector<size_type> &indices,
657  const std::vector<TrilinosScalar> &values);
658 
663  void add (const std::vector<size_type> &indices,
664  const ::Vector<TrilinosScalar> &values);
665 
671  void add (const size_type n_elements,
672  const size_type *indices,
673  const TrilinosScalar *values);
674 
678  VectorBase &operator *= (const TrilinosScalar factor);
679 
683  VectorBase &operator /= (const TrilinosScalar factor);
684 
688  VectorBase &operator += (const VectorBase &V);
689 
693  VectorBase &operator -= (const VectorBase &V);
694 
699  void add (const TrilinosScalar s);
700 
713  void add (const VectorBase &V,
714  const bool allow_different_maps = false);
715 
719  void add (const TrilinosScalar a,
720  const VectorBase &V);
721 
725  void add (const TrilinosScalar a,
726  const VectorBase &V,
727  const TrilinosScalar b,
728  const VectorBase &W);
729 
734  void sadd (const TrilinosScalar s,
735  const VectorBase &V);
736 
740  void sadd (const TrilinosScalar s,
741  const TrilinosScalar a,
742  const VectorBase &V);
743 
749  void sadd (const TrilinosScalar s,
750  const TrilinosScalar a,
751  const VectorBase &V,
752  const TrilinosScalar b,
753  const VectorBase &W) DEAL_II_DEPRECATED;
754 
761  void sadd (const TrilinosScalar s,
762  const TrilinosScalar a,
763  const VectorBase &V,
764  const TrilinosScalar b,
765  const VectorBase &W,
766  const TrilinosScalar c,
767  const VectorBase &X) DEAL_II_DEPRECATED;
768 
774  void scale (const VectorBase &scaling_factors);
775 
779  void equ (const TrilinosScalar a,
780  const VectorBase &V);
781 
787  void equ (const TrilinosScalar a,
788  const VectorBase &V,
789  const TrilinosScalar b,
790  const VectorBase &W) DEAL_II_DEPRECATED;
791 
802  void ratio (const VectorBase &a,
803  const VectorBase &b) DEAL_II_DEPRECATED;
805 
806 
811 
816  const Epetra_MultiVector &trilinos_vector () const;
817 
822  Epetra_FEVector &trilinos_vector ();
823 
828  const Epetra_Map &vector_partitioner () const;
829 
836  void print (const char *format = 0) const DEAL_II_DEPRECATED;
837 
845  void print (std::ostream &out,
846  const unsigned int precision = 3,
847  const bool scientific = true,
848  const bool across = true) const;
849 
863  void swap (VectorBase &v);
864 
868  std::size_t memory_consumption () const;
869 
874  const MPI_Comm &get_mpi_communicator () const;
876 
881 
886  int,
887  << "An error with error number " << arg1
888  << " occurred while calling a Trilinos function");
889 
894  size_type, size_type, size_type, size_type,
895  << "You tried to access element " << arg1
896  << " of a distributed vector, but this element is not stored "
897  << "on the current processor. Note: There are "
898  << arg2 << " elements stored "
899  << "on the current processor from within the range "
900  << arg3 << " through " << arg4
901  << " but Trilinos vectors need not store contiguous "
902  << "ranges on each processor, and not every element in "
903  << "this range may in fact be stored locally.");
904 
905 
906  private:
918  Epetra_CombineMode last_action;
919 
925 
931 
937  std_cxx11::shared_ptr<Epetra_FEVector> vector;
938 
944  std_cxx11::shared_ptr<Epetra_MultiVector> nonlocal_vector;
945 
950 
954  friend class internal::VectorReference;
955  friend class Vector;
956  friend class MPI::Vector;
957  };
958 
959 
960 
961 
962 // ------------------- inline and template functions --------------
963 
972  inline
973  void swap (VectorBase &u, VectorBase &v)
974  {
975  u.swap (v);
976  }
977 
978 
979 #ifndef DOXYGEN
980 
981  namespace internal
982  {
983  inline
984  VectorReference::VectorReference (VectorBase &vector,
985  const size_type index)
986  :
987  vector (vector),
988  index (index)
989  {}
990 
991 
992  inline
993  const VectorReference &
994  VectorReference::operator = (const VectorReference &r) const
995  {
996  // as explained in the class
997  // documentation, this is not the copy
998  // operator. so simply pass on to the
999  // "correct" assignment operator
1000  *this = static_cast<TrilinosScalar> (r);
1001 
1002  return *this;
1003  }
1004 
1005 
1006 
1007  inline
1008  VectorReference &
1009  VectorReference::operator = (const VectorReference &r)
1010  {
1011  // as above
1012  *this = static_cast<TrilinosScalar> (r);
1013 
1014  return *this;
1015  }
1016 
1017 
1018  inline
1019  const VectorReference &
1020  VectorReference::operator = (const TrilinosScalar &value) const
1021  {
1022  vector.set (1, &index, &value);
1023  return *this;
1024  }
1025 
1026 
1027 
1028  inline
1029  const VectorReference &
1030  VectorReference::operator += (const TrilinosScalar &value) const
1031  {
1032  vector.add (1, &index, &value);
1033  return *this;
1034  }
1035 
1036 
1037 
1038  inline
1039  const VectorReference &
1040  VectorReference::operator -= (const TrilinosScalar &value) const
1041  {
1042  TrilinosScalar new_value = -value;
1043  vector.add (1, &index, &new_value);
1044  return *this;
1045  }
1046 
1047 
1048 
1049  inline
1050  const VectorReference &
1051  VectorReference::operator *= (const TrilinosScalar &value) const
1052  {
1053  TrilinosScalar new_value = static_cast<TrilinosScalar>(*this) * value;
1054  vector.set (1, &index, &new_value);
1055  return *this;
1056  }
1057 
1058 
1059 
1060  inline
1061  const VectorReference &
1062  VectorReference::operator /= (const TrilinosScalar &value) const
1063  {
1064  TrilinosScalar new_value = static_cast<TrilinosScalar>(*this) / value;
1065  vector.set (1, &index, &new_value);
1066  return *this;
1067  }
1068  }
1069 
1070 
1071 
1072  inline
1073  bool
1074  VectorBase::is_compressed () const
1075  {
1076  return compressed;
1077  }
1078 
1079 
1080 
1081  inline
1082  bool
1083  VectorBase::in_local_range (const size_type index) const
1084  {
1085  std::pair<size_type, size_type> range = local_range();
1086 
1087  return ((index >= range.first) && (index < range.second));
1088  }
1089 
1090 
1091 
1092  inline
1093  IndexSet
1094  VectorBase::locally_owned_elements() const
1095  {
1096  Assert(owned_elements.size()==size(),
1097  ExcMessage("The locally owned elements have not been properly initialized!"
1098  " This happens for example if this object has been initialized"
1099  " with exactly one overlapping IndexSet."));
1100  return owned_elements;
1101  }
1102 
1103 
1104 
1105  inline
1106  bool
1107  VectorBase::has_ghost_elements() const
1108  {
1109  return has_ghosts;
1110  }
1111 
1112 
1113 
1114  inline
1115  internal::VectorReference
1116  VectorBase::operator () (const size_type index)
1117  {
1118  return internal::VectorReference (*this, index);
1119  }
1120 
1121 
1122 
1123  inline
1124  internal::VectorReference
1125  VectorBase::operator [] (const size_type index)
1126  {
1127  return operator() (index);
1128  }
1129 
1130 
1131  inline
1132  TrilinosScalar
1133  VectorBase::operator [] (const size_type index) const
1134  {
1135  return operator() (index);
1136  }
1137 
1138 
1139 
1140  inline
1141  void VectorBase::extract_subvector_to (const std::vector<size_type> &indices,
1142  std::vector<TrilinosScalar> &values) const
1143  {
1144  for (size_type i = 0; i < indices.size(); ++i)
1145  values[i] = operator()(indices[i]);
1146  }
1147 
1148 
1149 
1150  template <typename ForwardIterator, typename OutputIterator>
1151  inline
1152  void VectorBase::extract_subvector_to (ForwardIterator indices_begin,
1153  const ForwardIterator indices_end,
1154  OutputIterator values_begin) const
1155  {
1156  while (indices_begin != indices_end)
1157  {
1158  *values_begin = operator()(*indices_begin);
1159  indices_begin++;
1160  values_begin++;
1161  }
1162  }
1163 
1164 
1165 
1166  inline
1167  VectorBase::iterator
1168  VectorBase::begin()
1169  {
1170  return (*vector)[0];
1171  }
1172 
1173 
1174 
1175  inline
1176  VectorBase::iterator
1177  VectorBase::end()
1178  {
1179  return (*vector)[0]+local_size();
1180  }
1181 
1182 
1183 
1184  inline
1185  VectorBase::const_iterator
1186  VectorBase::begin() const
1187  {
1188  return (*vector)[0];
1189  }
1190 
1191 
1192 
1193  inline
1194  VectorBase::const_iterator
1195  VectorBase::end() const
1196  {
1197  return (*vector)[0]+local_size();
1198  }
1199 
1200 
1201 
1202  inline
1203  void
1204  VectorBase::reinit (const VectorBase &v,
1205  const bool omit_zeroing_entries)
1206  {
1207  Assert (vector.get() != 0,
1208  ExcMessage("Vector has not been constructed properly."));
1209 
1210  if (omit_zeroing_entries == false ||
1211  vector_partitioner().SameAs(v.vector_partitioner())==false)
1212  vector.reset (new Epetra_FEVector(*v.vector));
1213 
1214  if (v.nonlocal_vector.get() != 0)
1215  nonlocal_vector.reset(new Epetra_MultiVector(v.nonlocal_vector->Map(), 1));
1216 
1217  owned_elements = v.owned_elements;
1218  }
1219 
1220 
1221 
1222  inline
1223  VectorBase &
1224  VectorBase::operator = (const TrilinosScalar s)
1225  {
1226  AssertIsFinite(s);
1227 
1228  int ierr = vector->PutScalar(s);
1229  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1230 
1231  if (nonlocal_vector.get() != 0)
1232  {
1233  ierr = nonlocal_vector->PutScalar(0.);
1234  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1235  }
1236 
1237  return *this;
1238  }
1239 
1240 
1241 
1242  inline
1243  void
1244  VectorBase::set (const std::vector<size_type> &indices,
1245  const std::vector<TrilinosScalar> &values)
1246  {
1247  // if we have ghost values, do not allow
1248  // writing to this vector at all.
1249  Assert (!has_ghost_elements(), ExcGhostsPresent());
1250 
1251  Assert (indices.size() == values.size(),
1252  ExcDimensionMismatch(indices.size(),values.size()));
1253 
1254  set (indices.size(), &indices[0], &values[0]);
1255  }
1256 
1257 
1258 
1259  inline
1260  void
1261  VectorBase::set (const std::vector<size_type> &indices,
1262  const ::Vector<TrilinosScalar> &values)
1263  {
1264  // if we have ghost values, do not allow
1265  // writing to this vector at all.
1266  Assert (!has_ghost_elements(), ExcGhostsPresent());
1267 
1268  Assert (indices.size() == values.size(),
1269  ExcDimensionMismatch(indices.size(),values.size()));
1270 
1271  set (indices.size(), &indices[0], values.begin());
1272  }
1273 
1274 
1275 
1276  inline
1277  void
1278  VectorBase::set (const size_type n_elements,
1279  const size_type *indices,
1280  const TrilinosScalar *values)
1281  {
1282  // if we have ghost values, do not allow
1283  // writing to this vector at all.
1284  Assert (!has_ghost_elements(), ExcGhostsPresent());
1285 
1286  if (last_action == Add)
1287  {
1288  const int ierr = vector->GlobalAssemble(Add);
1289  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1290  }
1291 
1292  if (last_action != Insert)
1293  last_action = Insert;
1294 
1295  for (size_type i=0; i<n_elements; ++i)
1296  {
1297  const size_type row = indices[i];
1298  const TrilinosWrappers::types::int_type local_row = vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(row));
1299  if (local_row != -1)
1300  (*vector)[0][local_row] = values[i];
1301  else
1302  {
1303  const int ierr = vector->ReplaceGlobalValues (1,
1304  (const TrilinosWrappers::types::int_type *)(&row),
1305  &values[i]);
1306  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1307  compressed = false;
1308  }
1309  // in set operation, do not use the pre-allocated vector for nonlocal
1310  // entries even if it exists. This is to ensure that we really only
1311  // set the elements touched by the set() method and not all contained
1312  // in the nonlocal entries vector (there is no way to distinguish them
1313  // on the receiving processor)
1314  }
1315  }
1316 
1317 
1318 
1319  inline
1320  void
1321  VectorBase::add (const std::vector<size_type> &indices,
1322  const std::vector<TrilinosScalar> &values)
1323  {
1324  // if we have ghost values, do not allow
1325  // writing to this vector at all.
1326  Assert (!has_ghost_elements(), ExcGhostsPresent());
1327  Assert (indices.size() == values.size(),
1328  ExcDimensionMismatch(indices.size(),values.size()));
1329 
1330  add (indices.size(), &indices[0], &values[0]);
1331  }
1332 
1333 
1334 
1335  inline
1336  void
1337  VectorBase::add (const std::vector<size_type> &indices,
1338  const ::Vector<TrilinosScalar> &values)
1339  {
1340  // if we have ghost values, do not allow
1341  // writing to this vector at all.
1342  Assert (!has_ghost_elements(), ExcGhostsPresent());
1343  Assert (indices.size() == values.size(),
1344  ExcDimensionMismatch(indices.size(),values.size()));
1345 
1346  add (indices.size(), &indices[0], values.begin());
1347  }
1348 
1349 
1350 
1351  inline
1352  void
1353  VectorBase::add (const size_type n_elements,
1354  const size_type *indices,
1355  const TrilinosScalar *values)
1356  {
1357  // if we have ghost values, do not allow
1358  // writing to this vector at all.
1359  Assert (!has_ghost_elements(), ExcGhostsPresent());
1360 
1361  if (last_action != Add)
1362  {
1363  if (last_action == Insert)
1364  {
1365  const int ierr = vector->GlobalAssemble(Insert);
1366  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1367  }
1368  last_action = Add;
1369  }
1370 
1371  for (size_type i=0; i<n_elements; ++i)
1372  {
1373  const size_type row = indices[i];
1374  const TrilinosWrappers::types::int_type local_row = vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(row));
1375  if (local_row != -1)
1376  (*vector)[0][local_row] += values[i];
1377  else if (nonlocal_vector.get() == 0)
1378  {
1379  const int ierr = vector->SumIntoGlobalValues (1,
1380  (const TrilinosWrappers::types::int_type *)(&row),
1381  &values[i]);
1382  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1383  compressed = false;
1384  }
1385  else
1386  {
1387  // use pre-allocated vector for non-local entries if it exists for
1388  // addition operation
1389  const TrilinosWrappers::types::int_type my_row = nonlocal_vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(row));
1390  Assert(my_row != -1,
1391  ExcMessage("Attempted to write into off-processor vector entry "
1392  "that has not be specified as being writable upon "
1393  "initialization"));
1394  (*nonlocal_vector)[0][my_row] += values[i];
1395  compressed = false;
1396  }
1397  }
1398  }
1399 
1400 
1401 
1402  inline
1403  VectorBase::size_type
1404  VectorBase::size () const
1405  {
1406 #ifndef DEAL_II_WITH_64BIT_INDICES
1407  return (size_type) (vector->Map().MaxAllGID() + 1 - vector->Map().MinAllGID());
1408 #else
1409  return (size_type) (vector->Map().MaxAllGID64() + 1 - vector->Map().MinAllGID64());
1410 #endif
1411  }
1412 
1413 
1414 
1415  inline
1416  VectorBase::size_type
1417  VectorBase::local_size () const
1418  {
1419  return (size_type) vector->Map().NumMyElements();
1420  }
1421 
1422 
1423 
1424  inline
1425  std::pair<VectorBase::size_type, VectorBase::size_type>
1426  VectorBase::local_range () const
1427  {
1428 #ifndef DEAL_II_WITH_64BIT_INDICES
1429  const TrilinosWrappers::types::int_type begin = vector->Map().MinMyGID();
1430  const TrilinosWrappers::types::int_type end = vector->Map().MaxMyGID()+1;
1431 #else
1432  const TrilinosWrappers::types::int_type begin = vector->Map().MinMyGID64();
1433  const TrilinosWrappers::types::int_type end = vector->Map().MaxMyGID64()+1;
1434 #endif
1435 
1436  Assert (end-begin == vector->Map().NumMyElements(),
1437  ExcMessage ("This function only makes sense if the elements that this "
1438  "vector stores on the current processor form a contiguous range. "
1439  "This does not appear to be the case for the current vector."));
1440 
1441  return std::make_pair (begin, end);
1442  }
1443 
1444 
1445 
1446  inline
1447  TrilinosScalar
1448  VectorBase::operator * (const VectorBase &vec) const
1449  {
1450  Assert (vector->Map().SameAs(vec.vector->Map()),
1451  ExcDifferentParallelPartitioning());
1452  Assert (!has_ghost_elements(), ExcGhostsPresent());
1453 
1454  TrilinosScalar result;
1455 
1456  const int ierr = vector->Dot(*(vec.vector), &result);
1457  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1458 
1459  return result;
1460  }
1461 
1462 
1463 
1464  inline
1465  VectorBase::real_type
1466  VectorBase::norm_sqr () const
1467  {
1468  const TrilinosScalar d = l2_norm();
1469  return d*d;
1470  }
1471 
1472 
1473 
1474  inline
1475  TrilinosScalar
1476  VectorBase::mean_value () const
1477  {
1478  Assert (!has_ghost_elements(), ExcGhostsPresent());
1479 
1480  TrilinosScalar mean;
1481  const int ierr = vector->MeanValue (&mean);
1482  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1483 
1484  return mean;
1485  }
1486 
1487 
1488 
1489  inline
1490  TrilinosScalar
1491  VectorBase::minimal_value () const
1492  {
1493  return min();
1494  }
1495 
1496 
1497 
1498  inline
1499  TrilinosScalar
1500  VectorBase::min () const
1501  {
1502  TrilinosScalar min_value;
1503  const int ierr = vector->MinValue (&min_value);
1504  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1505 
1506  return min_value;
1507  }
1508 
1509 
1510 
1511  inline
1512  TrilinosScalar
1513  VectorBase::max () const
1514  {
1515  TrilinosScalar max_value;
1516  const int ierr = vector->MaxValue (&max_value);
1517  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1518 
1519  return max_value;
1520  }
1521 
1522 
1523 
1524  inline
1525  VectorBase::real_type
1526  VectorBase::l1_norm () const
1527  {
1528  Assert (!has_ghost_elements(), ExcGhostsPresent());
1529 
1530  TrilinosScalar d;
1531  const int ierr = vector->Norm1 (&d);
1532  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1533 
1534  return d;
1535  }
1536 
1537 
1538 
1539  inline
1540  VectorBase::real_type
1541  VectorBase::l2_norm () const
1542  {
1543  Assert (!has_ghost_elements(), ExcGhostsPresent());
1544 
1545  TrilinosScalar d;
1546  const int ierr = vector->Norm2 (&d);
1547  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1548 
1549  return d;
1550  }
1551 
1552 
1553 
1554  inline
1555  VectorBase::real_type
1556  VectorBase::lp_norm (const TrilinosScalar p) const
1557  {
1558  Assert (!has_ghost_elements(), ExcGhostsPresent());
1559 
1560  TrilinosScalar norm = 0;
1561  TrilinosScalar sum=0;
1562  const size_type n_local = local_size();
1563 
1564  // loop over all the elements because
1565  // Trilinos does not support lp norms
1566  for (size_type i=0; i<n_local; i++)
1567  sum += std::pow(std::fabs((*vector)[0][i]), p);
1568 
1569  norm = std::pow(sum, static_cast<TrilinosScalar>(1./p));
1570 
1571  return norm;
1572  }
1573 
1574 
1575 
1576  inline
1577  VectorBase::real_type
1578  VectorBase::linfty_norm () const
1579  {
1580  // while we disallow the other
1581  // norm operations on ghosted
1582  // vectors, this particular norm
1583  // is safe to run even in the
1584  // presence of ghost elements
1585  TrilinosScalar d;
1586  const int ierr = vector->NormInf (&d);
1587  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1588 
1589  return d;
1590  }
1591 
1592 
1593 
1594  inline
1595  TrilinosScalar
1596  VectorBase::add_and_dot (const TrilinosScalar a,
1597  const VectorBase &V,
1598  const VectorBase &W)
1599  {
1600  this->add(a, V);
1601  return *this * W;
1602  }
1603 
1604 
1605 
1606  // inline also scalar products, vector
1607  // additions etc. since they are all
1608  // representable by a single Trilinos
1609  // call. This reduces the overhead of the
1610  // wrapper class.
1611  inline
1612  VectorBase &
1613  VectorBase::operator *= (const TrilinosScalar a)
1614  {
1615  AssertIsFinite(a);
1616 
1617  const int ierr = vector->Scale(a);
1618  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1619 
1620  return *this;
1621  }
1622 
1623 
1624 
1625  inline
1626  VectorBase &
1627  VectorBase::operator /= (const TrilinosScalar a)
1628  {
1629  AssertIsFinite(a);
1630 
1631  const TrilinosScalar factor = 1./a;
1632 
1633  AssertIsFinite(factor);
1634 
1635  const int ierr = vector->Scale(factor);
1636  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1637 
1638  return *this;
1639  }
1640 
1641 
1642 
1643  inline
1644  VectorBase &
1645  VectorBase::operator += (const VectorBase &v)
1646  {
1647  Assert (size() == v.size(),
1648  ExcDimensionMismatch(size(), v.size()));
1649  Assert (vector->Map().SameAs(v.vector->Map()),
1650  ExcDifferentParallelPartitioning());
1651 
1652  const int ierr = vector->Update (1.0, *(v.vector), 1.0);
1653  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1654 
1655  return *this;
1656  }
1657 
1658 
1659 
1660  inline
1661  VectorBase &
1662  VectorBase::operator -= (const VectorBase &v)
1663  {
1664  Assert (size() == v.size(),
1665  ExcDimensionMismatch(size(), v.size()));
1666  Assert (vector->Map().SameAs(v.vector->Map()),
1667  ExcDifferentParallelPartitioning());
1668 
1669  const int ierr = vector->Update (-1.0, *(v.vector), 1.0);
1670  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1671 
1672  return *this;
1673  }
1674 
1675 
1676 
1677  inline
1678  void
1679  VectorBase::add (const TrilinosScalar s)
1680  {
1681  // if we have ghost values, do not allow
1682  // writing to this vector at all.
1683  Assert (!has_ghost_elements(), ExcGhostsPresent());
1684  AssertIsFinite(s);
1685 
1686  size_type n_local = local_size();
1687  for (size_type i=0; i<n_local; i++)
1688  (*vector)[0][i] += s;
1689  }
1690 
1691 
1692 
1693  inline
1694  void
1695  VectorBase::add (const TrilinosScalar a,
1696  const VectorBase &v)
1697  {
1698  // if we have ghost values, do not allow
1699  // writing to this vector at all.
1700  Assert (!has_ghost_elements(), ExcGhostsPresent());
1701  Assert (local_size() == v.local_size(),
1702  ExcDimensionMismatch(local_size(), v.local_size()));
1703 
1704  AssertIsFinite(a);
1705 
1706  const int ierr = vector->Update(a, *(v.vector), 1.);
1707  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1708  }
1709 
1710 
1711 
1712  inline
1713  void
1714  VectorBase::add (const TrilinosScalar a,
1715  const VectorBase &v,
1716  const TrilinosScalar b,
1717  const VectorBase &w)
1718  {
1719  // if we have ghost values, do not allow
1720  // writing to this vector at all.
1721  Assert (!has_ghost_elements(), ExcGhostsPresent());
1722  Assert (local_size() == v.local_size(),
1723  ExcDimensionMismatch(local_size(), v.local_size()));
1724  Assert (local_size() == w.local_size(),
1725  ExcDimensionMismatch(local_size(), w.local_size()));
1726 
1727  AssertIsFinite(a);
1728  AssertIsFinite(b);
1729 
1730  const int ierr = vector->Update(a, *(v.vector), b, *(w.vector), 1.);
1731 
1732  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1733  }
1734 
1735 
1736 
1737  inline
1738  void
1739  VectorBase::sadd (const TrilinosScalar s,
1740  const VectorBase &v)
1741  {
1742  // if we have ghost values, do not allow
1743  // writing to this vector at all.
1744  Assert (!has_ghost_elements(), ExcGhostsPresent());
1745  Assert (size() == v.size(),
1746  ExcDimensionMismatch (size(), v.size()));
1747 
1748  AssertIsFinite(s);
1749 
1750  // We assume that the vectors have the same Map
1751  // if the local size is the same and if the vectors are not ghosted
1752  if (local_size() == v.local_size() && !v.has_ghost_elements())
1753  {
1754  Assert (this->vector->Map().SameAs(v.vector->Map())==true,
1755  ExcDifferentParallelPartitioning());
1756  const int ierr = vector->Update(1., *(v.vector), s);
1757  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1758  }
1759  else
1760  {
1761  (*this) *= s;
1762  this->add(v, true);
1763  }
1764  }
1765 
1766 
1767 
1768  inline
1769  void
1770  VectorBase::sadd (const TrilinosScalar s,
1771  const TrilinosScalar a,
1772  const VectorBase &v)
1773  {
1774  // if we have ghost values, do not allow
1775  // writing to this vector at all.
1776  Assert (!has_ghost_elements(), ExcGhostsPresent());
1777  Assert (size() == v.size(),
1778  ExcDimensionMismatch (size(), v.size()));
1779  AssertIsFinite(s);
1780  AssertIsFinite(a);
1781 
1782  // We assume that the vectors have the same Map
1783  // if the local size is the same and if the vectors are not ghosted
1784  if (local_size() == v.local_size() && !v.has_ghost_elements())
1785  {
1786  Assert (this->vector->Map().SameAs(v.vector->Map())==true,
1787  ExcDifferentParallelPartitioning());
1788  const int ierr = vector->Update(a, *(v.vector), s);
1789  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1790  }
1791  else
1792  {
1793  (*this) *= s;
1794  VectorBase tmp = v;
1795  tmp *= a;
1796  this->add(tmp, true);
1797  }
1798  }
1799 
1800 
1801 
1802  inline
1803  void
1804  VectorBase::sadd (const TrilinosScalar s,
1805  const TrilinosScalar a,
1806  const VectorBase &v,
1807  const TrilinosScalar b,
1808  const VectorBase &w)
1809  {
1810  // if we have ghost values, do not allow
1811  // writing to this vector at all.
1812  Assert (!has_ghost_elements(), ExcGhostsPresent());
1813  Assert (size() == v.size(),
1814  ExcDimensionMismatch (size(), v.size()));
1815  Assert (size() == w.size(),
1816  ExcDimensionMismatch (size(), w.size()));
1817  AssertIsFinite(s);
1818  AssertIsFinite(a);
1819  AssertIsFinite(b);
1820 
1821  // We assume that the vectors have the same Map
1822  // if the local size is the same and if the vectors are not ghosted
1823  if (local_size() == v.local_size() && !v.has_ghost_elements() &&
1824  local_size() == w.local_size() && !w.has_ghost_elements())
1825  {
1826  Assert (this->vector->Map().SameAs(v.vector->Map())==true,
1827  ExcDifferentParallelPartitioning());
1828  Assert (this->vector->Map().SameAs(w.vector->Map())==true,
1829  ExcDifferentParallelPartitioning());
1830  const int ierr = vector->Update(a, *(v.vector), b, *(w.vector), s);
1831  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1832  }
1833  else
1834  {
1835  this->sadd( s, a, v);
1836  this->sadd(1., b, w);
1837  }
1838  }
1839 
1840 
1841 
1842  inline
1843  void
1844  VectorBase::sadd (const TrilinosScalar s,
1845  const TrilinosScalar a,
1846  const VectorBase &v,
1847  const TrilinosScalar b,
1848  const VectorBase &w,
1849  const TrilinosScalar c,
1850  const VectorBase &x)
1851  {
1852  // if we have ghost values, do not allow
1853  // writing to this vector at all.
1854  Assert (!has_ghost_elements(), ExcGhostsPresent());
1855  Assert (size() == v.size(),
1856  ExcDimensionMismatch (size(), v.size()));
1857  Assert (size() == w.size(),
1858  ExcDimensionMismatch (size(), w.size()));
1859  Assert (size() == x.size(),
1860  ExcDimensionMismatch (size(), x.size()));
1861  AssertIsFinite(s);
1862  AssertIsFinite(a);
1863  AssertIsFinite(b);
1864  AssertIsFinite(c);
1865 
1866  // We assume that the vectors have the same Map
1867  // if the local size is the same and if the vectors are not ghosted
1868  if (local_size() == v.local_size() && !v.has_ghost_elements() &&
1869  local_size() == w.local_size() && !w.has_ghost_elements() &&
1870  local_size() == x.local_size() && !x.has_ghost_elements())
1871  {
1872  Assert (this->vector->Map().SameAs(v.vector->Map())==true,
1873  ExcDifferentParallelPartitioning());
1874  Assert (this->vector->Map().SameAs(w.vector->Map())==true,
1875  ExcDifferentParallelPartitioning());
1876  Assert (this->vector->Map().SameAs(x.vector->Map())==true,
1877  ExcDifferentParallelPartitioning());
1878 
1879  // Update member can only
1880  // input two other vectors so
1881  // do it in two steps
1882  const int ierr = vector->Update(a, *(v.vector), b, *(w.vector), s);
1883  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1884 
1885  const int jerr = vector->Update(c, *(x.vector), 1.);
1886  Assert (jerr == 0, ExcTrilinosError(jerr));
1887  (void)jerr; // removes -Wunused-parameter warning in optimized mode
1888  }
1889  else
1890  {
1891  this->sadd( s, a, v);
1892  this->sadd(1., b, w);
1893  this->sadd(1., c, x);
1894  }
1895  }
1896 
1897 
1898 
1899  inline
1900  void
1901  VectorBase::scale (const VectorBase &factors)
1902  {
1903  // if we have ghost values, do not allow
1904  // writing to this vector at all.
1905  Assert (!has_ghost_elements(), ExcGhostsPresent());
1906  Assert (local_size() == factors.local_size(),
1907  ExcDimensionMismatch(local_size(), factors.local_size()));
1908 
1909  const int ierr = vector->Multiply (1.0, *(factors.vector), *vector, 0.0);
1910  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1911  }
1912 
1913 
1914 
1915  inline
1916  void
1917  VectorBase::equ (const TrilinosScalar a,
1918  const VectorBase &v)
1919  {
1920  // if we have ghost values, do not allow
1921  // writing to this vector at all.
1922  Assert (!has_ghost_elements(), ExcGhostsPresent());
1923  AssertIsFinite(a);
1924 
1925  // If we don't have the same map, copy.
1926  if (vector->Map().SameAs(v.vector->Map())==false)
1927  {
1928  this->sadd(0., a, v);
1929  }
1930  else
1931  {
1932  // Otherwise, just update
1933  int ierr = vector->Update(a, *v.vector, 0.0);
1934  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1935 
1936  last_action = Zero;
1937  }
1938 
1939  }
1940 
1941 
1942 
1943  inline
1944  void
1945  VectorBase::ratio (const VectorBase &v,
1946  const VectorBase &w)
1947  {
1948  Assert (v.local_size() == w.local_size(),
1949  ExcDimensionMismatch (v.local_size(), w.local_size()));
1950  Assert (local_size() == w.local_size(),
1951  ExcDimensionMismatch (local_size(), w.local_size()));
1952 
1953  const int ierr = vector->ReciprocalMultiply(1.0, *(w.vector),
1954  *(v.vector), 0.0);
1955 
1956  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1957  }
1958 
1959 
1960 
1961  inline
1962  const Epetra_MultiVector &
1963  VectorBase::trilinos_vector () const
1964  {
1965  return static_cast<const Epetra_MultiVector &>(*vector);
1966  }
1967 
1968 
1969 
1970  inline
1971  Epetra_FEVector &
1972  VectorBase::trilinos_vector ()
1973  {
1974  return *vector;
1975  }
1976 
1977 
1978 
1979  inline
1980  const Epetra_Map &
1981  VectorBase::vector_partitioner () const
1982  {
1983  return static_cast<const Epetra_Map &>(vector->Map());
1984  }
1985 
1986 
1987 
1988  inline
1989  const MPI_Comm &
1990  VectorBase::get_mpi_communicator () const
1991  {
1992  static MPI_Comm comm;
1993 
1994 #ifdef DEAL_II_WITH_MPI
1995 
1996  const Epetra_MpiComm *mpi_comm
1997  = dynamic_cast<const Epetra_MpiComm *>(&vector->Map().Comm());
1998  comm = mpi_comm->Comm();
1999 
2000 #else
2001 
2002  comm = MPI_COMM_SELF;
2003 
2004 #endif
2005 
2006  return comm;
2007  }
2008 
2009 
2010 
2011 #endif // DOXYGEN
2012 
2013 }
2014 
2017 DEAL_II_NAMESPACE_CLOSE
2018 
2019 #endif // DEAL_II_WITH_TRILINOS
2020 
2021 /*---------------------------- trilinos_vector_base.h ---------------------------*/
2022 
2023 #endif
2024 /*---------------------------- trilinos_vector_base.h ---------------------------*/
static ::ExceptionBase & ExcTrilinosError(int arg1)
void sadd(const TrilinosScalar s, const VectorBase &V)
real_type l1_norm() const
std_cxx11::shared_ptr< Epetra_MultiVector > nonlocal_vector
std::pair< size_type, size_type > local_range() const
void scale(const VectorBase &scaling_factors)
STL namespace.
void reinit(const VectorBase &v, const bool omit_zeroing_entries=false)
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
void set(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
static ::ExceptionBase & ExcMessage(std::string arg1)
real_type linfty_norm() const
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
VectorizedArray< Number > min(const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y)
std_cxx11::shared_ptr< Epetra_FEVector > vector
TrilinosScalar minimal_value() const 1
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const Epetra_Map & vector_partitioner() const
#define DeclException0(Exception0)
Definition: exceptions.h:541
std::size_t memory_consumption() const
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< TrilinosScalar > &values) const
static ::ExceptionBase & ExcTrilinosError(int arg1)
void compress(::VectorOperation::values operation)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcGhostsPresent()
TrilinosScalar el(const size_type index) const 1
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
TrilinosScalar mean_value() const
double norm(const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Du)
Definition: divergence.h:532
real_type norm_sqr() const
void ratio(const VectorBase &a, const VectorBase &b) 1
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:600
TrilinosScalar max() const
const MPI_Comm & get_mpi_communicator() const
const Epetra_MultiVector & trilinos_vector() const
bool in_local_range(const size_type index) const
friend class internal::VectorReference
TrilinosScalar add_and_dot(const TrilinosScalar a, const VectorBase &V, const VectorBase &W)
real_type lp_norm(const TrilinosScalar p) const
IndexSet locally_owned_elements() const
static ::ExceptionBase & ExcDifferentParallelPartitioning()
#define AssertIsFinite(number)
Definition: exceptions.h:1197
size_type local_size() const
TrilinosScalar min() const
void print(const char *format=0) const 1
real_type l2_norm() const
void equ(const TrilinosScalar a, const VectorBase &V)