16#ifndef dealii_trilinos_vector_h
17#define dealii_trilinos_vector_h
22#ifdef DEAL_II_WITH_TRILINOS
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>
50 template <
typename Number>
51 class ReadWriteVector;
102 VectorReference(MPI::Vector &vector,
const size_type index);
108 VectorReference(
const VectorReference &) =
default;
121 const VectorReference &
122 operator=(
const VectorReference &r)
const;
128 operator=(
const VectorReference &r);
133 const VectorReference &
139 const VectorReference &
145 const VectorReference &
151 const VectorReference &
157 const VectorReference &
171 <<
"An error with error number " << arg1
172 <<
" occurred while calling a Trilinos function");
187 friend class ::TrilinosWrappers::MPI::Vector;
194# ifndef DEAL_II_WITH_64BIT_INDICES
199 gid(
const Epetra_BlockMap &map,
int i)
208 gid(
const Epetra_BlockMap &map,
int i)
448 const MPI_Comm communicator = MPI_COMM_WORLD);
463 const MPI_Comm communicator = MPI_COMM_WORLD);
481 const MPI_Comm communicator = MPI_COMM_WORLD);
495 template <
typename Number>
497 const ::Vector<Number> &v,
498 const MPI_Comm communicator = MPI_COMM_WORLD);
547 const bool omit_zeroing_entries =
false,
548 const bool allow_different_maps =
false);
574 const MPI_Comm communicator = MPI_COMM_WORLD,
575 const bool omit_zeroing_entries =
false);
613 const IndexSet &locally_relevant_or_ghost_entries,
614 const MPI_Comm communicator = MPI_COMM_WORLD,
615 const bool vector_writable =
false);
629 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
630 const bool make_ghosted =
true,
631 const bool vector_writable =
false);
695 template <
typename Number>
718 const ::TrilinosWrappers::SparseMatrix &matrix,
734 DEAL_II_DEPRECATED_EARLY
810 std::pair<size_type, size_type>
1013 std::vector<TrilinosScalar> & values)
const;
1042 template <
typename ForwardIterator,
typename OutputIterator>
1045 const ForwardIterator indices_end,
1046 OutputIterator values_begin)
const;
1097 set(
const std::vector<size_type> & indices,
1098 const std::vector<TrilinosScalar> &values);
1105 set(
const std::vector<size_type> & indices,
1106 const ::Vector<TrilinosScalar> &values);
1123 add(
const std::vector<size_type> & indices,
1124 const std::vector<TrilinosScalar> &values);
1131 add(
const std::vector<size_type> & indices,
1132 const ::Vector<TrilinosScalar> &values);
1188 add(
const Vector &V,
const bool allow_different_maps =
false);
1242 const Epetra_MultiVector &
1256 const Epetra_BlockMap &
1267 print(std::ostream & out,
1268 const unsigned int precision = 3,
1269 const bool scientific =
true,
1270 const bool across =
true)
const;
1311 <<
"An error with error number " << arg1
1312 <<
" occurred while calling a Trilinos function");
1323 <<
"You are trying to access element " << arg1
1324 <<
" of a distributed vector, but this element is not stored "
1325 <<
"on the current processor. Note: There are " << arg2
1326 <<
" elements stored "
1327 <<
"on the current processor from within the range [" << arg3 <<
','
1328 << arg4 <<
"] but Trilinos vectors need not store contiguous "
1329 <<
"ranges on each processor, and not every element in "
1330 <<
"this range may in fact be stored locally."
1332 <<
"A common source for this kind of problem is that you "
1333 <<
"are passing a 'fully distributed' vector into a function "
1334 <<
"that needs read access to vector elements that correspond "
1335 <<
"to degrees of freedom on ghost cells (or at least to "
1336 <<
"'locally active' degrees of freedom that are not also "
1337 <<
"'locally owned'). You need to pass a vector that has these "
1338 <<
"elements as ghost entries.");
1412 inline VectorReference::VectorReference(
MPI::Vector & vector,
1413 const size_type index)
1419 inline const VectorReference &
1420 VectorReference::operator=(
const VectorReference &r)
const
1433 inline VectorReference &
1434 VectorReference::operator=(
const VectorReference &r)
1443 inline const VectorReference &
1446 vector.set(1, &index, &value);
1452 inline const VectorReference &
1455 vector.add(1, &index, &value);
1461 inline const VectorReference &
1465 vector.add(1, &index, &new_value);
1471 inline const VectorReference &
1475 vector.set(1, &index, &new_value);
1481 inline const VectorReference &
1485 vector.set(1, &index, &new_value);
1495 std::pair<size_type, size_type> range =
local_range();
1497 return ((index >= range.first) && (index < range.second));
1507 "The locally owned elements have not been properly initialized!"
1508 " This happens for example if this object has been initialized"
1509 " with exactly one overlapping IndexSet."));
1529 inline internal::VectorReference
1537 inline internal::VectorReference
1555 std::vector<TrilinosScalar> & values)
const
1557 for (
size_type i = 0; i < indices.size(); ++i)
1558 values[i] =
operator()(indices[i]);
1563 template <
typename ForwardIterator,
typename OutputIterator>
1566 const ForwardIterator indices_end,
1567 OutputIterator values_begin)
const
1569 while (indices_begin != indices_end)
1612 Vector::set(
const std::vector<size_type> & indices,
1613 const std::vector<TrilinosScalar> &values)
1621 set(indices.size(), indices.data(),
values.data());
1627 Vector::set(
const std::vector<size_type> & indices,
1628 const ::Vector<TrilinosScalar> &values)
1636 set(indices.size(), indices.data(),
values.begin());
1643 const size_type * indices,
1652 const int ierr =
vector->GlobalAssemble(Add);
1659 for (
size_type i = 0; i < n_elements; ++i)
1664 if (local_row != -1)
1665 (*vector)[0][local_row] =
values[i];
1668 const int ierr =
vector->ReplaceGlobalValues(1, &row, &values[i]);
1683 Vector::add(
const std::vector<size_type> & indices,
1684 const std::vector<TrilinosScalar> &values)
1691 add(indices.size(), indices.data(),
values.data());
1697 Vector::add(
const std::vector<size_type> & indices,
1698 const ::Vector<TrilinosScalar> &values)
1705 add(indices.size(), indices.data(),
values.begin());
1712 const size_type * indices,
1723 const int ierr =
vector->GlobalAssemble(Insert);
1729 for (
size_type i = 0; i < n_elements; ++i)
1734 if (local_row != -1)
1735 (*vector)[0][local_row] +=
values[i];
1738 const int ierr =
vector->SumIntoGlobalValues(
1755 "Attempted to write into off-processor vector entry "
1756 "that has not be specified as being writable upon "
1758 (*nonlocal_vector)[0][my_row] +=
values[i];
1769# ifndef DEAL_II_WITH_64BIT_INDICES
1770 return vector->Map().MaxAllGID() + 1 -
vector->Map().MinAllGID();
1772 return vector->Map().MaxAllGID64() + 1 -
vector->Map().MinAllGID64();
1781 return vector->Map().NumMyElements();
1794 inline std::pair<Vector::size_type, Vector::size_type>
1797# ifndef DEAL_II_WITH_64BIT_INDICES
1800 vector->Map().MaxMyGID() + 1;
1803 vector->Map().MinMyGID64();
1805 vector->Map().MaxMyGID64() + 1;
1811 "This function only makes sense if the elements that this "
1812 "vector stores on the current processor form a contiguous range. "
1813 "This does not appear to be the case for the current vector."));
1829 const int ierr =
vector->Dot(*(vec.vector), &result);
1852 const int ierr =
vector->MeanValue(&mean);
1864 const int ierr =
vector->MinValue(&min_value);
1876 const int ierr =
vector->MaxValue(&max_value);
1890 const int ierr =
vector->Norm1(&d);
1904 const int ierr =
vector->Norm2(&d);
1942 const int ierr =
vector->NormInf(&d);
1971 const int ierr =
vector->Scale(a);
1988 const int ierr =
vector->Scale(factor);
2003 const int ierr =
vector->Update(1.0, *(v.vector), 1.0);
2018 const int ierr =
vector->Update(-1.0, *(v.vector), 1.0);
2051 const int ierr =
vector->Update(a, *(v.vector), 1.);
2072 const int ierr =
vector->Update(a, *(v.vector), b, *(
w.vector), 1.);
2094 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2096 const int ierr =
vector->Update(1., *(v.vector), s);
2125 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2127 const int ierr =
vector->Update(a, *(v.vector), s);
2135 this->
add(tmp,
true);
2149 const int ierr =
vector->Multiply(1.0, *(factors.vector), *
vector, 0.0);
2164 if (
vector->Map().SameAs(v.vector->Map()) ==
false)
2166 this->
sadd(0., a, v);
2171 int ierr =
vector->Update(a, *v.vector, 0.0);
2180 inline const Epetra_MultiVector &
2183 return static_cast<const Epetra_MultiVector &
>(*vector);
2188 inline Epetra_FEVector &
2196 inline const Epetra_BlockMap &
2207 const Epetra_MpiComm *mpi_comm =
2208 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
2209 return mpi_comm->Comm();
2214 template <
typename number>
2216 const ::Vector<number> &v,
2231 int ierr =
vector->PutScalar(s);
2253 namespace LinearOperatorImplementation
2266 template <
typename Matrix>
2270 bool omit_zeroing_entries)
2272 v.
reinit(matrix.locally_owned_range_indices(),
2273 matrix.get_mpi_communicator(),
2274 omit_zeroing_entries);
2277 template <
typename Matrix>
2281 bool omit_zeroing_entries)
2283 v.
reinit(matrix.locally_owned_domain_indices(),
2284 matrix.get_mpi_communicator(),
2285 omit_zeroing_entries);
size_type n_elements() const
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
real_type l1_norm() const
Vector & operator/=(const TrilinosScalar factor)
void compress(VectorOperation::values operation)
Vector(const IndexSet ¶llel_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)
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)
void update_ghost_values() const
MPI_Comm get_mpi_communicator() const
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
bool in_local_range(const size_type index) const
size_type local_size() 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)
real_type l2_norm() const
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)
bool has_ghost_elements() const
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
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)
Epetra_CombineMode last_action
reference operator[](const size_type index)
const value_type * const_iterator
TrilinosScalar value_type
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)
bool is_non_negative() const
::types::global_dof_index size_type
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)
void swap(Vector &u, Vector &v)
TrilinosScalar max() const
Vector & operator=(const ::Vector< Number > &v)
const_iterator end() const
typename numbers::NumberTraits< Number >::real_type real_type
bool has_ghost_elements() const
const value_type * const_iterator
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
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#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)
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)
types::global_dof_index size_type
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
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