18#ifdef DEAL_II_WITH_TRILINOS
28# include <boost/io/ios_state.hpp>
30# include <Epetra_Export.h>
31# include <Epetra_Import.h>
32# include <Epetra_Vector.h>
54 vector.vector->Map().LID(
60 vector.vector->Map().MinMyGID(),
61 vector.vector->Map().MaxMyGID()));
64 return (*(vector.vector))[0][local_index];
76 , vector(new Epetra_FEVector(
77 Epetra_Map(0, 0, 0,
Utilities::Trilinos::comm_self())))
86 reinit(parallel_partitioning, communicator);
123 vector = std::make_unique<Epetra_FEVector>(
135 reinit(local, ghost, communicator,
false);
145 Epetra_Map map(0, 0, Epetra_MpiComm(MPI_COMM_SELF));
148 vector = std::make_unique<Epetra_FEVector>(map);
161 const bool overlapping =
167 vector = std::make_unique<Epetra_FEVector>(map);
197 const bool omit_zeroing_entries,
198 const bool allow_different_maps)
205 if (allow_different_maps ==
false)
211 const Epetra_MpiComm *my_comm =
212 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
213 const Epetra_MpiComm *v_comm =
214 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
215 const bool same_communicators =
216 my_comm !=
nullptr && v_comm !=
nullptr &&
217 my_comm->DataPtr() == v_comm->DataPtr();
218 if (!same_communicators ||
221 vector = std::make_unique<Epetra_FEVector>(v.
vector->Map());
226 else if (omit_zeroing_entries ==
false)
235 ierr =
vector->PutScalar(0.0);
248 Assert(omit_zeroing_entries ==
false,
250 "It is not possible to exchange data with the "
251 "option 'omit_zeroing_entries' set, which would not write "
257 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
259 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
265 const Epetra_MpiComm *comm_ptr =
266 dynamic_cast<const Epetra_MpiComm *
>(&(v.
vector->Comm()));
291 size_type n_elements = 0, added_elements = 0, block_offset = 0;
293 n_elements += v.
block(block).local_size();
294 std::vector<TrilinosWrappers::types::int_type> global_ids(n_elements, -1);
299 v.
block(block).trilinos_partitioner());
301 global_ids[added_elements++] = glob_elements[i] + block_offset;
304 block_offset += v.
block(block).size();
308 Epetra_Map new_map(v.
size(),
312 v.
block(0).trilinos_partitioner().Comm());
314 auto actual_vec = std::make_unique<Epetra_FEVector>(new_map);
319 v.
block(block).trilinos_vector().ExtractCopy(entries, 0);
320 entries += v.
block(block).local_size();
323 if (import_data ==
true)
326 *actual_vec)) == v.
size(),
331 Epetra_Import data_exchange(
vector->Map(), actual_vec->Map());
333 const int ierr =
vector->Import(*actual_vec, data_exchange, Insert);
339 vector = std::move(actual_vec);
341 const Epetra_MpiComm *comm_ptr =
342 dynamic_cast<const Epetra_MpiComm *
>(&(
vector->Comm()));
357 const bool vector_writable)
361 if (vector_writable ==
false)
363 IndexSet parallel_partitioner = locally_owned_entries;
367 vector = std::make_unique<Epetra_FEVector>(map);
374 ExcMessage(
"A writable vector must not have ghost entries in "
375 "its parallel partitioning"));
377 if (
vector->Map().SameAs(map) ==
false)
378 vector = std::make_unique<Epetra_FEVector>(map);
381 const int ierr =
vector->PutScalar(0.);
386 IndexSet nonlocal_entries(ghost_entries);
390 Epetra_Map nonlocal_map =
393 std::make_unique<Epetra_MultiVector>(nonlocal_map, 1);
413 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
414 const bool make_ghosted,
415 const bool vector_writable)
419 Assert(partitioner->ghost_indices_initialized(),
420 ExcMessage(
"You asked to create a ghosted vector, but the "
421 "partitioner does not provide ghost indices."));
423 this->
reinit(partitioner->locally_owned_range(),
424 partitioner->ghost_indices(),
425 partitioner->get_mpi_communicator(),
430 this->
reinit(partitioner->locally_owned_range(),
431 partitioner->get_mpi_communicator());
441 ExcMessage(
"Vector is not constructed properly."));
445 const Epetra_MpiComm *my_comm =
446 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
447 const Epetra_MpiComm *v_comm =
448 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
449 const bool same_communicators = my_comm !=
nullptr && v_comm !=
nullptr &&
450 my_comm->DataPtr() == v_comm->DataPtr();
477 if (same_communicators && v.
vector->Map().SameAs(
vector->Map()))
490 (v.
vector->Map().UniqueGIDs() ||
vector->Map().UniqueGIDs()))
523 template <
typename number>
549 "Cannot find exchange information!"));
551 ExcMessage(
"The input vector has overlapping data, "
552 "which is not allowed."));
558 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
559 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
574 "Both vectors need to have the same size for import_elements() to work!"));
584 (*
this)[idx] = rwv[idx];
589 (*
this)[idx] += rwv[idx];
617 "compress() can only be called with VectorOperation add, insert, or unknown"));
627 "The last operation on the Vector and the given last action in the compress() call do not agree!"));
635 const double double_mode = mode;
636 const Epetra_MpiComm *comm_ptr =
644 "Not all processors agree whether the last operation on "
645 "this vector was an addition or a set operation. This will "
646 "prevent the compress() operation from succeeding."));
653 ierr =
vector->GlobalAssemble(mode);
679 if (trilinos_i == -1)
685 vector->Map().MaxMyGID()));
688 value = (*vector)[0][trilinos_i];
698 if (allow_different_maps ==
false)
706# if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0)
707 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
709 vector->Import(*v.
vector, data_exchange, Epetra_AddLocalAlso);
716 Epetra_MultiVector dummy(
vector->Map(), 1,
false);
717 Epetra_Import data_exchange(dummy.Map(), v.
vector->Map());
719 int ierr = dummy.Import(*v.
vector, data_exchange, Insert);
722 ierr =
vector->Update(1.0, dummy, 1.0);
739 if ((*(v.
vector))[0][i] != (*vector)[0][i])
752 return (!(*
this == v));
764 unsigned int flag = 0;
777 const Epetra_MpiComm *mpi_comm =
778 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
781 return num_nonzero == 0;
793 unsigned int flag = 0;
806 const auto max_n_negative =
808 return max_n_negative == 0;
815 const unsigned int precision,
816 const bool scientific,
817 const bool across)
const
820 boost::io::ios_flags_saver restore_flags(out);
823 out.precision(precision);
825 out.setf(std::ios::scientific, std::ios::floatfield);
827 out.setf(std::ios::fixed, std::ios::floatfield);
831 auto global_id = [&](
const size_type index) {
834 out <<
"size:" <<
size() <<
" local_size:" <<
local_size() <<
" :"
837 out <<
"[" << global_id(i) <<
"]: " << (*(
vector))[0][i]
843 int leading_dimension;
844 int ierr =
vector->ExtractView(&val, &leading_dimension);
849 out <<
static_cast<double>(val[i]) <<
' ';
852 out <<
static_cast<double>(val[i]) << std::endl;
882 return sizeof(*this) +
889# include "trilinos_vector.inst"
unsigned int n_blocks() const
BlockType & block(const unsigned int i)
bool is_ascending_and_one_to_one(const MPI_Comm communicator) const
size_type n_elements() const
void set_size(const size_type size)
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
void subtract_set(const IndexSet &other)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
const IndexSet & get_stored_elements() const
void compress(VectorOperation::values operation)
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
std::unique_ptr< Epetra_MultiVector > nonlocal_vector
void import_elements(const LinearAlgebra::ReadWriteVector< double > &rwv, const VectorOperation::values operation)
MPI_Comm get_mpi_communicator() const
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 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 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
bool operator!=(const Vector &v) const
bool operator==(const Vector &v) const
Epetra_CombineMode last_action
Vector & operator=(const TrilinosScalar s)
bool is_non_negative() const
std::size_t memory_consumption() const
const Epetra_CrsMatrix & trilinos_matrix() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcGhostsPresent()
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void swap(BlockVector &u, BlockVector &v)
TrilinosWrappers::types::int_type global_length(const Epetra_MultiVector &vector)
int gid(const Epetra_BlockMap &map, int i)
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)
T max(const T &t, const MPI_Comm mpi_communicator)
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm mpi_communicator)