16 #ifndef dealii__la_parallel_vector_h 17 #define dealii__la_parallel_vector_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/mpi.h> 21 #include <deal.II/base/numbers.h> 22 #include <deal.II/base/partitioner.h> 23 #include <deal.II/base/std_cxx11/shared_ptr.h> 24 #include <deal.II/base/thread_management.h> 25 #include <deal.II/lac/vector_space_vector.h> 26 #include <deal.II/lac/vector_type_traits.h> 30 DEAL_II_NAMESPACE_OPEN
45 #ifdef DEAL_II_WITH_PETSC 55 #ifdef DEAL_II_WITH_TRILINOS 175 template <
typename Number>
180 typedef Number value_type;
181 typedef value_type *pointer;
182 typedef const value_type *const_pointer;
183 typedef value_type *iterator;
184 typedef const value_type *const_iterator;
185 typedef value_type &reference;
186 typedef const value_type &const_reference;
243 const MPI_Comm communicator);
249 const MPI_Comm communicator);
257 Vector (
const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &
partitioner);
269 const bool omit_zeroing_entries =
false);
281 template <
typename Number2>
283 const bool omit_zeroing_entries =
false);
303 const MPI_Comm communicator);
309 const MPI_Comm communicator);
317 void reinit (
const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &
partitioner);
361 template <
typename Number2>
365 #ifdef DEAL_II_WITH_PETSC 380 #ifdef DEAL_II_WITH_TRILINOS 462 void compress_start (
const unsigned int communication_channel = 0,
567 std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
568 std_cxx11::shared_ptr<const CommunicationPatternBase> ());
578 virtual void add(
const Number a);
595 virtual void add (
const std::vector<size_type> &indices,
596 const std::vector<Number> &values);
602 virtual void sadd(
const Number s,
const Number a,
621 virtual real_type
l1_norm()
const;
627 virtual real_type
l2_norm()
const;
659 virtual size_type
size()
const;
677 virtual void print(std::ostream &out,
678 const unsigned int precision=3,
679 const bool scientific=
true,
680 const bool across=
true)
const;
704 template <
typename OtherNumber>
705 void add (
const std::vector<size_type> &indices,
706 const ::Vector<OtherNumber> &values);
712 template <
typename OtherNumber>
713 void add (
const size_type n_elements,
714 const size_type *indices,
715 const OtherNumber *values);
721 void sadd (
const Number s,
729 void sadd (
const Number s,
764 std::pair<size_type, size_type>
local_range () const DEAL_II_DEPRECATED;
772 bool in_local_range (const size_type global_index) const DEAL_II_DEPRECATED;
812 const_iterator
begin () const;
824 const_iterator
end () const;
835 Number operator () (const size_type global_index) const;
846 Number &operator () (const size_type global_index);
855 Number operator [] (const size_type global_index) const;
863 Number &operator [] (const size_type global_index);
891 template <typename OtherNumber>
893 std::vector<OtherNumber> &values) const;
899 template <typename ForwardIterator, typename OutputIterator>
901 const ForwardIterator indices_end,
902 OutputIterator values_begin) const;
919 real_type
lp_norm (const real_type p) const;
982 Number, Number,
unsigned int,
984 << " the element received from a remote processor, value "
985 <<
std::setprecision(16) << arg1
986 << ", does not match with the value "
987 <<
std::setprecision(16) << arg2
988 << " on the owner processor " << arg3);
994 size_type, size_type, size_type, size_type,
995 << "You tried to access element " << arg1
996 << " of a distributed vector, but this element is not "
997 << "stored on the current processor. Note: The range of "
998 << "locally owned elements is " << arg2 << " to "
999 << arg3 << ", and there are " << arg4 << " ghost elements "
1000 << "that this vector can access.");
1011 template <typename Number2>
1045 const
Vector<Number> &W);
1086 #ifdef DEAL_II_WITH_MPI 1120 void resize_val (
const size_type new_allocated_size);
1125 template <
typename Number2>
friend class Vector;
1140 template <
typename Number>
1145 return vector_is_ghosted;
1150 template <
typename Number>
1152 typename Vector<Number>::size_type
1155 return partitioner->size();
1160 template <
typename Number>
1162 typename Vector<Number>::size_type
1165 return partitioner->local_size();
1170 template <
typename Number>
1172 std::pair<typename Vector<Number>::size_type,
1173 typename Vector<Number>::size_type>
1176 return partitioner->local_range();
1181 template <
typename Number>
1185 (
const size_type global_index)
const 1187 return partitioner->in_local_range (global_index);
1192 template <
typename Number>
1199 is.add_range (partitioner->local_range().first,
1200 partitioner->local_range().second);
1207 template <
typename Number>
1209 typename Vector<Number>::size_type
1212 return partitioner->n_ghost_indices();
1217 template <
typename Number>
1222 return partitioner->ghost_indices();
1227 template <
typename Number>
1232 return partitioner->is_ghost_entry (global_index);
1237 template <
typename Number>
1239 typename Vector<Number>::iterator
1247 template <
typename Number>
1249 typename Vector<Number>::const_iterator
1257 template <
typename Number>
1259 typename Vector<Number>::iterator
1262 return val+partitioner->local_size();
1267 template <
typename Number>
1269 typename Vector<Number>::const_iterator
1272 return val+partitioner->local_size();
1277 template <
typename Number>
1282 Assert (partitioner->in_local_range (global_index) ||
1283 partitioner->ghost_indices().is_element(global_index),
1284 ExcAccessToNonLocalElement(global_index, partitioner->local_range().first,
1285 partitioner->local_range().second,
1286 partitioner->ghost_indices().n_elements()));
1288 Assert (partitioner->in_local_range (global_index) || vector_is_ghosted ==
true,
1289 ExcMessage(
"You tried to read a ghost element of this vector, " 1290 "but it has not imported its ghost values."));
1291 return val[partitioner->global_to_local(global_index)];
1296 template <
typename Number>
1301 Assert (partitioner->in_local_range (global_index) ||
1302 partitioner->ghost_indices().is_element(global_index),
1303 ExcAccessToNonLocalElement(global_index, partitioner->local_range().first,
1304 partitioner->local_range().second,
1305 partitioner->ghost_indices().n_elements()));
1312 return val[partitioner->global_to_local (global_index)];
1317 template <
typename Number>
1327 template <
typename Number>
1337 template <
typename Number>
1343 partitioner->local_size()+
1344 partitioner->n_ghost_indices());
1346 Assert (local_index < local_size() || vector_is_ghosted ==
true,
1347 ExcMessage(
"You tried to read a ghost element of this vector, " 1348 "but it has not imported its ghost values."));
1349 return val[local_index];
1354 template <
typename Number>
1360 partitioner->local_size()+
1361 partitioner->n_ghost_indices());
1362 return val[local_index];
1367 template <
typename Number>
1368 template <
typename OtherNumber>
1371 std::vector<OtherNumber> &values)
const 1373 for (size_type i = 0; i < indices.size(); ++i)
1374 values[i] =
operator()(indices[i]);
1379 template <
typename Number>
1380 template <
typename ForwardIterator,
typename OutputIterator>
1383 const ForwardIterator indices_end,
1384 OutputIterator values_begin)
const 1386 while (indices_begin != indices_end)
1396 template <
typename Number>
1397 template <
typename OtherNumber>
1401 const ::Vector<OtherNumber> &values)
1404 for (size_type i=0; i<indices.size(); ++i)
1407 ExcMessage(
"The given value is not finite but either infinite or Not A Number (NaN)"));
1414 template <
typename Number>
1415 template <
typename OtherNumber>
1419 const size_type *indices,
1420 const OtherNumber *values)
1422 for (size_type i=0; i<n_elements; ++i, ++indices, ++values)
1425 ExcMessage(
"The given value is not finite but either infinite or Not A Number (NaN)"));
1432 template <
typename Number>
1437 return partitioner->get_mpi_communicator();
1442 template <
typename Number>
1444 const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &
1464 template <
typename Number>
1478 template <
typename Number>
1484 DEAL_II_NAMESPACE_CLOSE
virtual Number add_and_dot(const Number a, const VectorSpaceVector< Number > &V, const VectorSpaceVector< Number > &W)
void compress(::VectorOperation::values operation)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
#define AssertDimension(dim1, dim2)
virtual Number mean_value() const
std::vector< MPI_Request > compress_requests
virtual Vector< Number > & operator*=(const Number factor)
virtual void equ(const Number a, const VectorSpaceVector< Number > &V)
Number inner_product_local(const Vector< Number2 > &V) const
bool in_local_range(const size_type global_index) const 1
virtual size_type size() const
virtual void scale(const VectorSpaceVector< Number > &scaling_factors)
real_type lp_norm(const real_type p) const
real_type l1_norm_local() const
void swap(LinearAlgebra::distributed::Vector< Number > &u, LinearAlgebra::distributed::Vector< Number > &v)
virtual std::size_t memory_consumption() const
#define AssertIndexRange(index, range)
bool is_ghost_entry(const types::global_dof_index global_index) const 1
Number mean_value_local() const
std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > partitioner
virtual ::IndexSet locally_owned_elements() const
static ::ExceptionBase & ExcVectorTypeNotCompatible()
void clear_mpi_requests()
void update_ghost_values() const
bool is_finite(const double x)
virtual void sadd(const Number s, const Number a, const VectorSpaceVector< Number > &V)
std::vector< MPI_Request > update_ghost_values_requests
virtual Vector< Number > & operator/=(const Number factor)
void swap(Vector< Number > &v)
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual void add(const Number a)
bool has_ghost_elements() const
virtual real_type l2_norm() const
unsigned int global_dof_index
virtual Vector< Number > & operator+=(const VectorSpaceVector< Number > &V)
#define Assert(cond, exc)
void resize_val(const size_type new_allocated_size)
const IndexSet & ghost_elements() const 1
virtual Number operator*(const VectorSpaceVector< Number > &V) const
real_type linfty_norm_local() const
Number operator[](const size_type global_index) const
#define DeclException0(Exception0)
virtual void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
size_type local_size() const
virtual real_type l1_norm() const
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< OtherNumber > &values) const
const MPI_Comm & get_mpi_communicator() const
bool all_zero_local() const
Number add_and_dot_local(const Number a, const Vector< Number > &V, const Vector< Number > &W)
Vector< Number > & operator=(const Vector< Number > &in_vector)
void compress_start(const unsigned int communication_channel=0, ::VectorOperation::values operation=VectorOperation::add)
std_cxx11::shared_ptr< ::parallel::internal::TBBPartitioner > thread_loop_partitioner
std::pair< size_type, size_type > local_range() const 1
void compress_finish(::VectorOperation::values operation)
virtual real_type linfty_norm() const
Number operator()(const size_type global_index) const
bool partitioners_are_globally_compatible(const Utilities::MPI::Partitioner &part) const
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
static ::ExceptionBase & ExcNonMatchingElements(Number arg1, Number arg2, unsigned int arg3)
value_type operator()(const size_type i) const
bool partitioners_are_compatible(const Utilities::MPI::Partitioner &part) const
real_type lp_norm_local(const real_type p) const
Number local_element(const size_type local_index) const
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
real_type norm_sqr_local() const
const std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
size_type n_ghost_entries() const 1
virtual Vector< Number > & operator-=(const VectorSpaceVector< Number > &V)
static const bool supports_distributed_data
void update_ghost_values_start(const unsigned int communication_channel=0) const
void update_ghost_values_finish() const