deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
trilinos_vector.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2008 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_trilinos_vector_h
16#define dealii_trilinos_vector_h
17
18
19#include <deal.II/base/config.h>
20
21#ifdef DEAL_II_WITH_TRILINOS
25
28# include <deal.II/lac/vector.h>
31
33# include <Epetra_ConfigDefs.h>
34# include <Epetra_FEVector.h>
35# include <Epetra_LocalMap.h>
36# include <Epetra_Map.h>
37# include <Epetra_MpiComm.h>
39
40# include <memory>
41# include <utility>
42# include <vector>
43
45
46// Forward declarations
47# ifndef DOXYGEN
48namespace LinearAlgebra
49{
50 // Forward declaration
51 template <typename Number>
52 class ReadWriteVector;
53} // namespace LinearAlgebra
54# endif
55
60
66namespace TrilinosWrappers
67{
68 class SparseMatrix;
69
75 {
76 public:
78 };
79
83
90 namespace internal
91 {
101 class VectorReference
102 {
103 private:
104 using size_type = VectorTraits::size_type;
105
110 VectorReference(MPI::Vector &vector, const size_type index);
111
112 public:
116 VectorReference(const VectorReference &) = default;
117
129 const VectorReference &
130 operator=(const VectorReference &r) const;
131
135 VectorReference &
136 operator=(const VectorReference &r);
137
141 const VectorReference &
142 operator=(const TrilinosScalar &s) const;
143
147 const VectorReference &
148 operator+=(const TrilinosScalar &s) const;
149
153 const VectorReference &
154 operator-=(const TrilinosScalar &s) const;
155
159 const VectorReference &
160 operator*=(const TrilinosScalar &s) const;
161
165 const VectorReference &
166 operator/=(const TrilinosScalar &s) const;
167
172 operator TrilinosScalar() const;
173
178 int,
179 << "An error with error number " << arg1
180 << " occurred while calling a Trilinos function");
181
182 private:
186 MPI::Vector &vector;
187
191 const size_type index;
192
193 // Make the vector class a friend, so that it can create objects of the
194 // present type.
195 friend class ::TrilinosWrappers::MPI::Vector;
196 };
197 } // namespace internal
201
202# ifndef DEAL_II_WITH_64BIT_INDICES
203 // define a helper function that queries the global ID of local ID of
204 // an Epetra_BlockMap object by calling either the 32- or 64-bit
205 // function necessary.
206 inline int
207 gid(const Epetra_BlockMap &map, int i)
208 {
209 return map.GID(i);
210 }
211# else
212 // define a helper function that queries the global ID of local ID of
213 // an Epetra_BlockMap object by calling either the 32- or 64-bit
214 // function necessary.
215 inline long long int
216 gid(const Epetra_BlockMap &map, int i)
217 {
218 return map.GID64(i);
219 }
220# endif
221
227 namespace MPI
228 {
229 class BlockVector;
230
405 class Vector : public ReadVector<TrilinosScalar>
406 {
407 public:
417 using const_iterator = const value_type *;
420
430 Vector();
431
435 Vector(const Vector &v);
436
455 explicit Vector(const IndexSet &parallel_partitioning,
456 const MPI_Comm communicator = MPI_COMM_WORLD);
457
469 Vector(const IndexSet &local,
470 const IndexSet &ghost,
471 const MPI_Comm communicator = MPI_COMM_WORLD);
472
487 Vector(const IndexSet &parallel_partitioning,
488 const Vector &v,
489 const MPI_Comm communicator = MPI_COMM_WORLD);
490
503 template <typename Number>
504 Vector(const IndexSet &parallel_partitioning,
505 const ::Vector<Number> &v,
506 const MPI_Comm communicator = MPI_COMM_WORLD);
507
516 Vector(Vector &&v); // NOLINT
517
521 ~Vector() override = default;
522
527 void
528 clear();
529
553 void
554 reinit(const Vector &v,
555 const bool omit_zeroing_entries = false,
556 const bool allow_different_maps = false);
557
579 void
580 reinit(const IndexSet &parallel_partitioning,
581 const MPI_Comm communicator = MPI_COMM_WORLD,
582 const bool omit_zeroing_entries = false);
583
618 void
619 reinit(const IndexSet &locally_owned_entries,
620 const IndexSet &locally_relevant_or_ghost_entries,
621 const MPI_Comm communicator = MPI_COMM_WORLD,
622 const bool vector_writable = false);
623
634 void
635 reinit(
636 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
637 const bool make_ghosted = true,
638 const bool vector_writable = false);
639
643 void
644 reinit(const BlockVector &v, const bool import_data = false);
645
662 void
664
677 Vector &
679
709 Vector &
710 operator=(const Vector &v);
711
716 Vector &
717 operator=(Vector &&v) noexcept;
718
726 template <typename Number>
727 Vector &
728 operator=(const ::Vector<Number> &v);
729
747 void
749 const ::TrilinosWrappers::SparseMatrix &matrix,
750 const Vector &vector);
751
758 void
760 const VectorOperation::values operation);
761
766 void
768 const VectorOperation::values operation)
769 {
770 import_elements(rwv, operation);
771 }
772
773
779 bool
780 operator==(const Vector &v) const;
781
787 bool
788 operator!=(const Vector &v) const;
789
794 size() const override;
795
802
825 std::pair<size_type, size_type>
826 local_range() const;
827
835 bool
836 in_local_range(const size_type index) const;
837
853
861 bool
863
870 void
872
878 operator*(const Vector &vec) const;
879
884 norm_sqr() const;
885
890 mean_value() const;
891
896 min() const;
897
902 max() const;
903
908 l1_norm() const;
909
915 l2_norm() const;
916
922 lp_norm(const TrilinosScalar p) const;
923
928 linfty_norm() const;
929
950 add_and_dot(const TrilinosScalar a, const Vector &V, const Vector &W);
951
957 bool
958 all_zero() const;
959
965 bool
966 is_non_negative() const;
968
969
974
983 operator()(const size_type index);
984
993 operator()(const size_type index) const;
994
1000 reference
1001 operator[](const size_type index);
1002
1009 operator[](const size_type index) const;
1010
1026 void
1027 extract_subvector_to(const std::vector<size_type> &indices,
1028 std::vector<TrilinosScalar> &values) const;
1029
1033 virtual void
1035 const ArrayView<const size_type> &indices,
1036 const ArrayView<TrilinosScalar> &elements) const override;
1037
1065 template <typename ForwardIterator, typename OutputIterator>
1066 void
1067 extract_subvector_to(ForwardIterator indices_begin,
1068 const ForwardIterator indices_end,
1069 OutputIterator values_begin) const;
1070
1079 iterator
1081
1087 begin() const;
1088
1093 iterator
1095
1101 end() const;
1102
1104
1105
1110
1117 void
1118 set(const std::vector<size_type> &indices,
1119 const std::vector<TrilinosScalar> &values);
1120
1125 void
1126 set(const std::vector<size_type> &indices,
1127 const ::Vector<TrilinosScalar> &values);
1128
1134 void
1135 set(const size_type n_elements,
1136 const size_type *indices,
1137 const TrilinosScalar *values);
1138
1143 void
1144 add(const std::vector<size_type> &indices,
1145 const std::vector<TrilinosScalar> &values);
1146
1151 void
1152 add(const std::vector<size_type> &indices,
1153 const ::Vector<TrilinosScalar> &values);
1154
1160 void
1161 add(const size_type n_elements,
1162 const size_type *indices,
1163 const TrilinosScalar *values);
1164
1168 Vector &
1170
1174 Vector &
1176
1180 Vector &
1182
1186 Vector &
1188
1193 void
1195
1208 void
1209 add(const Vector &V, const bool allow_different_maps = false);
1210
1214 void
1215 add(const TrilinosScalar a, const Vector &V);
1216
1220 void
1222 const Vector &V,
1223 const TrilinosScalar b,
1224 const Vector &W);
1225
1230 void
1231 sadd(const TrilinosScalar s, const Vector &V);
1232
1236 void
1237 sadd(const TrilinosScalar s, const TrilinosScalar a, const Vector &V);
1238
1244 void
1245 scale(const Vector &scaling_factors);
1246
1250 void
1251 equ(const TrilinosScalar a, const Vector &V);
1253
1258
1263 const Epetra_MultiVector &
1265
1270 Epetra_FEVector &
1272
1277 const Epetra_BlockMap &
1279
1287 void
1288 print(std::ostream &out,
1289 const unsigned int precision = 3,
1290 const bool scientific = true,
1291 const bool across = true) const;
1292
1306 void
1307 swap(Vector &v) noexcept;
1308
1312 std::size_t
1313 memory_consumption() const;
1314
1318 MPI_Comm
1321
1326
1331 int,
1332 << "An error with error number " << arg1
1333 << " occurred while calling a Trilinos function");
1334
1340 size_type,
1341 size_type,
1342 size_type,
1343 size_type,
1344 << "You are trying to access element " << arg1
1345 << " of a distributed vector, but this element is not stored "
1346 << "on the current processor. Note: There are " << arg2
1347 << " elements stored "
1348 << "on the current processor from within the range [" << arg3 << ','
1349 << arg4 << "] but Trilinos vectors need not store contiguous "
1350 << "ranges on each processor, and not every element in "
1351 << "this range may in fact be stored locally."
1352 << "\n\n"
1353 << "A common source for this kind of problem is that you "
1354 << "are passing a 'fully distributed' vector into a function "
1355 << "that needs read access to vector elements that correspond "
1356 << "to degrees of freedom on ghost cells (or at least to "
1357 << "'locally active' degrees of freedom that are not also "
1358 << "'locally owned'). You need to pass a vector that has these "
1359 << "elements as ghost entries.");
1360
1361 private:
1373 Epetra_CombineMode last_action;
1374
1380
1386
1392 std::unique_ptr<Epetra_FEVector> vector;
1393
1399 std::unique_ptr<Epetra_MultiVector> nonlocal_vector;
1400
1405
1406 // Make the reference class a friend.
1408 };
1409
1410
1411
1412 // ------------------- inline and template functions --------------
1413
1414
1422 inline void
1423 swap(Vector &u, Vector &v) noexcept
1424 {
1425 u.swap(v);
1426 }
1427 } // namespace MPI
1428
1429# ifndef DOXYGEN
1430
1431 namespace internal
1432 {
1433 inline VectorReference::VectorReference(MPI::Vector &vector,
1434 const size_type index)
1435 : vector(vector)
1436 , index(index)
1437 {}
1438
1439
1440 inline const VectorReference &
1441 VectorReference::operator=(const VectorReference &r) const
1442 {
1443 // as explained in the class
1444 // documentation, this is not the copy
1445 // operator. so simply pass on to the
1446 // "correct" assignment operator
1447 *this = static_cast<TrilinosScalar>(r);
1448
1449 return *this;
1450 }
1451
1452
1453
1454 inline VectorReference &
1455 VectorReference::operator=(const VectorReference &r)
1456 {
1457 // as above
1458 *this = static_cast<TrilinosScalar>(r);
1459
1460 return *this;
1461 }
1462
1463
1464 inline const VectorReference &
1465 VectorReference::operator=(const TrilinosScalar &value) const
1466 {
1467 vector.set(1, &index, &value);
1468 return *this;
1469 }
1470
1471
1472
1473 inline const VectorReference &
1474 VectorReference::operator+=(const TrilinosScalar &value) const
1475 {
1476 vector.add(1, &index, &value);
1477 return *this;
1478 }
1479
1480
1481
1482 inline const VectorReference &
1483 VectorReference::operator-=(const TrilinosScalar &value) const
1484 {
1485 TrilinosScalar new_value = -value;
1486 vector.add(1, &index, &new_value);
1487 return *this;
1488 }
1489
1490
1491
1492 inline const VectorReference &
1493 VectorReference::operator*=(const TrilinosScalar &value) const
1494 {
1495 TrilinosScalar new_value = static_cast<TrilinosScalar>(*this) * value;
1496 vector.set(1, &index, &new_value);
1497 return *this;
1498 }
1499
1500
1501
1502 inline const VectorReference &
1503 VectorReference::operator/=(const TrilinosScalar &value) const
1504 {
1505 TrilinosScalar new_value = static_cast<TrilinosScalar>(*this) / value;
1506 vector.set(1, &index, &new_value);
1507 return *this;
1508 }
1509 } // namespace internal
1510
1511 namespace MPI
1512 {
1513 inline bool
1514 Vector::in_local_range(const size_type index) const
1515 {
1516 std::pair<size_type, size_type> range = local_range();
1517
1518 return ((index >= range.first) && (index < range.second));
1519 }
1520
1521
1522
1523 inline IndexSet
1525 {
1526 Assert(owned_elements.size() == size(),
1527 ExcMessage(
1528 "The locally owned elements have not been properly initialized!"
1529 " This happens for example if this object has been initialized"
1530 " with exactly one overlapping IndexSet."));
1531 return owned_elements;
1532 }
1533
1534
1535
1536 inline bool
1538 {
1539 return has_ghosts;
1540 }
1541
1542
1543
1544 inline void
1546 {}
1547
1548
1549
1550 inline internal::VectorReference
1551 Vector::operator()(const size_type index)
1552 {
1553 return internal::VectorReference(*this, index);
1554 }
1555
1556
1557
1558 inline internal::VectorReference
1559 Vector::operator[](const size_type index)
1560 {
1561 return operator()(index);
1562 }
1563
1564
1565
1566 inline TrilinosScalar
1567 Vector::operator[](const size_type index) const
1568 {
1569 return operator()(index);
1570 }
1571
1572
1573
1574 inline void
1575 Vector::extract_subvector_to(const std::vector<size_type> &indices,
1576 std::vector<TrilinosScalar> &values) const
1577 {
1578 for (size_type i = 0; i < indices.size(); ++i)
1579 values[i] = operator()(indices[i]);
1580 }
1581
1582
1583
1584 inline void
1586 const ArrayView<const size_type> &indices,
1587 const ArrayView<TrilinosScalar> &elements) const
1588 {
1589 AssertDimension(indices.size(), elements.size());
1590 for (unsigned int i = 0; i < indices.size(); ++i)
1591 {
1592 AssertIndexRange(indices[i], size());
1593 elements[i] = (*this)[indices[i]];
1594 }
1595 }
1596
1597
1598
1599 template <typename ForwardIterator, typename OutputIterator>
1600 inline void
1601 Vector::extract_subvector_to(ForwardIterator indices_begin,
1602 const ForwardIterator indices_end,
1603 OutputIterator values_begin) const
1604 {
1605 while (indices_begin != indices_end)
1606 {
1607 *values_begin = operator()(*indices_begin);
1608 ++indices_begin;
1609 ++values_begin;
1610 }
1611 }
1612
1613
1614
1615 inline Vector::iterator
1617 {
1618 return (*vector)[0];
1619 }
1620
1621
1622
1623 inline Vector::iterator
1624 Vector::end()
1625 {
1626 return (*vector)[0] + locally_owned_size();
1627 }
1628
1629
1630
1632 Vector::begin() const
1633 {
1634 return (*vector)[0];
1635 }
1636
1637
1638
1640 Vector::end() const
1641 {
1642 return (*vector)[0] + locally_owned_size();
1643 }
1644
1645
1646
1647 inline void
1648 Vector::set(const std::vector<size_type> &indices,
1649 const std::vector<TrilinosScalar> &values)
1650 {
1651 // if we have ghost values, do not allow
1652 // writing to this vector at all.
1654
1655 AssertDimension(indices.size(), values.size());
1656
1657 set(indices.size(), indices.data(), values.data());
1658 }
1659
1660
1661
1662 inline void
1663 Vector::set(const std::vector<size_type> &indices,
1664 const ::Vector<TrilinosScalar> &values)
1665 {
1666 // if we have ghost values, do not allow
1667 // writing to this vector at all.
1669
1670 AssertDimension(indices.size(), values.size());
1671
1672 set(indices.size(), indices.data(), values.begin());
1673 }
1674
1675
1676
1677 inline void
1678 Vector::set(const size_type n_elements,
1679 const size_type *indices,
1680 const TrilinosScalar *values)
1681 {
1682 // if we have ghost values, do not allow
1683 // writing to this vector at all.
1685
1686 if (last_action == Add)
1687 {
1688 const int ierr = vector->GlobalAssemble(Add);
1689 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1690 }
1691
1692 if (last_action != Insert)
1693 last_action = Insert;
1694
1695 for (size_type i = 0; i < n_elements; ++i)
1696 {
1697 const TrilinosWrappers::types::int_type row = indices[i];
1698 const TrilinosWrappers::types::int_type local_row =
1699 vector->Map().LID(row);
1700 if (local_row != -1)
1701 (*vector)[0][local_row] = values[i];
1702 else
1703 {
1704 const int ierr = vector->ReplaceGlobalValues(1, &row, &values[i]);
1705 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1706 compressed = false;
1707 }
1708 // in set operation, do not use the pre-allocated vector for nonlocal
1709 // entries even if it exists. This is to ensure that we really only
1710 // set the elements touched by the set() method and not all contained
1711 // in the nonlocal entries vector (there is no way to distinguish them
1712 // on the receiving processor)
1713 }
1714 }
1715
1716
1717
1718 inline void
1719 Vector::add(const std::vector<size_type> &indices,
1720 const std::vector<TrilinosScalar> &values)
1721 {
1722 // if we have ghost values, do not allow
1723 // writing to this vector at all.
1725 AssertDimension(indices.size(), values.size());
1726
1727 add(indices.size(), indices.data(), values.data());
1728 }
1729
1730
1731
1732 inline void
1733 Vector::add(const std::vector<size_type> &indices,
1734 const ::Vector<TrilinosScalar> &values)
1735 {
1736 // if we have ghost values, do not allow
1737 // writing to this vector at all.
1739 AssertDimension(indices.size(), values.size());
1740
1741 add(indices.size(), indices.data(), values.begin());
1742 }
1743
1744
1745
1746 inline void
1747 Vector::add(const size_type n_elements,
1748 const size_type *indices,
1749 const TrilinosScalar *values)
1750 {
1751 // if we have ghost values, do not allow
1752 // writing to this vector at all.
1754
1755 if (last_action != Add)
1756 {
1757 if (last_action == Insert)
1758 {
1759 const int ierr = vector->GlobalAssemble(Insert);
1760 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1761 }
1762 last_action = Add;
1763 }
1764
1765 for (size_type i = 0; i < n_elements; ++i)
1766 {
1767 const size_type row = indices[i];
1768 const TrilinosWrappers::types::int_type local_row = vector->Map().LID(
1769 static_cast<TrilinosWrappers::types::int_type>(row));
1770 if (local_row != -1)
1771 (*vector)[0][local_row] += values[i];
1772 else if (nonlocal_vector.get() == nullptr)
1773 {
1774 const int ierr = vector->SumIntoGlobalValues(
1775 1,
1776 reinterpret_cast<const TrilinosWrappers::types::int_type *>(
1777 &row),
1778 &values[i]);
1779 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1780 compressed = false;
1781 }
1782 else
1783 {
1784 // use pre-allocated vector for non-local entries if it exists for
1785 // addition operation
1787 nonlocal_vector->Map().LID(
1788 static_cast<TrilinosWrappers::types::int_type>(row));
1789 Assert(my_row != -1,
1790 ExcMessage(
1791 "Attempted to write into off-processor vector entry "
1792 "that has not be specified as being writable upon "
1793 "initialization"));
1794 (*nonlocal_vector)[0][my_row] += values[i];
1795 compressed = false;
1796 }
1797 }
1798 }
1799
1800
1801
1802 inline Vector::size_type
1803 Vector::size() const
1804 {
1805# ifndef DEAL_II_WITH_64BIT_INDICES
1806 return vector->Map().MaxAllGID() + 1 - vector->Map().MinAllGID();
1807# else
1808 return vector->Map().MaxAllGID64() + 1 - vector->Map().MinAllGID64();
1809# endif
1810 }
1811
1812
1813
1814 inline Vector::size_type
1816 {
1817 return owned_elements.n_elements();
1818 }
1819
1820
1821
1822 inline std::pair<Vector::size_type, Vector::size_type>
1823 Vector::local_range() const
1824 {
1825# ifndef DEAL_II_WITH_64BIT_INDICES
1826 const TrilinosWrappers::types::int_type begin = vector->Map().MinMyGID();
1828 vector->Map().MaxMyGID() + 1;
1829# else
1831 vector->Map().MinMyGID64();
1833 vector->Map().MaxMyGID64() + 1;
1834# endif
1835
1836 Assert(
1837 end - begin == vector->Map().NumMyElements(),
1838 ExcMessage(
1839 "This function only makes sense if the elements that this "
1840 "vector stores on the current processor form a contiguous range. "
1841 "This does not appear to be the case for the current vector."));
1842
1843 return std::make_pair(begin, end);
1844 }
1845
1846
1847
1848 inline TrilinosScalar
1849 Vector::operator*(const Vector &vec) const
1850 {
1851 Assert(vector->Map().SameAs(vec.vector->Map()),
1854
1855 TrilinosScalar result;
1856
1857 const int ierr = vector->Dot(*(vec.vector), &result);
1858 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1859
1860 return result;
1861 }
1862
1863
1864
1865 inline Vector::real_type
1866 Vector::norm_sqr() const
1867 {
1868 const TrilinosScalar d = l2_norm();
1869 return d * d;
1870 }
1871
1872
1873
1874 inline TrilinosScalar
1875 Vector::mean_value() const
1876 {
1878
1880 const int ierr = vector->MeanValue(&mean);
1881 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1882
1883 return mean;
1884 }
1885
1886
1887
1888 inline TrilinosScalar
1889 Vector::min() const
1890 {
1891 TrilinosScalar min_value;
1892 const int ierr = vector->MinValue(&min_value);
1893 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1894
1895 return min_value;
1896 }
1897
1898
1899
1900 inline TrilinosScalar
1901 Vector::max() const
1902 {
1903 TrilinosScalar max_value;
1904 const int ierr = vector->MaxValue(&max_value);
1905 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1906
1907 return max_value;
1908 }
1909
1910
1911
1912 inline Vector::real_type
1913 Vector::l1_norm() const
1914 {
1916
1918 const int ierr = vector->Norm1(&d);
1919 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1920
1921 return d;
1922 }
1923
1924
1925
1926 inline Vector::real_type
1927 Vector::l2_norm() const
1928 {
1930
1932 const int ierr = vector->Norm2(&d);
1933 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1934
1935 return d;
1936 }
1937
1938
1939
1940 inline Vector::real_type
1941 Vector::lp_norm(const TrilinosScalar p) const
1942 {
1944
1945 TrilinosScalar norm = 0;
1946 TrilinosScalar sum = 0;
1947 const size_type n_local = locally_owned_size();
1948
1949 // loop over all the elements because
1950 // Trilinos does not support lp norms
1951 for (size_type i = 0; i < n_local; ++i)
1952 sum += std::pow(std::fabs((*vector)[0][i]), p);
1953
1954 norm = std::pow(sum, static_cast<TrilinosScalar>(1. / p));
1955
1956 return norm;
1957 }
1958
1959
1960
1961 inline Vector::real_type
1962 Vector::linfty_norm() const
1963 {
1964 // while we disallow the other
1965 // norm operations on ghosted
1966 // vectors, this particular norm
1967 // is safe to run even in the
1968 // presence of ghost elements
1970 const int ierr = vector->NormInf(&d);
1971 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1972
1973 return d;
1974 }
1975
1976
1977
1978 inline TrilinosScalar
1980 const Vector &V,
1981 const Vector &W)
1982 {
1983 this->add(a, V);
1984 return *this * W;
1985 }
1986
1987
1988
1989 // inline also scalar products, vector
1990 // additions etc. since they are all
1991 // representable by a single Trilinos
1992 // call. This reduces the overhead of the
1993 // wrapper class.
1994 inline Vector &
1996 {
1997 AssertIsFinite(a);
1998
1999 const int ierr = vector->Scale(a);
2000 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2001
2002 return *this;
2003 }
2004
2005
2006
2007 inline Vector &
2009 {
2010 AssertIsFinite(a);
2011
2012 const TrilinosScalar factor = 1. / a;
2013
2014 AssertIsFinite(factor);
2015
2016 const int ierr = vector->Scale(factor);
2017 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2018
2019 return *this;
2020 }
2021
2022
2023
2024 inline Vector &
2025 Vector::operator+=(const Vector &v)
2026 {
2027 AssertDimension(size(), v.size());
2028 Assert(vector->Map().SameAs(v.vector->Map()),
2030
2031 const int ierr = vector->Update(1.0, *(v.vector), 1.0);
2032 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2033
2034 return *this;
2035 }
2036
2037
2038
2039 inline Vector &
2040 Vector::operator-=(const Vector &v)
2041 {
2042 AssertDimension(size(), v.size());
2043 Assert(vector->Map().SameAs(v.vector->Map()),
2045
2046 const int ierr = vector->Update(-1.0, *(v.vector), 1.0);
2047 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2048
2049 return *this;
2050 }
2051
2052
2053
2054 inline void
2056 {
2057 // if we have ghost values, do not allow
2058 // writing to this vector at all.
2060 AssertIsFinite(s);
2061
2062 size_type n_local = locally_owned_size();
2063 for (size_type i = 0; i < n_local; ++i)
2064 (*vector)[0][i] += s;
2065 }
2066
2067
2068
2069 inline void
2070 Vector::add(const TrilinosScalar a, const Vector &v)
2071 {
2072 // if we have ghost values, do not allow
2073 // writing to this vector at all.
2076
2077 AssertIsFinite(a);
2078
2079 const int ierr = vector->Update(a, *(v.vector), 1.);
2080 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2081 }
2082
2083
2084
2085 inline void
2087 const Vector &v,
2088 const TrilinosScalar b,
2089 const Vector &w)
2090 {
2091 // if we have ghost values, do not allow
2092 // writing to this vector at all.
2095 AssertDimension(locally_owned_size(), w.locally_owned_size());
2096
2097 AssertIsFinite(a);
2098 AssertIsFinite(b);
2099
2100 const int ierr = vector->Update(a, *(v.vector), b, *(w.vector), 1.);
2101
2102 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2103 }
2104
2105
2106
2107 inline void
2108 Vector::sadd(const TrilinosScalar s, const Vector &v)
2109 {
2110 // if we have ghost values, do not allow
2111 // writing to this vector at all.
2113 AssertDimension(size(), v.size());
2114
2115 AssertIsFinite(s);
2116
2117 // We assume that the vectors have the same Map
2118 // if the local size is the same and if the vectors are not ghosted
2120 !v.has_ghost_elements())
2121 {
2122 Assert(this->vector->Map().SameAs(v.vector->Map()) == true,
2124 const int ierr = vector->Update(1., *(v.vector), s);
2125 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2126 }
2127 else
2128 {
2129 (*this) *= s;
2130 this->add(v, true);
2131 }
2132 }
2133
2134
2135
2136 inline void
2138 const TrilinosScalar a,
2139 const Vector &v)
2140 {
2141 // if we have ghost values, do not allow
2142 // writing to this vector at all.
2144 AssertDimension(size(), v.size());
2145 AssertIsFinite(s);
2146 AssertIsFinite(a);
2147
2148 // We assume that the vectors have the same Map
2149 // if the local size is the same and if the vectors are not ghosted
2151 !v.has_ghost_elements())
2152 {
2153 Assert(this->vector->Map().SameAs(v.vector->Map()) == true,
2155 const int ierr = vector->Update(a, *(v.vector), s);
2156 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2157 }
2158 else
2159 {
2160 (*this) *= s;
2161 Vector tmp = v;
2162 tmp *= a;
2163 this->add(tmp, true);
2164 }
2165 }
2166
2167
2168
2169 inline void
2170 Vector::scale(const Vector &factors)
2171 {
2172 // if we have ghost values, do not allow
2173 // writing to this vector at all.
2176
2177 const int ierr = vector->Multiply(1.0, *(factors.vector), *vector, 0.0);
2178 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2179 }
2180
2181
2182
2183 inline void
2184 Vector::equ(const TrilinosScalar a, const Vector &v)
2185 {
2186 // if we have ghost values, do not allow
2187 // writing to this vector at all.
2189 AssertIsFinite(a);
2190
2191 // If we don't have the same map, copy.
2192 if (vector->Map().SameAs(v.vector->Map()) == false)
2193 {
2194 this->sadd(0., a, v);
2195 }
2196 else
2197 {
2198 // Otherwise, just update
2199 int ierr = vector->Update(a, *v.vector, 0.0);
2200 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2201
2202 last_action = Zero;
2203 }
2204 }
2205
2206
2207
2208 inline const Epetra_MultiVector &
2210 {
2211 return static_cast<const Epetra_MultiVector &>(*vector);
2212 }
2213
2214
2215
2216 inline Epetra_FEVector &
2218 {
2219 return *vector;
2220 }
2221
2222
2223
2224 inline const Epetra_BlockMap &
2226 {
2227 return vector->Map();
2228 }
2229
2230
2231
2232 inline MPI_Comm
2234 {
2235 const Epetra_MpiComm *mpi_comm =
2236 dynamic_cast<const Epetra_MpiComm *>(&vector->Map().Comm());
2237 return mpi_comm->Comm();
2238 }
2239
2240
2241
2242 template <typename number>
2243 Vector::Vector(const IndexSet &parallel_partitioner,
2244 const ::Vector<number> &v,
2245 const MPI_Comm communicator)
2246 {
2247 *this =
2248 Vector(parallel_partitioner.make_trilinos_map(communicator, true), v);
2249 owned_elements = parallel_partitioner;
2250 }
2251
2252
2253
2254 inline Vector &
2256 {
2257 AssertIsFinite(s);
2258
2259 int ierr = vector->PutScalar(s);
2260 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2261
2262 if (nonlocal_vector.get() != nullptr)
2263 {
2264 ierr = nonlocal_vector->PutScalar(0.);
2265 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2266 }
2267
2268 return *this;
2269 }
2270 } /* end of namespace MPI */
2271
2272# endif /* DOXYGEN */
2273
2274} /* end of namespace TrilinosWrappers */
2275
2277
2278
2279namespace internal
2280{
2282 {
2283 template <typename>
2284 class ReinitHelper;
2285
2290 template <>
2291 class ReinitHelper<TrilinosWrappers::MPI::Vector>
2292 {
2293 public:
2294 template <typename Matrix>
2295 static void
2296 reinit_range_vector(const Matrix &matrix,
2298 bool omit_zeroing_entries)
2299 {
2300 v.reinit(matrix.locally_owned_range_indices(),
2301 matrix.get_mpi_communicator(),
2302 omit_zeroing_entries);
2303 }
2304
2305 template <typename Matrix>
2306 static void
2307 reinit_domain_vector(const Matrix &matrix,
2309 bool omit_zeroing_entries)
2310 {
2311 v.reinit(matrix.locally_owned_domain_indices(),
2312 matrix.get_mpi_communicator(),
2313 omit_zeroing_entries);
2314 }
2315 };
2316
2317 } // namespace LinearOperatorImplementation
2318} /* namespace internal */
2319
2320
2321
2325template <>
2326struct is_serial_vector<TrilinosWrappers::MPI::Vector> : std::false_type
2327{};
2328
2329
2331
2332#else
2333
2334// Make sure the scripts that create the C++20 module input files have
2335// something to latch on if the preprocessor #ifdef above would
2336// otherwise lead to an empty content of the file.
2339
2340#endif
2341
2342#endif
std::size_t size() const
Definition array_view.h:689
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Vector & operator/=(const TrilinosScalar factor)
void compress(VectorOperation::values operation)
Vector(const IndexSet &parallel_partitioning, const ::Vector< Number > &v, const MPI_Comm communicator=MPI_COMM_WORLD)
void sadd(const TrilinosScalar s, const Vector &V)
TrilinosScalar mean_value() const
void add(const size_type n_elements, const size_type *indices, const TrilinosScalar *values)
VectorTraits::size_type size_type
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
std::unique_ptr< Epetra_MultiVector > nonlocal_vector
void add(const TrilinosScalar s)
void import_elements(const LinearAlgebra::ReadWriteVector< double > &rwv, const VectorOperation::values operation)
MPI_Comm get_mpi_communicator() const
void swap(Vector &v) noexcept
const Epetra_BlockMap & trilinos_partitioner() const
reference operator()(const size_type index)
void import_nonlocal_data_for_fe(const ::TrilinosWrappers::SparseMatrix &matrix, const Vector &vector)
std::unique_ptr< Epetra_FEVector > vector
size_type size() const override
bool in_local_range(const size_type index) const
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
friend class internal::VectorReference
real_type lp_norm(const TrilinosScalar p) const
IndexSet locally_owned_elements() const
Vector & operator-=(const Vector &V)
Vector & operator+=(const Vector &V)
const Epetra_MultiVector & trilinos_vector() const
const internal::VectorReference const_reference
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
std::pair< size_type, size_type > local_range() const
real_type norm_sqr() const
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< TrilinosScalar > &values) const
bool operator!=(const Vector &v) const
~Vector() override=default
virtual void extract_subvector_to(const ArrayView< const size_type > &indices, const ArrayView< TrilinosScalar > &elements) const override
Epetra_FEVector & trilinos_vector()
const_iterator begin() const
real_type linfty_norm() const
internal::VectorReference reference
TrilinosScalar min() const
void set(const std::vector< size_type > &indices, const ::Vector< TrilinosScalar > &values)
void set(const size_type n_elements, const size_type *indices, const TrilinosScalar *values)
void scale(const Vector &scaling_factors)
void equ(const TrilinosScalar a, const Vector &V)
bool operator==(const Vector &v) const
void sadd(const TrilinosScalar s, const TrilinosScalar a, const Vector &V)
reference operator[](const size_type index)
const value_type * const_iterator
void swap(Vector &u, Vector &v) noexcept
size_type locally_owned_size() const
TrilinosScalar operator[](const size_type index) const
TrilinosScalar add_and_dot(const TrilinosScalar a, const Vector &V, const Vector &W)
Vector & operator=(const TrilinosScalar s)
Vector & operator*=(const TrilinosScalar factor)
TrilinosScalar operator*(const Vector &vec) const
void add(const std::vector< size_type > &indices, const ::Vector< TrilinosScalar > &values)
void extract_subvector_to(ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const
std::size_t memory_consumption() const
void add(const TrilinosScalar a, const Vector &V, const TrilinosScalar b, const Vector &W)
void add(const TrilinosScalar a, const Vector &V)
void set(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
TrilinosScalar max() const
Vector & operator=(const ::Vector< Number > &v)
const_iterator end() const
::types::global_dof_index size_type
bool has_ghost_elements() const
virtual size_type size() const override
size_type locally_owned_size() const
static void reinit_domain_vector(const Matrix &matrix, TrilinosWrappers::MPI::Vector &v, bool omit_zeroing_entries)
static void reinit_range_vector(const Matrix &matrix, TrilinosWrappers::MPI::Vector &v, bool omit_zeroing_entries)
#define DEAL_II_DEPRECATED
Definition config.h:286
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
#define DeclException0(Exception0)
static ::ExceptionBase & ExcGhostsPresent()
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
static ::ExceptionBase & ExcDifferentParallelPartitioning()
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define Assert(cond, exc)
#define AssertIsFinite(number)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
const bool IsBlockVector< VectorType >::value
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
int gid(const Epetra_BlockMap &map, int i)
T sum(const T &t, const MPI_Comm mpi_communicator)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
unsigned int global_dof_index
Definition types.h:94
double TrilinosScalar
Definition types.h:198