16#ifndef dealii_trilinos_vector_h
17#define dealii_trilinos_vector_h
22#ifdef DEAL_II_WITH_TRILINOS
33# include <Epetra_ConfigDefs.h>
38# ifdef DEAL_II_WITH_MPI
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>
55 template <
typename Number>
56 class ReadWriteVector;
100 class VectorReference
113 VectorReference(
const VectorReference &) =
default;
126 const VectorReference &
127 operator=(
const VectorReference &r)
const;
133 operator=(
const VectorReference &r);
138 const VectorReference &
144 const VectorReference &
150 const VectorReference &
156 const VectorReference &
162 const VectorReference &
176 <<
"An error with error number " << arg1
177 <<
" occurred while calling a Trilinos function");
199# ifndef DEAL_II_WITH_64BIT_INDICES
204 gid(
const Epetra_BlockMap &map,
int i)
213 gid(
const Epetra_BlockMap &map,
int i)
454 const MPI_Comm &communicator = MPI_COMM_WORLD);
469 const MPI_Comm &communicator = MPI_COMM_WORLD);
487 const MPI_Comm &communicator = MPI_COMM_WORLD);
501 template <
typename Number>
503 const ::Vector<Number> &v,
504 const MPI_Comm & communicator = MPI_COMM_WORLD);
549 const bool omit_zeroing_entries =
false,
550 const bool allow_different_maps =
false);
576 const MPI_Comm &communicator = MPI_COMM_WORLD,
577 const bool omit_zeroing_entries =
false);
607 const MPI_Comm &communicator = MPI_COMM_WORLD,
608 const bool vector_writable =
false);
672 template <
typename Number>
776 std::pair<size_type, size_type>
976 std::vector<TrilinosScalar> &
values)
const;
1005 template <
typename ForwardIterator,
typename OutputIterator>
1008 const ForwardIterator indices_end,
1009 OutputIterator values_begin)
const;
1060 set(
const std::vector<size_type> & indices,
1061 const std::vector<TrilinosScalar> &
values);
1068 set(
const std::vector<size_type> & indices,
1069 const ::Vector<TrilinosScalar> &
values);
1086 add(
const std::vector<size_type> & indices,
1087 const std::vector<TrilinosScalar> &
values);
1094 add(
const std::vector<size_type> & indices,
1095 const ::Vector<TrilinosScalar> &
values);
1151 add(
const Vector &
V,
const bool allow_different_maps =
false);
1205 const Epetra_MultiVector &
1219 const Epetra_BlockMap &
1230 print(std::ostream & out,
1231 const unsigned int precision = 3,
1232 const bool scientific =
true,
1233 const bool across =
true)
const;
1275 <<
"An error with error number " << arg1
1276 <<
" occurred while calling a Trilinos function");
1287 <<
"You are trying to access element " << arg1
1288 <<
" of a distributed vector, but this element is not stored "
1289 <<
"on the current processor. Note: There are " << arg2
1290 <<
" elements stored "
1291 <<
"on the current processor from within the range [" << arg3 <<
","
1292 << arg4 <<
"] but Trilinos vectors need not store contiguous "
1293 <<
"ranges on each processor, and not every element in "
1294 <<
"this range may in fact be stored locally."
1296 <<
"A common source for this kind of problem is that you "
1297 <<
"are passing a 'fully distributed' vector into a function "
1298 <<
"that needs read access to vector elements that correspond "
1299 <<
"to degrees of freedom on ghost cells (or at least to "
1300 <<
"'locally active' degrees of freedom that are not also "
1301 <<
"'locally owned'). You need to pass a vector that has these "
1302 <<
"elements as ghost entries.");
1376 inline VectorReference::VectorReference(
MPI::Vector & vector,
1383 inline const VectorReference &
1384 VectorReference::operator=(
const VectorReference &r)
const
1397 inline VectorReference &
1398 VectorReference::operator=(
const VectorReference &r)
1407 inline const VectorReference &
1410 vector.
set(1, &index, &value);
1416 inline const VectorReference &
1419 vector.add(1, &index, &value);
1425 inline const VectorReference &
1429 vector.add(1, &index, &new_value);
1435 inline const VectorReference &
1439 vector.set(1, &index, &new_value);
1445 inline const VectorReference &
1449 vector.set(1, &index, &new_value);
1459 std::pair<size_type, size_type> range =
local_range();
1461 return ((index >= range.first) && (index < range.second));
1471 "The locally owned elements have not been properly initialized!"
1472 " This happens for example if this object has been initialized"
1473 " with exactly one overlapping IndexSet."));
1493 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)
1584 set(indices.size(), indices.data(),
values.data());
1590 Vector::set(
const std::vector<size_type> & indices,
1591 const ::Vector<TrilinosScalar> &
values)
1600 set(indices.size(), indices.data(),
values.begin());
1616 const int ierr =
vector->GlobalAssemble(Add);
1623 for (
size_type i = 0; i < n_elements; ++i)
1628 if (local_row != -1)
1629 (*vector)[0][local_row] =
values[i];
1632 const int ierr =
vector->ReplaceGlobalValues(1, &row, &
values[i]);
1647 Vector::add(
const std::vector<size_type> & indices,
1648 const std::vector<TrilinosScalar> &
values)
1656 add(indices.size(), indices.data(),
values.data());
1662 Vector::add(
const std::vector<size_type> & indices,
1663 const ::Vector<TrilinosScalar> &
values)
1671 add(indices.size(), indices.data(),
values.begin());
1689 const int ierr =
vector->GlobalAssemble(Insert);
1695 for (
size_type i = 0; i < n_elements; ++i)
1700 if (local_row != -1)
1701 (*vector)[0][local_row] +=
values[i];
1704 const int ierr =
vector->SumIntoGlobalValues(
1721 "Attempted to write into off-processor vector entry "
1722 "that has not be specified as being writable upon "
1724 (*nonlocal_vector)[0][my_row] +=
values[i];
1735# ifndef DEAL_II_WITH_64BIT_INDICES
1736 return vector->Map().MaxAllGID() + 1 -
vector->Map().MinAllGID();
1738 return vector->Map().MaxAllGID64() + 1 -
vector->Map().MinAllGID64();
1747 return vector->Map().NumMyElements();
1760 inline std::pair<Vector::size_type, Vector::size_type>
1763# ifndef DEAL_II_WITH_64BIT_INDICES
1766 vector->Map().MaxMyGID() + 1;
1769 vector->Map().MinMyGID64();
1771 vector->Map().MaxMyGID64() + 1;
1777 "This function only makes sense if the elements that this "
1778 "vector stores on the current processor form a contiguous range. "
1779 "This does not appear to be the case for the current vector."));
1794 const int ierr =
vector->Dot(*(vec.vector), &result);
1829 const int ierr =
vector->MinValue(&min_value);
1841 const int ierr =
vector->MaxValue(&max_value);
1855 const int ierr =
vector->Norm1(&
d);
1869 const int ierr =
vector->Norm2(&
d);
1907 const int ierr =
vector->NormInf(&
d);
1936 const int ierr =
vector->Scale(a);
1953 const int ierr =
vector->Scale(factor);
1968 const int ierr =
vector->Update(1.0, *(v.vector), 1.0);
1983 const int ierr =
vector->Update(-1.0, *(v.vector), 1.0);
2017 const int ierr =
vector->Update(a, *(v.vector), 1.);
2040 const int ierr =
vector->Update(a, *(v.vector),
b, *(
w.vector), 1.);
2061 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2063 const int ierr =
vector->Update(1., *(v.vector), s);
2091 Assert(this->
vector->Map().SameAs(v.vector->Map()) ==
true,
2093 const int ierr =
vector->Update(a, *(v.vector), s);
2101 this->
add(tmp,
true);
2116 const int ierr =
vector->Multiply(1.0, *(factors.vector), *
vector, 0.0);
2131 if (
vector->Map().SameAs(v.vector->Map()) ==
false)
2133 this->
sadd(0., a, v);
2138 int ierr =
vector->Update(a, *v.vector, 0.0);
2147 inline const Epetra_MultiVector &
2150 return static_cast<const Epetra_MultiVector &
>(*vector);
2155 inline Epetra_FEVector &
2163 inline const Epetra_BlockMap &
2176# ifdef DEAL_II_WITH_MPI
2178 const Epetra_MpiComm *mpi_comm =
2179 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
2180 comm = mpi_comm->Comm();
2184 comm = MPI_COMM_SELF;
2191 template <
typename number>
2193 const ::Vector<number> &v,
2208 int ierr =
vector->PutScalar(s);
2230 namespace LinearOperatorImplementation
2243 template <
typename Matrix>
2247 bool omit_zeroing_entries)
2250 matrix.get_mpi_communicator(),
2251 omit_zeroing_entries);
2254 template <
typename Matrix>
2258 bool omit_zeroing_entries)
2261 matrix.get_mpi_communicator(),
2262 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_NAMESPACE_OPEN
#define DEAL_II_DEPRECATED_EARLY
#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)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#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
Expression fabs(const Expression &x)
@ matrix
Contents is actually a matrix.
SparseMatrix< double > SparseMatrix
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 > b(const Tensor< 2, dim, Number > &F)
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
*braid_SplitCommworld & comm