18#ifdef DEAL_II_WITH_TRILINOS
20# ifdef DEAL_II_WITH_MPI
26# include <boost/io/ios_state.hpp>
28# include <Epetra_Import.h>
29# include <Epetra_Map.h>
30# include <Epetra_MpiComm.h>
39 namespace EpetraWrappers
43 , vector(new Epetra_FEVector(
51 , vector(new Epetra_FEVector(
V.trilinos_vector()))
59 , vector(new Epetra_FEVector(
60 parallel_partitioner.make_trilinos_map(communicator, false)))
68 const bool omit_zeroing_entries)
70 Epetra_Map input_map =
72 if (
vector->Map().SameAs(input_map) ==
false)
73 vector = std::make_unique<Epetra_FEVector>(input_map);
74 else if (omit_zeroing_entries ==
false)
76 const int ierr =
vector->PutScalar(0.);
86 const bool omit_zeroing_entries)
97 omit_zeroing_entries);
109 if (
vector->Map().SameAs(
V.trilinos_vector().Map()))
113 if (
size() ==
V.size())
115 Epetra_Import data_exchange(
vector->Map(),
116 V.trilinos_vector().Map());
119 vector->Import(
V.trilinos_vector(), data_exchange, Insert);
124 vector = std::make_unique<Epetra_FEVector>(
V.trilinos_vector());
137 const int ierr =
vector->PutScalar(s);
150 std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
151 communication_pattern)
155 if (communication_pattern ==
nullptr)
161 V.get_stored_elements().size()) ||
165 V.get_stored_elements(),
166 dynamic_cast<const Epetra_MpiComm &
>(
vector->Comm()).Comm());
172 std::dynamic_pointer_cast<const CommunicationPattern>(
173 communication_pattern);
177 std::string(
"The communication pattern is not of type ") +
178 "LinearAlgebra::EpetraWrappers::CommunicationPattern."));
184 Epetra_FEVector source_vector(
import.TargetMap());
185 double *
values = source_vector.Values();
189 vector->Export(source_vector,
import, Insert);
191 vector->Export(source_vector,
import, Add);
193 vector->Export(source_vector,
import, Epetra_Max);
195 vector->Export(source_vector,
import, Epetra_Min);
218 *
this *= 1. / factor;
246# if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0)
247 Epetra_Import data_exchange(
vector->Map(),
251 Epetra_AddLocalAlso);
258 Epetra_MultiVector dummy(
vector->Map(), 1,
false);
259 Epetra_Import data_exchange(dummy.Map(),
266 ierr =
vector->Update(1.0, dummy, 1.0);
314 const unsigned local_size(
vector->MyLength());
315 for (
unsigned int i = 0; i < local_size; ++i)
365 const int ierr =
vector->Update(
396 Assert(
dynamic_cast<const Vector *
>(&scaling_factors) !=
nullptr,
400 const Vector &down_scaling_factors =
401 dynamic_cast<const Vector &
>(scaling_factors);
405 const int ierr =
vector->Multiply(1.0,
426 this->
sadd(0., a,
V);
443 double * start_ptr = (*vector)[0];
444 const double *ptr = start_ptr, *eptr = start_ptr +
vector->MyLength();
445 unsigned int flag = 0;
457 const Epetra_MpiComm *mpi_comm =
458 dynamic_cast<const Epetra_MpiComm *
>(&
vector->Map().Comm());
462 return num_nonzero == 0;
535# ifndef DEAL_II_WITH_64BIT_INDICES
536 return vector->GlobalLength();
538 return vector->GlobalLength64();
547 return vector->MyLength();
555 const Epetra_MpiComm *epetra_comm =
556 dynamic_cast<const Epetra_MpiComm *
>(&(
vector->Comm()));
558 return epetra_comm->GetMpiComm();
569 if (
vector->Map().LinearMap())
571# ifndef DEAL_II_WITH_64BIT_INDICES
575 vector->Map().MaxMyGID64() + 1);
578 else if (
vector->Map().NumMyElements() > 0)
581# ifndef DEAL_II_WITH_64BIT_INDICES
582 unsigned int *vector_indices =
583 reinterpret_cast<unsigned int *
>(
vector->Map().MyGlobalElements());
588 is.
add_indices(vector_indices, vector_indices + n_indices);
597 const Epetra_FEVector &
615 const unsigned int precision,
616 const bool scientific,
617 const bool across)
const
620 boost::io::ios_flags_saver restore_flags(out);
625 int leading_dimension;
626 int ierr =
vector->ExtractView(&val, &leading_dimension);
630 out.precision(precision);
632 out.setf(std::ios::scientific, std::ios::floatfield);
634 out.setf(std::ios::fixed, std::ios::floatfield);
637 for (
int i = 0; i <
vector->MyLength(); ++i)
638 out << val[i] <<
' ';
640 for (
int i = 0; i <
vector->MyLength(); ++i)
641 out << val[i] << std::endl;
654 return sizeof(*this) +
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
void add_range(const size_type begin, const size_type end)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
virtual double mean_value() const override
virtual Vector & operator/=(const double factor) override
virtual double l2_norm() const override
Vector & operator=(const Vector &V)
MPI_Comm get_mpi_communicator() const
std::unique_ptr< Epetra_FEVector > vector
const Epetra_FEVector & trilinos_vector() const
virtual double l1_norm() const override
virtual ::IndexSet locally_owned_elements() const override
virtual double add_and_dot(const double a, const VectorSpaceVector< double > &V, const VectorSpaceVector< double > &W) override
virtual void scale(const VectorSpaceVector< double > &scaling_factors) override
void create_epetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm &mpi_comm)
virtual double operator*(const VectorSpaceVector< double > &V) const override
void reinit(const IndexSet ¶llel_partitioner, const MPI_Comm &communicator, const bool omit_zeroing_entries=false)
virtual size_type size() const override
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const override
virtual Vector & operator-=(const VectorSpaceVector< double > &V) override
size_type locally_owned_size() const
virtual double linfty_norm() const override
virtual bool all_zero() const override
virtual void import(const ReadWriteVector< double > &V, VectorOperation::values operation, std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > communication_pattern={}) override
virtual void add(const double a) override
::IndexSet source_stored_elements
virtual void equ(const double a, const VectorSpaceVector< double > &V) override
virtual std::size_t memory_consumption() const override
std::shared_ptr< const CommunicationPattern > epetra_comm_pattern
virtual Vector & operator*=(const double factor) override
virtual Vector & operator+=(const VectorSpaceVector< double > &V) override
virtual void sadd(const double s, const double a, const VectorSpaceVector< double > &V) override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcVectorTypeNotCompatible()
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcZero()
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIsFinite(number)
static ::ExceptionBase & ExcDifferentParallelPartitioning()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
T sum(const T &t, const MPI_Comm &mpi_communicator)
const Epetra_Comm & comm_self()
void copy(const T *begin, const T *end, U *dest)