16 #ifndef dealii_trilinos_vector_h
17 #define dealii_trilinos_vector_h
20 #include <deal.II/base/config.h>
22 #ifdef DEAL_II_WITH_TRILINOS
23 # include <deal.II/base/index_set.h>
24 # include <deal.II/base/mpi.h>
25 # include <deal.II/base/subscriptor.h>
26 # include <deal.II/base/utilities.h>
28 # include <deal.II/lac/exceptions.h>
29 # include <deal.II/lac/vector.h>
30 # include <deal.II/lac/vector_operation.h>
31 # include <deal.II/lac/vector_type_traits.h>
33 # include <Epetra_ConfigDefs.h>
38 # ifdef DEAL_II_WITH_MPI // only if MPI is installed
39 # include <Epetra_MpiComm.h>
42 # include <Epetra_SerialComm.h>
44 # include <Epetra_FEVector.h>
45 # include <Epetra_LocalMap.h>
46 # include <Epetra_Map.h>
48 DEAL_II_NAMESPACE_OPEN
53 template <
typename Number>
54 class ReadWriteVector;
118 const VectorReference &
119 operator=(
const VectorReference &r)
const;
125 operator=(
const VectorReference &r);
130 const VectorReference &
131 operator=(
const TrilinosScalar &s)
const;
136 const VectorReference &
137 operator+=(
const TrilinosScalar &s)
const;
142 const VectorReference &
143 operator-=(
const TrilinosScalar &s)
const;
148 const VectorReference &
149 operator*=(
const TrilinosScalar &s)
const;
154 const VectorReference &
155 operator/=(
const TrilinosScalar &s)
const;
161 operator TrilinosScalar()
const;
168 <<
"An error with error number " << arg1
169 <<
" occurred while calling a Trilinos function");
193 # ifndef DEAL_II_WITH_64BIT_INDICES
198 gid(
const Epetra_BlockMap &map,
int i)
207 gid(
const Epetra_BlockMap &map,
int i)
408 using real_type = TrilinosScalar;
450 const MPI_Comm &communicator = MPI_COMM_WORLD);
465 const MPI_Comm &communicator = MPI_COMM_WORLD);
483 const MPI_Comm &communicator = MPI_COMM_WORLD);
497 template <
typename Number>
499 const ::Vector<Number> &v,
500 const MPI_Comm & communicator = MPI_COMM_WORLD);
545 const bool omit_zeroing_entries =
false,
546 const bool allow_different_maps =
false);
572 const MPI_Comm &communicator = MPI_COMM_WORLD,
573 const bool omit_zeroing_entries =
false);
603 const MPI_Comm &communicator = MPI_COMM_WORLD,
604 const bool vector_writable =
false);
668 template <
typename Number>
691 const ::TrilinosWrappers::SparseMatrix &matrix,
762 std::pair<size_type, size_type>
858 lp_norm(
const TrilinosScalar p)
const;
943 TrilinosScalar
operator[](
const size_type index)
const;
962 std::vector<TrilinosScalar> & values)
const;
991 template <
typename ForwardIterator,
typename OutputIterator>
994 const ForwardIterator indices_end,
995 OutputIterator values_begin)
const;
1046 set(
const std::vector<size_type> & indices,
1047 const std::vector<TrilinosScalar> &values);
1054 set(
const std::vector<size_type> & indices,
1055 const ::Vector<TrilinosScalar> &values);
1063 set(
const size_type n_elements,
1064 const size_type * indices,
1065 const TrilinosScalar *values);
1072 add(
const std::vector<size_type> & indices,
1073 const std::vector<TrilinosScalar> &values);
1080 add(
const std::vector<size_type> & indices,
1081 const ::Vector<TrilinosScalar> &values);
1089 add(
const size_type n_elements,
1090 const size_type * indices,
1091 const TrilinosScalar *values);
1122 add(
const TrilinosScalar s);
1137 add(
const Vector &V,
const bool allow_different_maps =
false);
1143 add(
const TrilinosScalar a,
const Vector &V);
1149 add(
const TrilinosScalar a,
1151 const TrilinosScalar b,
1159 sadd(
const TrilinosScalar s,
const Vector &V);
1165 sadd(
const TrilinosScalar s,
const TrilinosScalar a,
const Vector &V);
1179 equ(
const TrilinosScalar a,
const Vector &V);
1191 const Epetra_MultiVector &
1215 const Epetra_BlockMap &
1226 print(std::ostream & out,
1227 const unsigned int precision = 3,
1228 const bool scientific =
true,
1229 const bool across =
true)
const;
1271 <<
"An error with error number " << arg1
1272 <<
" occurred while calling a Trilinos function");
1283 <<
"You tried to access element " << arg1
1284 <<
" of a distributed vector, but this element is not stored "
1285 <<
"on the current processor. Note: There are " << arg2
1286 <<
" elements stored "
1287 <<
"on the current processor from within the range " << arg3
1288 <<
" through " << arg4
1289 <<
" but Trilinos vectors need not store contiguous "
1290 <<
"ranges on each processor, and not every element in "
1291 <<
"this range may in fact be stored locally.");
1368 inline VectorReference::VectorReference(
MPI::Vector & vector,
1375 inline const VectorReference &
1376 VectorReference::operator=(
const VectorReference &r)
const
1382 *
this = static_cast<TrilinosScalar>(r);
1389 inline VectorReference &
1390 VectorReference::operator=(
const VectorReference &r)
1393 *
this = static_cast<TrilinosScalar>(r);
1399 inline const VectorReference &
1400 VectorReference::operator=(
const TrilinosScalar &value)
const
1402 vector.set(1, &index, &value);
1408 inline const VectorReference &
1409 VectorReference::operator+=(
const TrilinosScalar &value)
const
1411 vector.add(1, &index, &value);
1417 inline const VectorReference &
1418 VectorReference::operator-=(
const TrilinosScalar &value)
const
1420 TrilinosScalar new_value = -value;
1421 vector.add(1, &index, &new_value);
1427 inline const VectorReference &
1428 VectorReference::operator*=(
const TrilinosScalar &value)
const
1430 TrilinosScalar new_value = static_cast<TrilinosScalar>(*
this) * value;
1431 vector.set(1, &index, &new_value);
1437 inline const VectorReference &
1438 VectorReference::operator/=(
const TrilinosScalar &value)
const
1440 TrilinosScalar new_value = static_cast<TrilinosScalar>(*
this) / value;
1441 vector.set(1, &index, &new_value);
1451 std::pair<size_type, size_type> range = local_range();
1453 return ((index >= range.first) && (index < range.second));
1461 Assert(owned_elements.size() == size(),
1463 "The locally owned elements have not been properly initialized!"
1464 " This happens for example if this object has been initialized"
1465 " with exactly one overlapping IndexSet."));
1466 return owned_elements;
1472 Vector::has_ghost_elements()
const
1485 inline internal::VectorReference
1488 return internal::VectorReference(*
this, index);
1495 return operator()(index);
1502 return operator()(index);
1509 std::vector<TrilinosScalar> & values)
const
1511 for (
size_type i = 0; i < indices.size(); ++i)
1512 values[i] =
operator()(indices[i]);
1517 template <
typename ForwardIterator,
typename OutputIterator>
1520 const ForwardIterator indices_end,
1521 OutputIterator values_begin)
const
1523 while (indices_begin != indices_end)
1525 *values_begin = operator()(*indices_begin);
1533 inline Vector::iterator
1536 return (*vector)[0];
1541 inline Vector::iterator
1544 return (*vector)[0] + local_size();
1549 inline Vector::const_iterator
1552 return (*vector)[0];
1557 inline Vector::const_iterator
1560 return (*vector)[0] + local_size();
1566 Vector::set(
const std::vector<size_type> & indices,
1567 const std::vector<TrilinosScalar> &values)
1573 Assert(indices.size() == values.size(),
1576 set(indices.size(), indices.data(), values.data());
1582 Vector::set(
const std::vector<size_type> & indices,
1583 const ::Vector<TrilinosScalar> &values)
1589 Assert(indices.size() == values.size(),
1592 set(indices.size(), indices.data(), values.begin());
1600 const TrilinosScalar *values)
1606 if (last_action == Add)
1608 const int ierr = vector->GlobalAssemble(Add);
1612 if (last_action != Insert)
1613 last_action = Insert;
1615 for (
size_type i = 0; i < n_elements; ++i)
1617 const TrilinosWrappers::types::int_type row = indices[i];
1618 const TrilinosWrappers::types::int_type local_row =
1619 vector->Map().LID(row);
1620 if (local_row != -1)
1621 (*vector)[0][local_row] = values[i];
1624 const int ierr = vector->ReplaceGlobalValues(1, &row, &values[i]);
1639 Vector::add(
const std::vector<size_type> & indices,
1640 const std::vector<TrilinosScalar> &values)
1645 Assert(indices.size() == values.size(),
1648 add(indices.size(), indices.data(), values.data());
1654 Vector::add(
const std::vector<size_type> & indices,
1655 const ::Vector<TrilinosScalar> &values)
1660 Assert(indices.size() == values.size(),
1663 add(indices.size(), indices.data(), values.begin());
1671 const TrilinosScalar *values)
1677 if (last_action != Add)
1679 if (last_action == Insert)
1681 const int ierr = vector->GlobalAssemble(Insert);
1687 for (
size_type i = 0; i < n_elements; ++i)
1690 const TrilinosWrappers::types::int_type local_row = vector->Map().LID(
1691 static_cast<TrilinosWrappers::types::int_type>(row));
1692 if (local_row != -1)
1693 (*vector)[0][local_row] += values[i];
1694 else if (nonlocal_vector.get() ==
nullptr)
1696 const int ierr = vector->SumIntoGlobalValues(
1698 reinterpret_cast<const TrilinosWrappers::types::int_type *>(
1708 const TrilinosWrappers::types::int_type my_row =
1709 nonlocal_vector->Map().LID(
1710 static_cast<TrilinosWrappers::types::int_type>(row));
1713 "Attempted to write into off-processor vector entry "
1714 "that has not be specified as being writable upon "
1716 (*nonlocal_vector)[0][my_row] += values[i];
1724 inline Vector::size_type
1727 # ifndef DEAL_II_WITH_64BIT_INDICES
1728 return vector->Map().MaxAllGID() + 1 - vector->Map().MinAllGID();
1730 return vector->Map().MaxAllGID64() + 1 - vector->Map().MinAllGID64();
1736 inline Vector::size_type
1737 Vector::local_size()
const
1739 return vector->Map().NumMyElements();
1744 inline std::pair<Vector::size_type, Vector::size_type>
1745 Vector::local_range()
const
1747 # ifndef DEAL_II_WITH_64BIT_INDICES
1748 const TrilinosWrappers::types::int_type begin = vector->Map().MinMyGID();
1749 const TrilinosWrappers::types::int_type end =
1750 vector->Map().MaxMyGID() + 1;
1752 const TrilinosWrappers::types::int_type begin =
1753 vector->Map().MinMyGID64();
1754 const TrilinosWrappers::types::int_type end =
1755 vector->Map().MaxMyGID64() + 1;
1759 end - begin == vector->Map().NumMyElements(),
1761 "This function only makes sense if the elements that this "
1762 "vector stores on the current processor form a contiguous range. "
1763 "This does not appear to be the case for the current vector."));
1765 return std::make_pair(begin, end);
1772 Assert(vector->Map().SameAs(vec.vector->Map()),
1773 ExcDifferentParallelPartitioning());
1776 TrilinosScalar result;
1778 const int ierr = vector->Dot(*(vec.vector), &result);
1789 const TrilinosScalar
d = l2_norm();
1795 inline TrilinosScalar
1800 TrilinosScalar
mean;
1801 const int ierr = vector->MeanValue(&mean);
1809 inline TrilinosScalar
1812 TrilinosScalar min_value;
1813 const int ierr = vector->MinValue(&min_value);
1821 inline TrilinosScalar
1824 TrilinosScalar max_value;
1825 const int ierr = vector->MaxValue(&max_value);
1839 const int ierr = vector->Norm1(&d);
1853 const int ierr = vector->Norm2(&d);
1866 TrilinosScalar
norm = 0;
1867 TrilinosScalar
sum = 0;
1873 sum += std::pow(std::fabs((*vector)[0][i]), p);
1875 norm = std::pow(sum, static_cast<TrilinosScalar>(1. / p));
1891 const int ierr = vector->NormInf(&d);
1899 inline TrilinosScalar
1920 const int ierr = vector->Scale(a);
1933 const TrilinosScalar factor = 1. / a;
1937 const int ierr = vector->Scale(factor);
1949 Assert(vector->Map().SameAs(v.vector->Map()),
1950 ExcDifferentParallelPartitioning());
1952 const int ierr = vector->Update(1.0, *(v.vector), 1.0);
1964 Assert(vector->Map().SameAs(v.vector->Map()),
1965 ExcDifferentParallelPartitioning());
1967 const int ierr = vector->Update(-1.0, *(v.vector), 1.0);
1985 (*vector)[0][i] += s;
1996 Assert(local_size() == v.local_size(),
2001 const int ierr = vector->Update(a, *(v.vector), 1.);
2010 const TrilinosScalar b,
2016 Assert(local_size() == v.local_size(),
2018 Assert(local_size() ==
w.local_size(),
2024 const int ierr = vector->Update(a, *(v.vector), b, *(
w.vector), 1.);
2043 if (local_size() == v.local_size() && !v.has_ghost_elements())
2045 Assert(this->vector->Map().SameAs(v.vector->Map()) ==
true,
2046 ExcDifferentParallelPartitioning());
2047 const int ierr = vector->Update(1., *(v.vector), s);
2061 const TrilinosScalar a,
2073 if (local_size() == v.local_size() && !v.has_ghost_elements())
2075 Assert(this->vector->Map().SameAs(v.vector->Map()) ==
true,
2076 ExcDifferentParallelPartitioning());
2077 const int ierr = vector->Update(a, *(v.vector), s);
2085 this->add(tmp,
true);
2097 Assert(local_size() == factors.local_size(),
2100 const int ierr = vector->Multiply(1.0, *(factors.vector), *vector, 0.0);
2115 if (vector->Map().SameAs(v.vector->Map()) ==
false)
2117 this->
sadd(0., a, v);
2122 int ierr = vector->Update(a, *v.vector, 0.0);
2131 inline const Epetra_MultiVector &
2132 Vector::trilinos_vector()
const
2134 return static_cast<const Epetra_MultiVector &>(*vector);
2139 inline Epetra_FEVector &
2140 Vector::trilinos_vector()
2147 inline const Epetra_Map &
2148 Vector::vector_partitioner()
const
2151 return static_cast<const Epetra_Map &>(vector->Map());
2156 inline const Epetra_BlockMap &
2157 Vector::trilinos_partitioner()
const
2159 return vector->Map();
2164 inline const MPI_Comm &
2165 Vector::get_mpi_communicator()
const
2167 static MPI_Comm comm;
2169 # ifdef DEAL_II_WITH_MPI
2171 const Epetra_MpiComm *mpi_comm =
2172 dynamic_cast<const Epetra_MpiComm *>(&vector->Map().Comm());
2173 comm = mpi_comm->Comm();
2177 comm = MPI_COMM_SELF;
2184 template <
typename number>
2186 const ::Vector<number> &v,
2187 const MPI_Comm & communicator)
2191 owned_elements = parallel_partitioner;
2201 int ierr = vector->PutScalar(s);
2204 if (nonlocal_vector.get() !=
nullptr)
2206 ierr = nonlocal_vector->PutScalar(0.);
2223 namespace LinearOperatorImplementation
2236 template <
typename Matrix>
2240 bool omit_zeroing_entries)
2242 v.
reinit(matrix.locally_owned_range_indices(),
2243 matrix.get_mpi_communicator(),
2244 omit_zeroing_entries);
2247 template <
typename Matrix>
2251 bool omit_zeroing_entries)
2253 v.
reinit(matrix.locally_owned_domain_indices(),
2254 matrix.get_mpi_communicator(),
2255 omit_zeroing_entries);
2274 DEAL_II_NAMESPACE_CLOSE
2276 #endif // DEAL_II_WITH_TRILINOS
2280 #endif // dealii_trilinos_vector_h