15#ifndef dealii_trilinos_vector_h
16#define dealii_trilinos_vector_h
21#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;
100 class VectorReference
109 VectorReference(
MPI::Vector &vector,
const size_type index);
115 VectorReference(
const VectorReference &) =
default;
128 const VectorReference &
129 operator=(
const VectorReference &r)
const;
135 operator=(
const VectorReference &r);
140 const VectorReference &
146 const VectorReference &
152 const VectorReference &
158 const VectorReference &
164 const VectorReference &
178 <<
"An error with error number " << arg1
179 <<
" occurred while calling a Trilinos function");
190 const size_type index;
194 friend class ::TrilinosWrappers::MPI::Vector;
201# ifndef DEAL_II_WITH_64BIT_INDICES
206 gid(
const Epetra_BlockMap &map,
int i)
215 gid(
const Epetra_BlockMap &map,
int i)
503 template <
typename Number>
505 const ::Vector<Number> &v,
637 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
727 template <
typename Number>
750 const ::TrilinosWrappers::SparseMatrix &matrix,
826 std::pair<size_type, size_type>
1029 std::vector<TrilinosScalar> &values)
const;
1065 template <
typename ForwardIterator,
typename OutputIterator>
1118 set(
const std::vector<size_type> &indices,
1119 const std::vector<TrilinosScalar> &values);
1126 set(
const std::vector<size_type> &indices,
1127 const ::Vector<TrilinosScalar> &values);
1144 add(
const std::vector<size_type> &indices,
1145 const std::vector<TrilinosScalar> &values);
1152 add(
const std::vector<size_type> &indices,
1153 const ::Vector<TrilinosScalar> &values);
1263 const Epetra_MultiVector &
1288 print(std::ostream &out,
1289 const unsigned int precision = 3,
1290 const bool scientific =
true,
1291 const bool across =
true)
const;
1332 <<
"An error with error number " <<
arg1
1333 <<
" occurred while calling a Trilinos function");
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."
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.");
1433 inline VectorReference::VectorReference(
MPI::Vector &vector,
1434 const size_type index)
1440 inline const VectorReference &
1441 VectorReference::operator=(
const VectorReference &r)
const
1454 inline VectorReference &
1455 VectorReference::operator=(
const VectorReference &r)
1464 inline const VectorReference &
1467 vector.set(1, &index, &value);
1473 inline const VectorReference &
1476 vector.add(1, &index, &value);
1482 inline const VectorReference &
1486 vector.add(1, &index, &new_value);
1492 inline const VectorReference &
1496 vector.set(1, &index, &new_value);
1502 inline const VectorReference &
1506 vector.set(1, &index, &new_value);
1518 return ((index >=
range.first) && (index <
range.second));
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."));
1550 inline internal::VectorReference
1558 inline internal::VectorReference
1576 std::vector<TrilinosScalar> &values)
const
1578 for (
size_type i = 0; i < indices.size(); ++i)
1579 values[i] =
operator()(indices[i]);
1589 for (
unsigned int i = 0; i < indices.
size(); ++i)
1592 elements[i] = (*this)[indices[i]];
1598 template <
typename ForwardIterator,
typename OutputIterator>
1601 const ForwardIterator indices_end,
1602 OutputIterator values_begin)
const
1647 Vector::set(
const std::vector<size_type> &indices,
1648 const std::vector<TrilinosScalar> &values)
1656 set(indices.size(), indices.data(),
values.data());
1662 Vector::set(
const std::vector<size_type> &indices,
1663 const ::Vector<TrilinosScalar> &values)
1671 set(indices.size(), indices.data(),
values.begin());
1678 const size_type *indices,
1694 for (
size_type i = 0; i < n_elements; ++i)
1703 const int ierr =
vector->ReplaceGlobalValues(1, &row, &values[i]);
1718 Vector::add(
const std::vector<size_type> &indices,
1719 const std::vector<TrilinosScalar> &values)
1726 add(indices.size(), indices.data(),
values.data());
1732 Vector::add(
const std::vector<size_type> &indices,
1733 const ::Vector<TrilinosScalar> &values)
1740 add(indices.size(), indices.data(),
values.begin());
1747 const size_type *indices,
1764 for (
size_type i = 0; i < n_elements; ++i)
1773 const int ierr =
vector->SumIntoGlobalValues(
1790 "Attempted to write into off-processor vector entry "
1791 "that has not be specified as being writable upon "
1804# ifndef DEAL_II_WITH_64BIT_INDICES
1805 return vector->Map().MaxAllGID() + 1 -
vector->Map().MinAllGID();
1807 return vector->Map().MaxAllGID64() + 1 -
vector->Map().MinAllGID64();
1821 inline std::pair<Vector::size_type, Vector::size_type>
1824# ifndef DEAL_II_WITH_64BIT_INDICES
1827 vector->Map().MaxMyGID() + 1;
1830 vector->Map().MinMyGID64();
1832 vector->Map().MaxMyGID64() + 1;
1838 "This function only makes sense if the elements that this "
1839 "vector stores on the current processor form a contiguous range. "
1840 "This does not appear to be the case for the current vector."));
2030 const int ierr =
vector->Update(1.0, *(v.vector), 1.0);
2045 const int ierr =
vector->Update(-1.0, *(v.vector), 1.0);
2078 const int ierr =
vector->Update(a, *(v.vector), 1.);
2099 const int ierr =
vector->Update(a, *(v.vector), b, *(
w.vector), 1.);
2121 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2123 const int ierr =
vector->Update(1., *(v.vector), s);
2152 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2154 const int ierr =
vector->Update(a, *(v.vector), s);
2162 this->
add(tmp,
true);
2191 if (
vector->Map().SameAs(v.vector->Map()) ==
false)
2193 this->
sadd(0., a, v);
2198 int ierr =
vector->Update(a, *v.vector, 0.0);
2207 inline const Epetra_MultiVector &
2210 return static_cast<const Epetra_MultiVector &
>(*vector);
2215 inline Epetra_FEVector &
2223 inline const Epetra_BlockMap &
2234 const Epetra_MpiComm *mpi_comm =
2235 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
2236 return mpi_comm->Comm();
2241 template <
typename number>
2243 const ::Vector<number> &v,
2280 namespace LinearOperatorImplementation
2293 template <
typename Matrix>
2297 bool omit_zeroing_entries)
2299 v.
reinit(matrix.locally_owned_range_indices(),
2300 matrix.get_mpi_communicator(),
2301 omit_zeroing_entries);
2304 template <
typename Matrix>
2308 bool omit_zeroing_entries)
2310 v.
reinit(matrix.locally_owned_domain_indices(),
2311 matrix.get_mpi_communicator(),
2312 omit_zeroing_entries);
size_type n_elements() 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)
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)
void update_ghost_values() const
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
virtual void extract_subvector_to(const ArrayView< const size_type > &indices, ArrayView< TrilinosScalar > &elements) const override
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
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)
bool is_non_negative() const
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
typename numbers::NumberTraits< Number >::real_type real_type
bool has_ghost_elements() const
const value_type * const_iterator
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
#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)
#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)
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