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>
51 template <
typename Number>
52 class ReadWriteVector;
103 VectorReference(MPI::Vector &vector,
const size_type index);
109 VectorReference(
const VectorReference &) =
default;
122 const VectorReference &
123 operator=(
const VectorReference &r)
const;
129 operator=(
const VectorReference &r);
134 const VectorReference &
140 const VectorReference &
146 const VectorReference &
152 const VectorReference &
158 const VectorReference &
172 <<
"An error with error number " << arg1
173 <<
" occurred while calling a Trilinos function");
188 friend class ::TrilinosWrappers::MPI::Vector;
195# ifndef DEAL_II_WITH_64BIT_INDICES
200 gid(
const Epetra_BlockMap &map,
int i)
209 gid(
const Epetra_BlockMap &map,
int i)
449 const MPI_Comm &communicator = MPI_COMM_WORLD);
464 const MPI_Comm &communicator = MPI_COMM_WORLD);
482 const MPI_Comm &communicator = MPI_COMM_WORLD);
496 template <
typename Number>
498 const ::Vector<Number> &v,
499 const MPI_Comm & communicator = MPI_COMM_WORLD);
544 const bool omit_zeroing_entries =
false,
545 const bool allow_different_maps =
false);
571 const MPI_Comm &communicator = MPI_COMM_WORLD,
572 const bool omit_zeroing_entries =
false);
602 const MPI_Comm &communicator = MPI_COMM_WORLD,
603 const bool vector_writable =
false);
667 template <
typename Number>
690 const ::TrilinosWrappers::SparseMatrix &matrix,
771 std::pair<size_type, size_type>
974 std::vector<TrilinosScalar> & values)
const;
1003 template <
typename ForwardIterator,
typename OutputIterator>
1006 const ForwardIterator indices_end,
1007 OutputIterator values_begin)
const;
1058 set(
const std::vector<size_type> & indices,
1059 const std::vector<TrilinosScalar> &values);
1066 set(
const std::vector<size_type> & indices,
1067 const ::Vector<TrilinosScalar> &values);
1084 add(
const std::vector<size_type> & indices,
1085 const std::vector<TrilinosScalar> &values);
1092 add(
const std::vector<size_type> & indices,
1093 const ::Vector<TrilinosScalar> &values);
1149 add(
const Vector &V,
const bool allow_different_maps =
false);
1203 const Epetra_MultiVector &
1217 const Epetra_BlockMap &
1228 print(std::ostream & out,
1229 const unsigned int precision = 3,
1230 const bool scientific =
true,
1231 const bool across =
true)
const;
1273 <<
"An error with error number " << arg1
1274 <<
" occurred while calling a Trilinos function");
1285 <<
"You are trying to access element " << arg1
1286 <<
" of a distributed vector, but this element is not stored "
1287 <<
"on the current processor. Note: There are " << arg2
1288 <<
" elements stored "
1289 <<
"on the current processor from within the range [" << arg3 <<
','
1290 << arg4 <<
"] but Trilinos vectors need not store contiguous "
1291 <<
"ranges on each processor, and not every element in "
1292 <<
"this range may in fact be stored locally."
1294 <<
"A common source for this kind of problem is that you "
1295 <<
"are passing a 'fully distributed' vector into a function "
1296 <<
"that needs read access to vector elements that correspond "
1297 <<
"to degrees of freedom on ghost cells (or at least to "
1298 <<
"'locally active' degrees of freedom that are not also "
1299 <<
"'locally owned'). You need to pass a vector that has these "
1300 <<
"elements as ghost entries.");
1374 inline VectorReference::VectorReference(
MPI::Vector & vector,
1375 const size_type index)
1381 inline const VectorReference &
1382 VectorReference::operator=(
const VectorReference &r)
const
1395 inline VectorReference &
1396 VectorReference::operator=(
const VectorReference &r)
1405 inline const VectorReference &
1408 vector.set(1, &index, &value);
1414 inline const VectorReference &
1417 vector.add(1, &index, &value);
1423 inline const VectorReference &
1427 vector.add(1, &index, &new_value);
1433 inline const VectorReference &
1437 vector.set(1, &index, &new_value);
1443 inline const VectorReference &
1447 vector.set(1, &index, &new_value);
1457 std::pair<size_type, size_type> range =
local_range();
1459 return ((index >= range.first) && (index < range.second));
1469 "The locally owned elements have not been properly initialized!"
1470 " This happens for example if this object has been initialized"
1471 " with exactly one overlapping IndexSet."));
1491 inline internal::VectorReference
1499 inline internal::VectorReference
1517 std::vector<TrilinosScalar> & values)
const
1519 for (
size_type i = 0; i < indices.size(); ++i)
1520 values[i] =
operator()(indices[i]);
1525 template <
typename ForwardIterator,
typename OutputIterator>
1528 const ForwardIterator indices_end,
1529 OutputIterator values_begin)
const
1531 while (indices_begin != indices_end)
1574 Vector::set(
const std::vector<size_type> & indices,
1575 const std::vector<TrilinosScalar> &values)
1583 set(indices.size(), indices.data(),
values.data());
1589 Vector::set(
const std::vector<size_type> & indices,
1590 const ::Vector<TrilinosScalar> &values)
1598 set(indices.size(), indices.data(),
values.begin());
1605 const size_type * indices,
1614 const int ierr =
vector->GlobalAssemble(Add);
1621 for (
size_type i = 0; i < n_elements; ++i)
1626 if (local_row != -1)
1627 (*vector)[0][local_row] =
values[i];
1630 const int ierr =
vector->ReplaceGlobalValues(1, &row, &values[i]);
1645 Vector::add(
const std::vector<size_type> & indices,
1646 const std::vector<TrilinosScalar> &values)
1653 add(indices.size(), indices.data(),
values.data());
1659 Vector::add(
const std::vector<size_type> & indices,
1660 const ::Vector<TrilinosScalar> &values)
1667 add(indices.size(), indices.data(),
values.begin());
1674 const size_type * indices,
1685 const int ierr =
vector->GlobalAssemble(Insert);
1691 for (
size_type i = 0; i < n_elements; ++i)
1696 if (local_row != -1)
1697 (*vector)[0][local_row] +=
values[i];
1700 const int ierr =
vector->SumIntoGlobalValues(
1717 "Attempted to write into off-processor vector entry "
1718 "that has not be specified as being writable upon "
1720 (*nonlocal_vector)[0][my_row] +=
values[i];
1731# ifndef DEAL_II_WITH_64BIT_INDICES
1732 return vector->Map().MaxAllGID() + 1 -
vector->Map().MinAllGID();
1734 return vector->Map().MaxAllGID64() + 1 -
vector->Map().MinAllGID64();
1743 return vector->Map().NumMyElements();
1756 inline std::pair<Vector::size_type, Vector::size_type>
1759# ifndef DEAL_II_WITH_64BIT_INDICES
1762 vector->Map().MaxMyGID() + 1;
1765 vector->Map().MinMyGID64();
1767 vector->Map().MaxMyGID64() + 1;
1773 "This function only makes sense if the elements that this "
1774 "vector stores on the current processor form a contiguous range. "
1775 "This does not appear to be the case for the current vector."));
1791 const int ierr =
vector->Dot(*(vec.vector), &result);
1814 const int ierr =
vector->MeanValue(&mean);
1826 const int ierr =
vector->MinValue(&min_value);
1838 const int ierr =
vector->MaxValue(&max_value);
1852 const int ierr =
vector->Norm1(&d);
1866 const int ierr =
vector->Norm2(&d);
1904 const int ierr =
vector->NormInf(&d);
1933 const int ierr =
vector->Scale(a);
1950 const int ierr =
vector->Scale(factor);
1965 const int ierr =
vector->Update(1.0, *(v.vector), 1.0);
1980 const int ierr =
vector->Update(-1.0, *(v.vector), 1.0);
2013 const int ierr =
vector->Update(a, *(v.vector), 1.);
2034 const int ierr =
vector->Update(a, *(v.vector), b, *(
w.vector), 1.);
2056 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2058 const int ierr =
vector->Update(1., *(v.vector), s);
2087 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2089 const int ierr =
vector->Update(a, *(v.vector), s);
2097 this->
add(tmp,
true);
2111 const int ierr =
vector->Multiply(1.0, *(factors.vector), *
vector, 0.0);
2126 if (
vector->Map().SameAs(v.vector->Map()) ==
false)
2128 this->
sadd(0., a, v);
2133 int ierr =
vector->Update(a, *v.vector, 0.0);
2142 inline const Epetra_MultiVector &
2145 return static_cast<const Epetra_MultiVector &
>(*vector);
2150 inline Epetra_FEVector &
2158 inline const Epetra_BlockMap &
2171 const Epetra_MpiComm *mpi_comm =
2172 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
2173 comm = mpi_comm->Comm();
2178 template <
typename number>
2180 const ::Vector<number> &v,
2195 int ierr =
vector->PutScalar(s);
2217 namespace LinearOperatorImplementation
2230 template <
typename Matrix>
2234 bool omit_zeroing_entries)
2236 v.
reinit(matrix.locally_owned_range_indices(),
2237 matrix.get_mpi_communicator(),
2238 omit_zeroing_entries);
2241 template <
typename Matrix>
2245 bool omit_zeroing_entries)
2247 v.
reinit(matrix.locally_owned_domain_indices(),
2248 matrix.get_mpi_communicator(),
2249 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
#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)
#define Assert(cond, exc)
#define AssertIsFinite(number)
#define AssertDimension(dim1, dim2)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
real_type l1_norm() const
Vector & operator/=(const TrilinosScalar factor)
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)
const MPI_Comm & get_mpi_communicator() const
static void reinit_domain_vector(const Matrix &matrix, TrilinosWrappers::MPI::Vector &v, bool omit_zeroing_entries)
std::unique_ptr< Epetra_MultiVector > nonlocal_vector
void add(const TrilinosScalar s)
void update_ghost_values() const
static ::ExceptionBase & ExcDifferentParallelPartitioning()
const Epetra_BlockMap & trilinos_partitioner() const
int gid(const Epetra_BlockMap &map, int i)
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
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
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
void compress(::VectorOperation::values operation)
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)
static void reinit_range_vector(const Matrix &matrix, TrilinosWrappers::MPI::Vector &v, bool omit_zeroing_entries)
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
static ::ExceptionBase & ExcTrilinosError(int arg1)
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
Vector(const IndexSet ¶llel_partitioning, const ::Vector< Number > &v, const MPI_Comm &communicator=MPI_COMM_WORLD)
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
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)
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