16 #include <deal.II/lac/trilinos_vector.h> 18 #ifdef DEAL_II_WITH_TRILINOS 20 # include <deal.II/base/mpi.h> 21 # include <deal.II/base/std_cxx14/memory.h> 23 # include <deal.II/lac/read_write_vector.h> 24 # include <deal.II/lac/trilinos_index_access.h> 25 # include <deal.II/lac/trilinos_parallel_block_vector.h> 26 # include <deal.II/lac/trilinos_sparse_matrix.h> 28 # include <boost/io/ios_state.hpp> 30 # include <Epetra_Export.h> 31 # include <Epetra_Import.h> 32 # include <Epetra_Vector.h> 37 DEAL_II_NAMESPACE_OPEN
43 VectorReference::operator TrilinosScalar()
const 51 const TrilinosWrappers::types::int_type local_index =
52 vector.vector->Map().LID(
53 static_cast<TrilinosWrappers::types::int_type>(index));
58 vector.vector->Map().MinMyGID(),
59 vector.vector->Map().MaxMyGID()));
62 return (*(vector.vector))[0][local_index];
73 , vector(new Epetra_FEVector(
74 Epetra_Map(0, 0, 0,
Utilities::Trilinos::comm_self())))
80 const MPI_Comm &communicator)
83 reinit(parallel_partitioning, communicator);
92 vector = std_cxx14::make_unique<Epetra_FEVector>(*v.
vector);
109 const MPI_Comm &communicator)
113 static_cast<size_type
>(
119 vector = std_cxx14::make_unique<Epetra_FEVector>(
128 const MPI_Comm &communicator)
131 reinit(local, ghost, communicator,
false);
141 # ifdef DEAL_II_WITH_MPI 142 Epetra_Map map(0, 0, Epetra_MpiComm(MPI_COMM_SELF));
144 Epetra_Map map(0, 0, Epetra_SerialComm());
148 vector = std_cxx14::make_unique<Epetra_FEVector>(map);
156 const MPI_Comm &communicator,
161 const bool overlapping =
167 vector = std_cxx14::make_unique<Epetra_FEVector>(map);
184 const size_type n_elements_global =
197 const bool omit_zeroing_entries,
198 const bool allow_different_maps)
205 if (allow_different_maps ==
false)
211 # ifdef DEAL_II_WITH_MPI 212 const Epetra_MpiComm *my_comm =
213 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
214 const Epetra_MpiComm *v_comm =
215 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
216 const bool same_communicators =
217 my_comm !=
nullptr && v_comm !=
nullptr &&
218 my_comm->DataPtr() == v_comm->DataPtr();
220 const bool same_communicators =
true;
222 if (!same_communicators ||
225 vector = std_cxx14::make_unique<Epetra_FEVector>(v.
vector->Map());
230 else if (omit_zeroing_entries ==
false)
239 ierr =
vector->PutScalar(0.0);
252 Assert(omit_zeroing_entries ==
false,
254 "It is not possible to exchange data with the " 255 "option 'omit_zeroing_entries' set, which would not write " 261 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
263 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
268 # if defined(DEBUG) && defined(DEAL_II_WITH_MPI) 269 const Epetra_MpiComm *comm_ptr =
270 dynamic_cast<const Epetra_MpiComm *
>(&(v.
vector->Comm()));
272 const size_type n_elements_global =
295 size_type n_elements = 0, added_elements = 0, block_offset = 0;
296 for (size_type block = 0; block < v.
n_blocks(); ++block)
297 n_elements += v.
block(block).local_size();
298 std::vector<TrilinosWrappers::types::int_type> global_ids(n_elements, -1);
299 for (size_type block = 0; block < v.
n_blocks(); ++block)
301 TrilinosWrappers::types::int_type *glob_elements =
303 v.
block(block).trilinos_partitioner());
304 for (size_type i = 0; i < v.
block(block).local_size(); ++i)
305 global_ids[added_elements++] = glob_elements[i] + block_offset;
308 block_offset += v.
block(block).size();
312 Epetra_Map new_map(v.
size(),
316 v.
block(0).trilinos_partitioner().Comm());
318 auto actual_vec = std_cxx14::make_unique<Epetra_FEVector>(new_map);
320 TrilinosScalar *entries = (*actual_vec)[0];
321 for (size_type block = 0; block < v.
n_blocks(); ++block)
323 v.
block(block).trilinos_vector().ExtractCopy(entries, 0);
324 entries += v.
block(block).local_size();
327 if (import_data ==
true)
330 *actual_vec)) == v.
size(),
335 Epetra_Import data_exchange(
vector->Map(), actual_vec->Map());
337 const int ierr =
vector->Import(*actual_vec, data_exchange, Insert);
343 vector = std::move(actual_vec);
344 # if defined(DEBUG) && defined(DEAL_II_WITH_MPI) 345 const Epetra_MpiComm *comm_ptr =
346 dynamic_cast<const Epetra_MpiComm *
>(&(
vector->Comm()));
348 const size_type n_elements_global =
360 const MPI_Comm &communicator,
361 const bool vector_writable)
365 if (vector_writable ==
false)
367 IndexSet parallel_partitioner = locally_owned_entries;
371 vector = std_cxx14::make_unique<Epetra_FEVector>(map);
378 ExcMessage(
"A writable vector must not have ghost entries in " 379 "its parallel partitioning"));
381 if (
vector->Map().SameAs(map) ==
false)
382 vector = std_cxx14::make_unique<Epetra_FEVector>(map);
385 const int ierr =
vector->PutScalar(0.);
390 IndexSet nonlocal_entries(ghost_entries);
394 Epetra_Map nonlocal_map =
397 std_cxx14::make_unique<Epetra_MultiVector>(nonlocal_map, 1);
406 const size_type n_elements_global =
419 ExcMessage(
"Vector is not constructed properly."));
423 # ifdef DEAL_II_WITH_MPI 424 const Epetra_MpiComm *my_comm =
425 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Comm());
426 const Epetra_MpiComm *v_comm =
427 dynamic_cast<const Epetra_MpiComm *
>(&v.
vector->Comm());
428 const bool same_communicators = my_comm !=
nullptr && v_comm !=
nullptr &&
429 my_comm->DataPtr() == v_comm->DataPtr();
452 const bool same_communicators =
true;
459 if (same_communicators && v.
vector->Map().SameAs(
vector->Map()))
472 (v.
vector->Map().UniqueGIDs() ||
vector->Map().UniqueGIDs()))
480 vector = std_cxx14::make_unique<Epetra_FEVector>(*v.
vector);
505 template <
typename number>
531 "Cannot find exchange information!"));
533 ExcMessage(
"The input vector has overlapping data, " 534 "which is not allowed."));
540 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
541 const int ierr =
vector->Import(*v.
vector, data_exchange, Insert);
556 "Both vectors need to have the same size for import() to work!"));
566 (*
this)[idx] = rwv[idx];
571 (*
this)[idx] += rwv[idx];
599 "compress() can only be called with VectorOperation add, insert, or unknown"));
609 "The last operation on the Vector and the given last action in the compress() call do not agree!"));
614 # ifdef DEAL_II_WITH_MPI 618 double double_mode = mode;
619 const Epetra_MpiComm *comm_ptr =
626 "Not all processors agree whether the last operation on " 627 "this vector was an addition or a set operation. This will " 628 "prevent the compress() operation from succeeding."));
636 ierr =
vector->GlobalAssemble(mode);
656 TrilinosWrappers::types::int_type trilinos_i =
vector->Map().LID(
657 static_cast<TrilinosWrappers::types::int_type>(index));
658 TrilinosScalar value = 0.;
662 if (trilinos_i == -1)
668 vector->Map().MaxMyGID()));
671 value = (*vector)[0][trilinos_i];
681 if (allow_different_maps ==
false)
689 # if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0) 690 Epetra_Import data_exchange(
vector->Map(), v.
vector->Map());
692 vector->Import(*v.
vector, data_exchange, Epetra_AddLocalAlso);
699 Epetra_MultiVector dummy(
vector->Map(), 1,
false);
700 Epetra_Import data_exchange(dummy.Map(), v.
vector->Map());
702 int ierr = dummy.Import(*v.
vector, data_exchange, Insert);
705 ierr =
vector->Update(1.0, dummy, 1.0);
722 if ((*(v.
vector))[0][i] != (*vector)[0][i])
735 return (!(*
this == v));
745 TrilinosScalar * start_ptr = (*vector)[0];
746 const TrilinosScalar *ptr = start_ptr, *eptr = start_ptr +
local_size();
747 unsigned int flag = 0;
758 # ifdef DEAL_II_WITH_MPI 761 const Epetra_MpiComm *mpi_comm =
762 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
765 return num_nonzero == 0;
776 # ifdef DEAL_II_WITH_MPI 786 TrilinosScalar *start_ptr;
787 int leading_dimension;
788 int ierr =
vector->ExtractView(&start_ptr, &leading_dimension);
795 const TrilinosScalar *ptr = start_ptr, *eptr = start_ptr +
size();
814 const unsigned int precision,
815 const bool scientific,
816 const bool across)
const 819 boost::io::ios_flags_saver restore_flags(out);
822 out.precision(precision);
824 out.setf(std::ios::scientific, std::ios::floatfield);
826 out.setf(std::ios::fixed, std::ios::floatfield);
830 auto global_id = [&](
const size_type index) {
831 return gid(
vector->Map(), index);
833 out <<
"size:" <<
size() <<
" local_size:" <<
local_size() <<
" :" 836 out <<
"[" << global_id(i) <<
"]: " << (*(
vector))[0][i]
842 int leading_dimension;
843 int ierr =
vector->ExtractView(&val, &leading_dimension);
847 for (size_type i = 0; i <
size(); ++i)
848 out << static_cast<double>(val[i]) <<
' ';
850 for (size_type i = 0; i <
size(); ++i)
851 out << static_cast<double>(val[i]) << std::endl;
881 return sizeof(*this) +
883 (
sizeof(double) +
sizeof(TrilinosWrappers::types::int_type));
892 # include "trilinos_vector.inst" 896 DEAL_II_NAMESPACE_CLOSE
898 #endif // DEAL_II_WITH_TRILINOS Vector & operator=(const TrilinosScalar s)
std::unique_ptr< Epetra_MultiVector > nonlocal_vector
static ::ExceptionBase & ExcIO()
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)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
std::pair< size_type, size_type > local_range() const
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
LinearAlgebra::distributed::Vector< Number > Vector
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)
void swap(BlockVector &u, BlockVector &v)
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
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
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 swap(Vector< Number > &u, Vector< Number > &v)
void set_size(const size_type size)
size_type local_size() const
unsigned int n_blocks() const
void import(const LinearAlgebra::ReadWriteVector< double > &rwv, const VectorOperation::values operation)
const Epetra_BlockMap & trilinos_partitioner() const
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()