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> 53 vector.vector->Map().LID(
54 static_cast<TrilinosWrappers::types::int_type>(index));
59 vector.vector->Map().MinMyGID(),
60 vector.vector->Map().MaxMyGID()));
63 return (*(vector.vector))[0][local_index];
75 , vector(new Epetra_FEVector(
82 const MPI_Comm &communicator)
85 reinit(parallel_partitioning, communicator);
94 vector = std_cxx14::make_unique<Epetra_FEVector>(*v.
vector);
111 const MPI_Comm &communicator)
121 vector = std_cxx14::make_unique<Epetra_FEVector>(
130 const MPI_Comm &communicator)
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_cxx14::make_unique<Epetra_FEVector>(map);
158 const MPI_Comm &communicator,
163 const bool overlapping =
169 vector = std_cxx14::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_cxx14::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_cxx14::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()));
362 const MPI_Comm &communicator,
363 const bool vector_writable)
367 if (vector_writable ==
false)
369 IndexSet parallel_partitioner = locally_owned_entries;
373 vector = std_cxx14::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_cxx14::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_cxx14::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()))
482 vector = std_cxx14::make_unique<Epetra_FEVector>(*v.
vector);
507 template <
typename number>
519 for (
size_type i = local_range.first; i < local_range.second; ++i)
520 (*vector)[0][i - local_range.first] = v(i);
533 "Cannot find exchange information!"));
535 ExcMessage(
"The input vector has overlapping data, " 536 "which is not allowed."));
542 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
543 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
558 "Both vectors need to have the same size for import() to work!"));
568 (*
this)[idx] = rwv[idx];
573 (*
this)[idx] += rwv[idx];
601 "compress() can only be called with VectorOperation add, insert, or unknown"));
611 "The last operation on the Vector and the given last action in the compress() call do not agree!"));
616 # ifdef DEAL_II_WITH_MPI 620 double double_mode = mode;
621 const Epetra_MpiComm *comm_ptr =
628 "Not all processors agree whether the last operation on " 629 "this vector was an addition or a set operation. This will " 630 "prevent the compress() operation from succeeding."));
638 ierr =
vector->GlobalAssemble(mode);
659 static_cast<TrilinosWrappers::types::int_type>(index));
664 if (trilinos_i == -1)
670 vector->Map().MaxMyGID()));
673 value = (*vector)[0][trilinos_i];
683 if (allow_different_maps ==
false)
691 # if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0) 692 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
694 vector->Import(*v.
vector, data_exchange, Epetra_AddLocalAlso);
702 Epetra_Import data_exchange(
dummy.Map(), v.
vector->Map());
704 int ierr =
dummy.Import(*v.
vector, data_exchange, Insert);
724 if ((*(v.
vector))[0][i] != (*vector)[0][i])
737 return (!(*
this == v));
749 unsigned int flag = 0;
760 # ifdef DEAL_II_WITH_MPI 763 const Epetra_MpiComm *mpi_comm =
764 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
767 return num_nonzero == 0;
778 # ifdef DEAL_II_WITH_MPI 789 int leading_dimension;
790 int ierr =
vector->ExtractView(&start_ptr, &leading_dimension);
816 const unsigned int precision,
817 const bool scientific,
818 const bool across)
const 821 boost::io::ios_flags_saver restore_flags(out);
824 out.precision(precision);
826 out.setf(std::ios::scientific, std::ios::floatfield);
828 out.setf(std::ios::fixed, std::ios::floatfield);
832 auto global_id = [&](
const size_type index) {
835 out <<
"size:" <<
size() <<
" local_size:" <<
local_size() <<
" :" 838 out <<
"[" << global_id(i) <<
"]: " << (*(
vector))[0][i]
844 int leading_dimension;
845 int ierr =
vector->ExtractView(&val, &leading_dimension);
850 out << static_cast<double>(val[i]) <<
' ';
853 out << static_cast<double>(val[i]) << std::endl;
883 return sizeof(*this) +
890 # include "trilinos_vector.inst" 897 #endif // DEAL_II_WITH_TRILINOS
Vector & operator=(const TrilinosScalar s)
std::unique_ptr< Epetra_MultiVector > nonlocal_vector
int gid(const Epetra_BlockMap &map, int i)
static ::ExceptionBase & ExcIO()
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
TrilinosWrappers::types::int_type global_length(const Epetra_MultiVector &vector)
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
#define AssertIndexRange(index, range)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
std::pair< size_type, size_type > local_range() const
#define AssertThrow(cond, exc)
const Epetra_Comm & comm_self()
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
bool operator==(const Vector &v) const
static ::ExceptionBase & ExcMessage(std::string arg1)
reference operator()(const size_type index)
T sum(const T &t, const MPI_Comm &mpi_communicator)
TrilinosWrappers::types::int_type * my_global_elements(const Epetra_BlockMap &map)
void subtract_set(const IndexSet &other)
#define Assert(cond, exc)
bool is_non_negative() const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::unique_ptr< Epetra_FEVector > vector
void compress(::VectorOperation::values operation)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
#define DEAL_II_NAMESPACE_CLOSE
std::size_t memory_consumption() const
const Epetra_CrsMatrix & trilinos_matrix() const
void import_nonlocal_data_for_fe(const ::TrilinosWrappers::SparseMatrix &matrix, const Vector &vector)
bool has_ghost_elements() const
const types::global_dof_index * dummy()
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
IndexSet locally_owned_elements() const
static ::ExceptionBase & ExcGhostsPresent()
TrilinosWrappers::types::int_type n_global_elements(const Epetra_BlockMap &map)
void set_size(const size_type size)
size_type local_size() const
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
unsigned int n_blocks() const
void import(const LinearAlgebra::ReadWriteVector< double > &rwv, const VectorOperation::values operation)
const Epetra_BlockMap & trilinos_partitioner() const
#define DEAL_II_NAMESPACE_OPEN
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
static ::ExceptionBase & ExcNotImplemented()
Epetra_CombineMode last_action
bool operator!=(const Vector &v) const
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcTrilinosError(int arg1)
BlockType & block(const unsigned int i)
size_type n_elements() const
const IndexSet & get_stored_elements() const
static ::ExceptionBase & ExcInternalError()