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)
455 const MPI_Comm communicator = MPI_COMM_WORLD);
470 const MPI_Comm communicator = MPI_COMM_WORLD);
488 const MPI_Comm communicator = MPI_COMM_WORLD);
502 template <
typename Number>
504 const ::Vector<Number> &v,
505 const MPI_Comm communicator = MPI_COMM_WORLD);
554 const bool omit_zeroing_entries =
false,
555 const bool allow_different_maps =
false);
581 const MPI_Comm communicator = MPI_COMM_WORLD,
582 const bool omit_zeroing_entries =
false);
620 const IndexSet &locally_relevant_or_ghost_entries,
621 const MPI_Comm communicator = MPI_COMM_WORLD,
622 const bool vector_writable =
false);
636 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
637 const bool make_ghosted =
true,
638 const bool vector_writable =
false);
726 template <
typename Number>
749 const ::TrilinosWrappers::SparseMatrix &matrix,
825 std::pair<size_type, size_type>
1028 std::vector<TrilinosScalar> &values)
const;
1064 template <
typename ForwardIterator,
typename OutputIterator>
1067 const ForwardIterator indices_end,
1068 OutputIterator values_begin)
const;
1117 set(
const std::vector<size_type> &indices,
1118 const std::vector<TrilinosScalar> &values);
1125 set(
const std::vector<size_type> &indices,
1126 const ::Vector<TrilinosScalar> &values);
1143 add(
const std::vector<size_type> &indices,
1144 const std::vector<TrilinosScalar> &values);
1151 add(
const std::vector<size_type> &indices,
1152 const ::Vector<TrilinosScalar> &values);
1208 add(
const Vector &V,
const bool allow_different_maps =
false);
1262 const Epetra_MultiVector &
1276 const Epetra_BlockMap &
1287 print(std::ostream &out,
1288 const unsigned int precision = 3,
1289 const bool scientific =
true,
1290 const bool across =
true)
const;
1331 <<
"An error with error number " << arg1
1332 <<
" occurred while calling a Trilinos function");
1343 <<
"You are trying to access element " << arg1
1344 <<
" of a distributed vector, but this element is not stored "
1345 <<
"on the current processor. Note: There are " << arg2
1346 <<
" elements stored "
1347 <<
"on the current processor from within the range [" << arg3 <<
','
1348 << arg4 <<
"] but Trilinos vectors need not store contiguous "
1349 <<
"ranges on each processor, and not every element in "
1350 <<
"this range may in fact be stored locally."
1352 <<
"A common source for this kind of problem is that you "
1353 <<
"are passing a 'fully distributed' vector into a function "
1354 <<
"that needs read access to vector elements that correspond "
1355 <<
"to degrees of freedom on ghost cells (or at least to "
1356 <<
"'locally active' degrees of freedom that are not also "
1357 <<
"'locally owned'). You need to pass a vector that has these "
1358 <<
"elements as ghost entries.");
1432 inline VectorReference::VectorReference(
MPI::Vector &vector,
1433 const size_type index)
1439 inline const VectorReference &
1440 VectorReference::operator=(
const VectorReference &r)
const
1453 inline VectorReference &
1454 VectorReference::operator=(
const VectorReference &r)
1463 inline const VectorReference &
1466 vector.set(1, &index, &value);
1472 inline const VectorReference &
1475 vector.add(1, &index, &value);
1481 inline const VectorReference &
1485 vector.add(1, &index, &new_value);
1491 inline const VectorReference &
1495 vector.set(1, &index, &new_value);
1501 inline const VectorReference &
1505 vector.set(1, &index, &new_value);
1515 std::pair<size_type, size_type> range =
local_range();
1517 return ((index >= range.first) && (index < range.second));
1527 "The locally owned elements have not been properly initialized!"
1528 " This happens for example if this object has been initialized"
1529 " with exactly one overlapping IndexSet."));
1549 inline internal::VectorReference
1557 inline internal::VectorReference
1575 std::vector<TrilinosScalar> &values)
const
1577 for (size_type i = 0; i < indices.size(); ++i)
1578 values[i] =
operator()(indices[i]);
1588 for (
unsigned int i = 0; i < indices.
size(); ++i)
1591 elements[i] = (*this)[indices[i]];
1597 template <
typename ForwardIterator,
typename OutputIterator>
1600 const ForwardIterator indices_end,
1601 OutputIterator values_begin)
const
1603 while (indices_begin != indices_end)
1646 Vector::set(
const std::vector<size_type> &indices,
1647 const std::vector<TrilinosScalar> &values)
1655 set(indices.size(), indices.data(), values.data());
1661 Vector::set(
const std::vector<size_type> &indices,
1662 const ::Vector<TrilinosScalar> &values)
1670 set(indices.size(), indices.data(), values.begin());
1677 const size_type *indices,
1686 const int ierr =
vector->GlobalAssemble(Add);
1693 for (size_type i = 0; i < n_elements; ++i)
1698 if (local_row != -1)
1699 (*vector)[0][local_row] = values[i];
1702 const int ierr =
vector->ReplaceGlobalValues(1, &row, &values[i]);
1717 Vector::add(
const std::vector<size_type> &indices,
1718 const std::vector<TrilinosScalar> &values)
1725 add(indices.size(), indices.data(), values.data());
1731 Vector::add(
const std::vector<size_type> &indices,
1732 const ::Vector<TrilinosScalar> &values)
1739 add(indices.size(), indices.data(), values.begin());
1746 const size_type *indices,
1757 const int ierr =
vector->GlobalAssemble(Insert);
1763 for (size_type i = 0; i < n_elements; ++i)
1765 const size_type row = indices[i];
1768 if (local_row != -1)
1769 (*vector)[0][local_row] += values[i];
1772 const int ierr =
vector->SumIntoGlobalValues(
1789 "Attempted to write into off-processor vector entry "
1790 "that has not be specified as being writable upon "
1792 (*nonlocal_vector)[0][my_row] += values[i];
1803# ifndef DEAL_II_WITH_64BIT_INDICES
1804 return vector->Map().MaxAllGID() + 1 -
vector->Map().MinAllGID();
1806 return vector->Map().MaxAllGID64() + 1 -
vector->Map().MinAllGID64();
1820 inline std::pair<Vector::size_type, Vector::size_type>
1823# ifndef DEAL_II_WITH_64BIT_INDICES
1826 vector->Map().MaxMyGID() + 1;
1829 vector->Map().MinMyGID64();
1831 vector->Map().MaxMyGID64() + 1;
1837 "This function only makes sense if the elements that this "
1838 "vector stores on the current processor form a contiguous range. "
1839 "This does not appear to be the case for the current vector."));
1855 const int ierr =
vector->Dot(*(vec.vector), &result);
1878 const int ierr =
vector->MeanValue(&mean);
1890 const int ierr =
vector->MinValue(&min_value);
1902 const int ierr =
vector->MaxValue(&max_value);
1916 const int ierr =
vector->Norm1(&d);
1930 const int ierr =
vector->Norm2(&d);
1949 for (size_type i = 0; i < n_local; ++i)
1968 const int ierr =
vector->NormInf(&d);
1997 const int ierr =
vector->Scale(a);
2014 const int ierr =
vector->Scale(factor);
2029 const int ierr =
vector->Update(1.0, *(v.vector), 1.0);
2044 const int ierr =
vector->Update(-1.0, *(v.vector), 1.0);
2061 for (size_type i = 0; i < n_local; ++i)
2077 const int ierr =
vector->Update(a, *(v.vector), 1.);
2098 const int ierr =
vector->Update(a, *(v.vector), b, *(
w.vector), 1.);
2120 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2122 const int ierr =
vector->Update(1., *(v.vector), s);
2151 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2153 const int ierr =
vector->Update(a, *(v.vector), s);
2161 this->
add(tmp,
true);
2175 const int ierr =
vector->Multiply(1.0, *(factors.vector), *
vector, 0.0);
2190 if (
vector->Map().SameAs(v.vector->Map()) ==
false)
2192 this->
sadd(0., a, v);
2197 int ierr =
vector->Update(a, *v.vector, 0.0);
2206 inline const Epetra_MultiVector &
2209 return static_cast<const Epetra_MultiVector &
>(*vector);
2214 inline Epetra_FEVector &
2222 inline const Epetra_BlockMap &
2233 const Epetra_MpiComm *mpi_comm =
2234 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
2235 return mpi_comm->Comm();
2240 template <
typename number>
2242 const ::Vector<number> &v,
2257 int ierr =
vector->PutScalar(s);
2279 namespace LinearOperatorImplementation
2292 template <
typename Matrix>
2296 bool omit_zeroing_entries)
2298 v.
reinit(matrix.locally_owned_range_indices(),
2299 matrix.get_mpi_communicator(),
2300 omit_zeroing_entries);
2303 template <
typename Matrix>
2307 bool omit_zeroing_entries)
2309 v.
reinit(matrix.locally_owned_domain_indices(),
2310 matrix.get_mpi_communicator(),
2311 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)
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