16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/lac/trilinos_vector_base.h> 19 #ifdef DEAL_II_WITH_TRILINOS 23 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
24 # include <Epetra_Import.h> 25 # include <Epetra_Export.h> 29 DEAL_II_NAMESPACE_OPEN
35 #ifndef DEAL_II_WITH_64BIT_INDICES 39 int global_length(
const Epetra_FEVector &vector)
41 return vector.GlobalLength();
47 long long int global_length(
const Epetra_FEVector &vector)
49 return vector.GlobalLength64();
56 VectorReference::operator TrilinosScalar ()
const 58 Assert (index < vector.size(),
59 ExcIndexRange (index, 0, vector.size()));
65 const TrilinosWrappers::types::int_type local_index =
66 vector.vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(index));
68 VectorBase::ExcAccessToNonLocalElement (index, vector.local_size(),
69 vector.vector->Map().MinMyGID(),
70 vector.vector->Map().MaxMyGID()));
73 return (*(vector.vector))[0][local_index];
84 #ifdef DEAL_II_WITH_MPI
85 vector(new Epetra_FEVector(
86 Epetra_Map(0,0,Epetra_MpiComm(MPI_COMM_SELF))))
88 vector(new Epetra_FEVector(
89 Epetra_Map(0,0,Epetra_SerialComm())))
101 vector(new Epetra_FEVector(*v.vector))
116 #ifdef DEAL_II_WITH_MPI 117 Epetra_Map map (0, 0, Epetra_MpiComm(MPI_COMM_SELF));
119 Epetra_Map map (0, 0, Epetra_SerialComm());
123 vector.reset (
new Epetra_FEVector(map));
132 Assert (vector.get() != 0,
133 ExcMessage(
"Vector is not constructed properly."));
138 vector.reset (
new Epetra_FEVector(*v.
vector));
144 ExcMessage (
"The Epetra maps in the assignment operator =" 145 " do not match, even though the local_range " 146 " seems to be the same. Check vector setup!"));
151 ierr = vector->Update(1.0, *v.
vector, 0.0);
162 template <
typename number>
167 ExcDimensionMismatch(
size(), v.size()));
177 std::pair<size_type, size_type>
179 for (size_type i=local_range.first; i<local_range.second; ++i)
180 (*vector)[0][i-local_range.first] = v(i);
197 if (given_last_action==::VectorOperation::add)
199 else if (given_last_action==::VectorOperation::insert)
204 # ifdef DEAL_II_WITH_MPI 208 double double_mode = mode;
211 dynamic_cast<const Epetra_MpiComm *>
213 Assert(result.max-result.min<1e-5,
214 ExcMessage (
"Not all processors agree whether the last operation on " 215 "this vector was an addition or a set operation. This will " 216 "prevent the compress() operation from succeeding."));
224 ierr = vector->GlobalAssemble(mode);
243 TrilinosWrappers::types::int_type trilinos_i =
244 vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(index));
248 if (trilinos_i == -1)
251 return (*vector)[0][trilinos_i];
260 TrilinosWrappers::types::int_type trilinos_i =
261 vector->Map().LID(static_cast<TrilinosWrappers::types::int_type>(index));
262 TrilinosScalar value = 0.;
266 if (trilinos_i == -1)
269 vector->Map().MinMyGID(),
270 vector->Map().MaxMyGID()));
273 value = (*vector)[0][trilinos_i];
282 const bool allow_different_maps)
284 if (allow_different_maps ==
false)
290 ExcDimensionMismatch (
size(), v.
size()));
292 #if DEAL_II_TRILINOS_VERSION_GTE(11,11,0) 293 Epetra_Import data_exchange (vector->Map(), v.
vector->Map());
294 int ierr = vector->Import(*v.
vector, data_exchange, Epetra_AddLocalAlso);
301 Epetra_MultiVector dummy(vector->Map(), 1,
false);
302 Epetra_Import data_exchange (dummy.Map(), v.
vector->Map());
304 int ierr = dummy.Import(*v.
vector, data_exchange, Insert);
307 ierr = vector->Update (1.0, dummy, 1.0);
319 ExcDimensionMismatch(
size(), v.
size()));
325 if ((*(v.
vector))[0][i]!=(*vector)[0][i])
return false;
336 ExcDimensionMismatch(
size(), v.
size()));
338 return (!(*
this==v));
348 TrilinosScalar *start_ptr = (*vector)[0];
349 const TrilinosScalar *ptr = start_ptr,
351 unsigned int flag = 0;
362 #ifdef DEAL_II_WITH_MPI 365 const Epetra_MpiComm *mpi_comm
366 =
dynamic_cast<const Epetra_MpiComm *
>(&vector->Map().Comm());
368 return num_nonzero == 0;
380 #ifdef DEAL_II_WITH_MPI 390 TrilinosScalar *start_ptr;
391 int leading_dimension;
392 int ierr = vector->ExtractView (&start_ptr, &leading_dimension);
399 const TrilinosScalar *ptr = start_ptr,
400 *eptr = start_ptr +
size();
420 const TrilinosScalar b,
433 if (vector->Map().SameAs(v.
vector->Map())==
false)
435 sadd(0., a, v, b, w);
445 ExcDifferentParallelPartitioning());
446 int ierr = vector->Update(a, *v.
vector, b, *w.
vector, 0.0);
462 Assert (global_length(*vector)!=0, ExcEmptyObject());
465 for (size_type j=0; j<
size(); ++j)
467 double t = (*vector)[0][j];
470 std::printf (format, t);
472 std::printf (
" %5.2f",
double(t));
481 const unsigned int precision,
482 const bool scientific,
483 const bool across)
const 494 int leading_dimension;
495 int ierr = vector->ExtractView (&val, &leading_dimension);
498 out.precision (precision);
500 out.setf (std::ios::scientific, std::ios::floatfield);
502 out.setf (std::ios::fixed, std::ios::floatfield);
505 for (size_type i=0; i<
size(); ++i)
506 out << static_cast<double>(val[i]) <<
' ';
508 for (size_type i=0; i<
size(); ++i)
509 out << static_cast<double>(val[i]) << std::endl;
524 std::swap(vector, v.
vector);
539 sizeof(TrilinosWrappers::types::int_type) );
547 #include "trilinos_vector_base.inst" 550 DEAL_II_NAMESPACE_CLOSE
552 #endif // DEAL_II_WITH_TRILINOS
TrilinosScalar el(const size_type index) const DEAL_II_DEPRECATED
bool is_non_negative() const
reference operator()(const size_type index)
void sadd(const TrilinosScalar s, const VectorBase &V)
::ExceptionBase & ExcMessage(std::string arg1)
std_cxx11::shared_ptr< Epetra_MultiVector > nonlocal_vector
std::pair< size_type, size_type > local_range() const
void print(const char *format=0) const DEAL_II_DEPRECATED
#define AssertThrow(cond, exc)
void add(const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values)
T sum(const T &t, const MPI_Comm &mpi_communicator)
std_cxx11::shared_ptr< Epetra_FEVector > vector
#define Assert(cond, exc)
const Epetra_Map & vector_partitioner() const
std::size_t memory_consumption() const
bool has_ghost_elements() const
void compress(::VectorOperation::values operation)
VectorBase & operator=(const TrilinosScalar s)
bool operator!=(const VectorBase &v) const
Epetra_CombineMode last_action
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
bool operator==(const VectorBase &v) const
#define AssertIsFinite(number)
size_type local_size() const
void equ(const TrilinosScalar a, const VectorBase &V)