18#ifdef DEAL_II_WITH_TRILINOS
27# include <boost/io/ios_state.hpp>
29# include <Epetra_Export.h>
30# include <Epetra_Import.h>
31# include <Epetra_Vector.h>
53 vector.vector->Map().LID(
59 vector.vector->Map().MinMyGID(),
60 vector.vector->Map().MaxMyGID()));
63 return (*(vector.vector))[0][local_index];
75 , vector(new Epetra_FEVector(
76 Epetra_Map(0, 0, 0,
Utilities::Trilinos::comm_self())))
85 reinit(parallel_partitioning, communicator);
121 vector = std::make_unique<Epetra_FEVector>(
133 reinit(local, ghost, communicator,
false);
143 Epetra_Map map(0, 0, Epetra_MpiComm(MPI_COMM_SELF));
146 vector = std::make_unique<Epetra_FEVector>(map);
159 const bool overlapping =
165 vector = std::make_unique<Epetra_FEVector>(map);
195 const bool omit_zeroing_entries,
196 const bool allow_different_maps)
203 if (allow_different_maps ==
false)
209 const Epetra_MpiComm *my_comm =
210 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
211 const Epetra_MpiComm *v_comm =
212 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
213 const bool same_communicators =
214 my_comm !=
nullptr && v_comm !=
nullptr &&
215 my_comm->DataPtr() == v_comm->DataPtr();
216 if (!same_communicators ||
219 vector = std::make_unique<Epetra_FEVector>(v.
vector->Map());
224 else if (omit_zeroing_entries ==
false)
233 ierr =
vector->PutScalar(0.0);
246 Assert(omit_zeroing_entries ==
false,
248 "It is not possible to exchange data with the "
249 "option 'omit_zeroing_entries' set, which would not write "
255 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
257 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
263 const Epetra_MpiComm *comm_ptr =
264 dynamic_cast<const Epetra_MpiComm *
>(&(v.
vector->Comm()));
289 size_type n_elements = 0, added_elements = 0, block_offset = 0;
291 n_elements += v.
block(block).local_size();
292 std::vector<TrilinosWrappers::types::int_type> global_ids(n_elements, -1);
297 v.
block(block).trilinos_partitioner());
299 global_ids[added_elements++] = glob_elements[i] + block_offset;
302 block_offset += v.
block(block).size();
306 Epetra_Map new_map(v.
size(),
310 v.
block(0).trilinos_partitioner().Comm());
312 auto actual_vec = std::make_unique<Epetra_FEVector>(new_map);
317 v.
block(block).trilinos_vector().ExtractCopy(entries, 0);
318 entries += v.
block(block).local_size();
321 if (import_data ==
true)
324 *actual_vec)) == v.
size(),
329 Epetra_Import data_exchange(
vector->Map(), actual_vec->Map());
331 const int ierr =
vector->Import(*actual_vec, data_exchange, Insert);
337 vector = std::move(actual_vec);
339 const Epetra_MpiComm *comm_ptr =
340 dynamic_cast<const Epetra_MpiComm *
>(&(
vector->Comm()));
355 const bool vector_writable)
359 if (vector_writable ==
false)
361 IndexSet parallel_partitioner = locally_owned_entries;
365 vector = std::make_unique<Epetra_FEVector>(map);
372 ExcMessage(
"A writable vector must not have ghost entries in "
373 "its parallel partitioning"));
375 if (
vector->Map().SameAs(map) ==
false)
376 vector = std::make_unique<Epetra_FEVector>(map);
379 const int ierr =
vector->PutScalar(0.);
384 IndexSet nonlocal_entries(ghost_entries);
388 Epetra_Map nonlocal_map =
391 std::make_unique<Epetra_MultiVector>(nonlocal_map, 1);
413 ExcMessage(
"Vector is not constructed properly."));
417 const Epetra_MpiComm *my_comm =
418 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
419 const Epetra_MpiComm *v_comm =
420 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
421 const bool same_communicators = my_comm !=
nullptr && v_comm !=
nullptr &&
422 my_comm->DataPtr() == v_comm->DataPtr();
449 if (same_communicators && v.
vector->Map().SameAs(
vector->Map()))
462 (v.
vector->Map().UniqueGIDs() ||
vector->Map().UniqueGIDs()))
494 template <
typename number>
520 "Cannot find exchange information!"));
522 ExcMessage(
"The input vector has overlapping data, "
523 "which is not allowed."));
529 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
530 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
545 "Both vectors need to have the same size for import() to work!"));
555 (*this)[idx] = rwv[idx];
560 (*this)[idx] += rwv[idx];
588 "compress() can only be called with VectorOperation add, insert, or unknown"));
598 "The last operation on the Vector and the given last action in the compress() call do not agree!"));
606 const double double_mode = mode;
607 const Epetra_MpiComm *comm_ptr =
615 "Not all processors agree whether the last operation on "
616 "this vector was an addition or a set operation. This will "
617 "prevent the compress() operation from succeeding."));
624 ierr =
vector->GlobalAssemble(mode);
650 if (trilinos_i == -1)
656 vector->Map().MaxMyGID()));
659 value = (*vector)[0][trilinos_i];
669 if (allow_different_maps ==
false)
677# if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0)
678 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
680 vector->Import(*v.
vector, data_exchange, Epetra_AddLocalAlso);
687 Epetra_MultiVector dummy(
vector->Map(), 1,
false);
688 Epetra_Import data_exchange(dummy.Map(), v.
vector->Map());
690 int ierr = dummy.Import(*v.
vector, data_exchange, Insert);
693 ierr =
vector->Update(1.0, dummy, 1.0);
710 if ((*(v.
vector))[0][i] != (*vector)[0][i])
723 return (!(*
this == v));
735 unsigned int flag = 0;
748 const Epetra_MpiComm *mpi_comm =
749 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
752 return num_nonzero == 0;
764 unsigned int flag = 0;
777 const auto max_n_negative =
779 return max_n_negative == 0;
786 const unsigned int precision,
787 const bool scientific,
788 const bool across)
const
791 boost::io::ios_flags_saver restore_flags(out);
794 out.precision(precision);
796 out.setf(std::ios::scientific, std::ios::floatfield);
798 out.setf(std::ios::fixed, std::ios::floatfield);
802 auto global_id = [&](
const size_type index) {
805 out <<
"size:" <<
size() <<
" local_size:" <<
local_size() <<
" :"
808 out <<
"[" << global_id(i) <<
"]: " << (*(
vector))[0][i]
814 int leading_dimension;
815 int ierr =
vector->ExtractView(&val, &leading_dimension);
820 out <<
static_cast<double>(val[i]) <<
' ';
823 out <<
static_cast<double>(val[i]) << std::endl;
853 return sizeof(*this) +
860# include "trilinos_vector.inst"
size_type n_elements() const
void set_size(const size_type size)
void subtract_set(const IndexSet &other)
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
const Epetra_CrsMatrix & trilinos_matrix() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcGhostsPresent()
unsigned int n_blocks() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
BlockType & block(const unsigned int i)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
const MPI_Comm & get_mpi_communicator() const
std::unique_ptr< Epetra_MultiVector > nonlocal_vector
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
size_type local_size() const
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
IndexSet locally_owned_elements() const
void swap(BlockVector &u, BlockVector &v)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
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)
bool operator!=(const Vector &v) const
bool operator==(const Vector &v) const
Epetra_CombineMode last_action
static ::ExceptionBase & ExcTrilinosError(int arg1)
Vector & operator=(const TrilinosScalar s)
void import(const LinearAlgebra::ReadWriteVector< double > &rwv, const VectorOperation::values operation)
bool is_non_negative() const
std::size_t memory_consumption() const
const IndexSet & get_stored_elements() const
TrilinosWrappers::types::int_type global_length(const Epetra_MultiVector &vector)
TrilinosWrappers::types::int_type * my_global_elements(const Epetra_BlockMap &map)
TrilinosWrappers::types::int64_type n_global_elements(const Epetra_BlockMap &map)
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
T max(const T &t, const MPI_Comm &mpi_communicator)