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(
85 reinit(parallel_partitioning, communicator);
121 vector = std::make_unique<Epetra_FEVector>(
133 reinit(local, ghost, communicator,
false);
143# ifdef DEAL_II_WITH_MPI
144 Epetra_Map map(0, 0, Epetra_MpiComm(MPI_COMM_SELF));
146 Epetra_Map map(0, 0, Epetra_SerialComm());
150 vector = std::make_unique<Epetra_FEVector>(map);
163 const bool overlapping =
169 vector = std::make_unique<Epetra_FEVector>(map);
199 const bool omit_zeroing_entries,
200 const bool allow_different_maps)
207 if (allow_different_maps ==
false)
213# ifdef DEAL_II_WITH_MPI
214 const Epetra_MpiComm *my_comm =
215 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
216 const Epetra_MpiComm *v_comm =
217 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
218 const bool same_communicators =
219 my_comm !=
nullptr && v_comm !=
nullptr &&
220 my_comm->DataPtr() == v_comm->DataPtr();
222 const bool same_communicators =
true;
224 if (!same_communicators ||
227 vector = std::make_unique<Epetra_FEVector>(v.
vector->Map());
232 else if (omit_zeroing_entries ==
false)
241 ierr =
vector->PutScalar(0.0);
254 Assert(omit_zeroing_entries ==
false,
256 "It is not possible to exchange data with the "
257 "option 'omit_zeroing_entries' set, which would not write "
263 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
265 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
270# if defined(DEBUG) && defined(DEAL_II_WITH_MPI)
271 const Epetra_MpiComm *comm_ptr =
272 dynamic_cast<const Epetra_MpiComm *
>(&(v.
vector->Comm()));
297 size_type n_elements = 0, added_elements = 0, block_offset = 0;
299 n_elements += v.
block(block).local_size();
300 std::vector<TrilinosWrappers::types::int_type> global_ids(n_elements, -1);
305 v.
block(block).trilinos_partitioner());
307 global_ids[added_elements++] = glob_elements[i] + block_offset;
310 block_offset += v.
block(block).size();
314 Epetra_Map new_map(v.
size(),
318 v.
block(0).trilinos_partitioner().Comm());
320 auto actual_vec = std::make_unique<Epetra_FEVector>(new_map);
325 v.
block(block).trilinos_vector().ExtractCopy(entries, 0);
326 entries += v.
block(block).local_size();
329 if (import_data ==
true)
332 *actual_vec)) == v.
size(),
337 Epetra_Import data_exchange(
vector->Map(), actual_vec->Map());
339 const int ierr =
vector->Import(*actual_vec, data_exchange, Insert);
345 vector = std::move(actual_vec);
346# if defined(DEBUG) && defined(DEAL_II_WITH_MPI)
347 const Epetra_MpiComm *comm_ptr =
348 dynamic_cast<const Epetra_MpiComm *
>(&(
vector->Comm()));
363 const bool vector_writable)
367 if (vector_writable ==
false)
369 IndexSet parallel_partitioner = locally_owned_entries;
373 vector = std::make_unique<Epetra_FEVector>(map);
380 ExcMessage(
"A writable vector must not have ghost entries in "
381 "its parallel partitioning"));
383 if (
vector->Map().SameAs(map) ==
false)
384 vector = std::make_unique<Epetra_FEVector>(map);
387 const int ierr =
vector->PutScalar(0.);
392 IndexSet nonlocal_entries(ghost_entries);
396 Epetra_Map nonlocal_map =
399 std::make_unique<Epetra_MultiVector>(nonlocal_map, 1);
421 ExcMessage(
"Vector is not constructed properly."));
425# ifdef DEAL_II_WITH_MPI
426 const Epetra_MpiComm *my_comm =
427 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
428 const Epetra_MpiComm *v_comm =
429 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
430 const bool same_communicators = my_comm !=
nullptr && v_comm !=
nullptr &&
431 my_comm->DataPtr() == v_comm->DataPtr();
454 const bool same_communicators =
true;
461 if (same_communicators && v.
vector->Map().SameAs(
vector->Map()))
474 (v.
vector->Map().UniqueGIDs() ||
vector->Map().UniqueGIDs()))
506 template <
typename number>
532 "Cannot find exchange information!"));
534 ExcMessage(
"The input vector has overlapping data, "
535 "which is not allowed."));
541 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
542 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
557 "Both vectors need to have the same size for import() to work!"));
567 (*this)[idx] = rwv[idx];
572 (*this)[idx] += rwv[idx];
600 "compress() can only be called with VectorOperation add, insert, or unknown"));
610 "The last operation on the Vector and the given last action in the compress() call do not agree!"));
615# ifdef DEAL_II_WITH_MPI
619 double double_mode = mode;
620 const Epetra_MpiComm *comm_ptr =
627 "Not all processors agree whether the last operation on "
628 "this vector was an addition or a set operation. This will "
629 "prevent the compress() operation from succeeding."));
637 ierr =
vector->GlobalAssemble(mode);
663 if (trilinos_i == -1)
669 vector->Map().MaxMyGID()));
672 value = (*vector)[0][trilinos_i];
682 if (allow_different_maps ==
false)
690# if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0)
691 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
693 vector->Import(*v.
vector, data_exchange, Epetra_AddLocalAlso);
700 Epetra_MultiVector dummy(
vector->Map(), 1,
false);
701 Epetra_Import data_exchange(dummy.Map(), v.
vector->Map());
703 int ierr = dummy.Import(*v.
vector, data_exchange, Insert);
706 ierr =
vector->Update(1.0, dummy, 1.0);
723 if ((*(v.
vector))[0][i] != (*vector)[0][i])
736 return (!(*
this == v));
748 unsigned int flag = 0;
759# ifdef DEAL_II_WITH_MPI
762 const Epetra_MpiComm *mpi_comm =
763 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
766 return num_nonzero == 0;
777# ifdef DEAL_II_WITH_MPI
788 int leading_dimension;
789 int ierr =
vector->ExtractView(&start_ptr, &leading_dimension);
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"
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)
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
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
TrilinosWrappers::types::int_type global_length(const Epetra_MultiVector &vector)
TrilinosWrappers::types::int_type n_global_elements(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type * my_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)
const Epetra_Comm & comm_self()