16 #include <deal.II/base/std_cxx14/memory.h> 18 #include <deal.II/lac/trilinos_epetra_vector.h> 20 #ifdef DEAL_II_WITH_TRILINOS 22 # ifdef DEAL_II_WITH_MPI 24 # include <deal.II/base/index_set.h> 26 # include <deal.II/lac/read_write_vector.h> 28 # include <boost/io/ios_state.hpp> 30 # include <Epetra_Import.h> 31 # include <Epetra_Map.h> 32 # include <Epetra_MpiComm.h> 37 DEAL_II_NAMESPACE_OPEN
41 namespace EpetraWrappers
45 , vector(new Epetra_FEVector(
46 Epetra_Map(0, 0, 0,
Utilities::Trilinos::comm_self())))
53 , vector(new Epetra_FEVector(V.trilinos_vector()))
59 const MPI_Comm &communicator)
61 , vector(new Epetra_FEVector(
62 parallel_partitioner.make_trilinos_map(communicator, false)))
69 const MPI_Comm &communicator,
70 const bool omit_zeroing_entries)
72 Epetra_Map input_map =
74 if (
vector->Map().SameAs(input_map) ==
false)
75 vector = std_cxx14::make_unique<Epetra_FEVector>(input_map);
76 else if (omit_zeroing_entries ==
false)
78 const int ierr =
vector->PutScalar(0.);
88 const bool omit_zeroing_entries)
91 Assert(dynamic_cast<const Vector *>(&V) !=
nullptr,
99 omit_zeroing_entries);
117 Epetra_Import data_exchange(
vector->Map(),
140 const int ierr =
vector->PutScalar(s);
153 std::shared_ptr<const CommunicationPatternBase> communication_pattern)
157 if (communication_pattern ==
nullptr)
168 dynamic_cast<const Epetra_MpiComm &
>(
vector->Comm()).Comm());
175 communication_pattern);
179 std::string(
"The communication pattern is not of type ") +
180 "LinearAlgebra::EpetraWrappers::CommunicationPattern."));
186 Epetra_FEVector source_vector(
import.TargetMap());
187 double * values = source_vector.Values();
188 std::copy(V.
begin(), V.
end(), values);
191 vector->Export(source_vector,
import, Insert);
193 vector->Export(source_vector,
import, Add);
195 vector->Export(source_vector,
import, Epetra_Max);
197 vector->Export(source_vector,
import, Epetra_Min);
220 *
this *= 1. / factor;
231 Assert(dynamic_cast<const Vector *>(&V) !=
nullptr,
248 # if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0) 249 Epetra_Import data_exchange(
vector->Map(),
253 Epetra_AddLocalAlso);
260 Epetra_MultiVector dummy(
vector->Map(), 1,
false);
261 Epetra_Import data_exchange(dummy.Map(),
268 ierr =
vector->Update(1.0, dummy, 1.0);
292 Assert(dynamic_cast<const Vector *>(&V) !=
nullptr,
316 const unsigned local_size(
vector->MyLength());
317 for (
unsigned int i = 0; i < local_size; ++i)
327 Assert(dynamic_cast<const Vector *>(&V) !=
nullptr,
350 Assert(dynamic_cast<const Vector *>(&V) !=
nullptr,
353 Assert(dynamic_cast<const Vector *>(&W) !=
nullptr,
367 const int ierr =
vector->Update(
381 Assert(dynamic_cast<const Vector *>(&V) !=
nullptr,
398 Assert(dynamic_cast<const Vector *>(&scaling_factors) !=
nullptr,
402 const Vector &down_scaling_factors =
403 dynamic_cast<const Vector &
>(scaling_factors);
407 const int ierr =
vector->Multiply(1.0,
421 Assert(dynamic_cast<const Vector *>(&V) !=
nullptr,
428 this->
sadd(0., a, V);
445 double * start_ptr = (*vector)[0];
446 const double *ptr = start_ptr, *eptr = start_ptr +
vector->MyLength();
447 unsigned int flag = 0;
459 const Epetra_MpiComm *mpi_comm =
460 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
464 return num_nonzero == 0;
487 int ierr =
vector->Norm1(&norm);
500 int ierr =
vector->Norm2(&norm);
513 int ierr =
vector->NormInf(&norm);
537 # ifndef DEAL_II_WITH_64BIT_INDICES 538 return vector->GlobalLength();
540 return vector->GlobalLength64();
549 const Epetra_MpiComm *epetra_comm =
550 dynamic_cast<const Epetra_MpiComm *
>(&(
vector->Comm()));
552 return epetra_comm->GetMpiComm();
563 if (
vector->Map().LinearMap())
565 # ifndef DEAL_II_WITH_64BIT_INDICES 569 vector->Map().MaxMyGID64() + 1);
572 else if (
vector->Map().NumMyElements() > 0)
574 const size_type n_indices =
vector->Map().NumMyElements();
575 # ifndef DEAL_II_WITH_64BIT_INDICES 576 unsigned int *vector_indices =
577 reinterpret_cast<unsigned int *
>(
vector->Map().MyGlobalElements());
579 size_type *vector_indices =
580 reinterpret_cast<size_type *
>(
vector->Map().MyGlobalElements64());
582 is.
add_indices(vector_indices, vector_indices + n_indices);
591 const Epetra_FEVector &
609 const unsigned int precision,
610 const bool scientific,
611 const bool across)
const 614 boost::io::ios_flags_saver restore_flags(out);
619 int leading_dimension;
620 int ierr =
vector->ExtractView(&val, &leading_dimension);
624 out.precision(precision);
626 out.setf(std::ios::scientific, std::ios::floatfield);
628 out.setf(std::ios::fixed, std::ios::floatfield);
631 for (
int i = 0; i <
vector->MyLength(); ++i)
632 out << val[i] <<
' ';
634 for (
int i = 0; i <
vector->MyLength(); ++i)
635 out << val[i] << std::endl;
648 return sizeof(*this) +
650 (
sizeof(double) +
sizeof(TrilinosWrappers::types::int_type));
657 const MPI_Comm &mpi_comm)
668 DEAL_II_NAMESPACE_CLOSE
virtual double l1_norm() const override
static ::ExceptionBase & ExcDifferentParallelPartitioning()
virtual bool all_zero() const override
static ::ExceptionBase & ExcIO()
virtual ::IndexSet locally_owned_elements() const override
virtual double l2_norm() const override
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
virtual size_type size() const override
virtual void scale(const VectorSpaceVector< double > &scaling_factors) override
#define AssertThrow(cond, exc)
virtual std::size_t memory_consumption() const override
std::unique_ptr< Epetra_FEVector > vector
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
virtual void equ(const double a, const VectorSpaceVector< double > &V) override
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm &mpi_communicator)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void create_epetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm &mpi_comm)
virtual Vector & operator*=(const double factor) override
virtual void import(const ReadWriteVector< double > &V, VectorOperation::values operation, std::shared_ptr< const CommunicationPatternBase > communication_pattern=std::shared_ptr< const CommunicationPatternBase >()) override
const Epetra_FEVector & trilinos_vector() const
Vector & operator=(const Vector &V)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcVectorTypeNotCompatible()
void reinit(const IndexSet ¶llel_partitioner, const MPI_Comm &communicator, const bool omit_zeroing_entries=false)
std::shared_ptr< const CommunicationPattern > epetra_comm_pattern
virtual void sadd(const double s, const double a, const VectorSpaceVector< double > &V) override
void add_range(const size_type begin, const size_type end)
virtual Vector & operator+=(const VectorSpaceVector< double > &V) override
virtual Vector & operator/=(const double factor) override
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const override
virtual double mean_value() const override
virtual Vector & operator-=(const VectorSpaceVector< double > &V) override
MPI_Comm get_mpi_communicator() const
static ::ExceptionBase & ExcNotImplemented()
virtual double operator*(const VectorSpaceVector< double > &V) const override
virtual void add(const double a) override
static ::ExceptionBase & ExcZero()
#define AssertIsFinite(number)
virtual double linfty_norm() const override
const IndexSet & get_stored_elements() const
static ::ExceptionBase & ExcInternalError()
virtual double add_and_dot(const double a, const VectorSpaceVector< double > &V, const VectorSpaceVector< double > &W) override
::IndexSet source_stored_elements