15#ifndef dealii_block_vector_base_h
16#define dealii_block_vector_base_h
49 using has_block_t =
decltype(std::declval<const T>().block(0));
79template <
typename VectorType>
93template <
typename VectorType>
103 namespace BlockVectorIterators
120 template <
typename BlockVectorType,
bool Constness>
135 std::conditional_t<Constness,
136 const typename BlockVectorType::value_type,
137 typename BlockVectorType::value_type>;
158 std::conditional_t<Constness,
160 typename BlockVectorType::BlockType::reference>;
167 std::conditional_t<Constness, const BlockVectorType, BlockVectorType>;
262 template <
bool OtherConstness>
270 template <
bool OtherConstness>
279 template <
bool OtherConstness>
286 template <
bool OtherConstness>
293 template <
bool OtherConstness>
300 template <
bool OtherConstness>
307 template <
bool OtherConstness>
349 "Your program tried to compare iterators pointing to "
350 "different block vectors. There is no reasonable way "
396 template <
typename,
bool>
440template <
typename VectorType>
442 public ReadVector<typename VectorType::value_type>
650 template <typename OtherNumber>
653 std::vector<OtherNumber> &values) const;
686 template <typename ForwardIterator, typename OutputIterator>
689 const ForwardIterator indices_end,
690 OutputIterator values_begin) const;
716 template <typename VectorType2>
724 operator=(const VectorType &v);
730 template <typename VectorType2>
839 template <typename Number>
847 template <typename Number>
856 template <typename Number>
860 const Number *values);
934 template <class BlockVector2>
941 template <class BlockVector2>
972 template <typename N,
bool C>
973 friend class ::
internal::BlockVectorIterators::Iterator;
988 namespace BlockVectorIterators
990 template <
typename BlockVectorType,
bool Constness>
991 inline Iterator<BlockVectorType, Constness>::Iterator(
992 const Iterator<BlockVectorType, Constness> &c)
994 , global_index(c.global_index)
995 , current_block(c.current_block)
996 , index_within_block(c.index_within_block)
997 , next_break_forward(c.next_break_forward)
998 , next_break_backward(c.next_break_backward)
1003 template <
typename BlockVectorType,
bool Constness>
1004 inline Iterator<BlockVectorType, Constness>::Iterator(
1005 const Iterator<BlockVectorType, !Constness> &c)
1007 , global_index(c.global_index)
1008 , current_block(c.current_block)
1009 , index_within_block(c.index_within_block)
1010 , next_break_forward(c.next_break_forward)
1011 , next_break_backward(c.next_break_backward)
1016 static_assert(Constness ==
true,
1017 "Constructing a non-const iterator from a const iterator "
1018 "does not make sense.");
1023 template <
typename BlockVectorType,
bool Constness>
1024 inline Iterator<BlockVectorType, Constness>::Iterator(
1032 , global_index(global_index)
1033 , current_block(current_block)
1034 , index_within_block(index_within_block)
1035 , next_break_forward(next_break_forward)
1036 , next_break_backward(next_break_backward)
1041 template <
typename BlockVectorType,
bool Constness>
1042 inline Iterator<BlockVectorType, Constness> &
1043 Iterator<BlockVectorType, Constness>::operator=(
const Iterator &c)
1046 global_index = c.global_index;
1047 index_within_block = c.index_within_block;
1048 current_block = c.current_block;
1049 next_break_forward = c.next_break_forward;
1050 next_break_backward = c.next_break_backward;
1057 template <
typename BlockVectorType,
bool Constness>
1058 inline typename Iterator<BlockVectorType, Constness>::dereference_type
1059 Iterator<BlockVectorType, Constness>::operator*()
const
1061 return parent->
block(current_block)(index_within_block);
1066 template <
typename BlockVectorType,
bool Constness>
1067 inline typename Iterator<BlockVectorType, Constness>::dereference_type
1068 Iterator<BlockVectorType, Constness>::operator[](
1069 const difference_type d)
const
1076 if ((global_index + d >= next_break_backward) &&
1077 (global_index + d <= next_break_forward))
1078 return parent->
block(current_block)(index_within_block + d);
1087 return (*parent)(global_index + d);
1092 template <
typename BlockVectorType,
bool Constness>
1093 inline Iterator<BlockVectorType, Constness> &
1094 Iterator<BlockVectorType, Constness>::operator++()
1102 template <
typename BlockVectorType,
bool Constness>
1103 inline Iterator<BlockVectorType, Constness>
1104 Iterator<BlockVectorType, Constness>::operator++(
int)
1106 const Iterator old_value = *
this;
1113 template <
typename BlockVectorType,
bool Constness>
1114 inline Iterator<BlockVectorType, Constness> &
1115 Iterator<BlockVectorType, Constness>::operator--()
1123 template <
typename BlockVectorType,
bool Constness>
1124 inline Iterator<BlockVectorType, Constness>
1125 Iterator<BlockVectorType, Constness>::operator--(
int)
1127 const Iterator old_value = *
this;
1134 template <
typename BlockVectorType,
bool Constness>
1135 template <
bool OtherConstness>
1137 Iterator<BlockVectorType, Constness>::operator==(
1138 const Iterator<BlockVectorType, OtherConstness> &i)
const
1140 Assert(parent == i.parent, ExcPointerToDifferentVectors());
1142 return (global_index == i.global_index);
1147 template <
typename BlockVectorType,
bool Constness>
1148 template <
bool OtherConstness>
1150 Iterator<BlockVectorType, Constness>::operator!=(
1151 const Iterator<BlockVectorType, OtherConstness> &i)
const
1153 Assert(parent == i.parent, ExcPointerToDifferentVectors());
1155 return (global_index != i.global_index);
1160 template <
typename BlockVectorType,
bool Constness>
1161 template <
bool OtherConstness>
1163 Iterator<BlockVectorType, Constness>::operator<(
1164 const Iterator<BlockVectorType, OtherConstness> &i)
const
1166 Assert(parent == i.parent, ExcPointerToDifferentVectors());
1168 return (global_index < i.global_index);
1173 template <
typename BlockVectorType,
bool Constness>
1174 template <
bool OtherConstness>
1176 Iterator<BlockVectorType, Constness>::operator<=(
1177 const Iterator<BlockVectorType, OtherConstness> &i)
const
1179 Assert(parent == i.parent, ExcPointerToDifferentVectors());
1181 return (global_index <= i.global_index);
1186 template <
typename BlockVectorType,
bool Constness>
1187 template <
bool OtherConstness>
1189 Iterator<BlockVectorType, Constness>::operator>(
1190 const Iterator<BlockVectorType, OtherConstness> &i)
const
1192 Assert(parent == i.parent, ExcPointerToDifferentVectors());
1194 return (global_index > i.global_index);
1199 template <
typename BlockVectorType,
bool Constness>
1200 template <
bool OtherConstness>
1202 Iterator<BlockVectorType, Constness>::operator>=(
1203 const Iterator<BlockVectorType, OtherConstness> &i)
const
1205 Assert(parent == i.parent, ExcPointerToDifferentVectors());
1207 return (global_index >= i.global_index);
1212 template <
typename BlockVectorType,
bool Constness>
1213 template <
bool OtherConstness>
1214 inline typename Iterator<BlockVectorType, Constness>::difference_type
1215 Iterator<BlockVectorType, Constness>::operator-(
1216 const Iterator<BlockVectorType, OtherConstness> &i)
const
1218 Assert(parent == i.parent, ExcPointerToDifferentVectors());
1220 return (
static_cast<signed int>(global_index) -
1221 static_cast<signed int>(i.global_index));
1226 template <
typename BlockVectorType,
bool Constness>
1227 inline Iterator<BlockVectorType, Constness>
1228 Iterator<BlockVectorType, Constness>::operator+(
1229 const difference_type &d)
const
1236 if ((global_index + d >= next_break_backward) &&
1237 (global_index + d <= next_break_forward))
1238 return Iterator(*parent,
1241 index_within_block + d,
1243 next_break_backward);
1248 return Iterator(*parent, global_index + d);
1253 template <
typename BlockVectorType,
bool Constness>
1254 inline Iterator<BlockVectorType, Constness>
1255 Iterator<BlockVectorType, Constness>::operator-(
1256 const difference_type &d)
const
1263 if ((global_index - d >= next_break_backward) &&
1264 (global_index - d <= next_break_forward))
1265 return Iterator(*parent,
1268 index_within_block - d,
1270 next_break_backward);
1275 return Iterator(*parent, global_index - d);
1280 template <
typename BlockVectorType,
bool Constness>
1281 inline Iterator<BlockVectorType, Constness> &
1282 Iterator<BlockVectorType, Constness>::operator+=(
const difference_type &d)
1289 if ((global_index + d >= next_break_backward) &&
1290 (global_index + d <= next_break_forward))
1293 index_within_block += d;
1299 *
this = Iterator(*parent, global_index + d);
1306 template <
typename BlockVectorType,
bool Constness>
1307 inline Iterator<BlockVectorType, Constness> &
1308 Iterator<BlockVectorType, Constness>::operator-=(
const difference_type &d)
1315 if ((global_index - d >= next_break_backward) &&
1316 (global_index - d <= next_break_forward))
1319 index_within_block -= d;
1325 *
this = Iterator(*parent, global_index - d);
1331 template <
typename BlockVectorType,
bool Constness>
1332 Iterator<BlockVectorType, Constness>::Iterator(
BlockVector &parent,
1335 , global_index(global_index)
1343 if (global_index < parent.
size())
1345 const std::pair<size_type, size_type> indices =
1347 current_block = indices.first;
1348 index_within_block = indices.second;
1350 next_break_backward =
1352 next_break_forward =
1360 this->global_index = parent.
size();
1362 index_within_block = 0;
1363 next_break_backward = global_index;
1370 template <
typename BlockVectorType,
bool Constness>
1372 Iterator<BlockVectorType, Constness>::move_forward()
1374 if (global_index != next_break_forward)
1375 ++index_within_block;
1380 index_within_block = 0;
1385 next_break_backward = next_break_forward + 1;
1388 if (current_block < parent->block_indices.
size())
1389 next_break_forward +=
1404 template <
typename BlockVectorType,
bool Constness>
1406 Iterator<BlockVectorType, Constness>::move_backward()
1408 if (global_index != next_break_backward)
1409 --index_within_block;
1410 else if (current_block != 0)
1415 index_within_block =
1420 next_break_forward = next_break_backward - 1;
1423 next_break_backward -=
1432 next_break_forward = 0;
1433 next_break_backward = 0;
1446template <
typename VectorType>
1455template <
typename VectorType>
1459 std::size_t local_size = 0;
1460 for (
unsigned int b = 0;
b <
n_blocks(); ++
b)
1461 local_size +=
block(b).locally_owned_size();
1467template <
typename VectorType>
1475 for (
unsigned int b = 0;
b <
n_blocks(); ++
b)
1488template <
typename VectorType>
1496template <
typename VectorType>
1507template <
typename VectorType>
1518template <
typename VectorType>
1526template <
typename VectorType>
1530 std::vector<size_type> sizes(
n_blocks());
1532 for (size_type i = 0; i <
n_blocks(); ++i)
1533 sizes[i] =
block(i).size();
1540template <
typename VectorType>
1544 for (
unsigned int i = 0; i <
n_blocks(); ++i)
1545 block(i).compress(operation);
1550template <
typename VectorType>
1559template <
typename VectorType>
1567template <
typename VectorType>
1576template <
typename VectorType>
1584template <
typename VectorType>
1588 const std::pair<size_type, size_type> local_index =
1591 return components[local_index.first].in_local_range(global_index);
1595template <
typename VectorType>
1599 for (size_type i = 0; i <
n_blocks(); ++i)
1608template <
typename VectorType>
1612 for (size_type i = 0; i <
n_blocks(); ++i)
1621template <
typename VectorType>
1629 value_type
sum = 0.;
1630 for (size_type i = 0; i <
n_blocks(); ++i)
1637template <
typename VectorType>
1642 for (size_type i = 0; i <
n_blocks(); ++i)
1650template <
typename VectorType>
1654 value_type
sum = 0.;
1657 for (size_type i = 0; i <
n_blocks(); ++i)
1667template <
typename VectorType>
1672 for (size_type i = 0; i <
n_blocks(); ++i)
1680template <
typename VectorType>
1689template <
typename VectorType>
1694 for (size_type i = 0; i <
n_blocks(); ++i)
1696 value_type newval =
components[i].linfty_norm();
1705template <
typename VectorType>
1715 value_type
sum = 0.;
1716 for (size_type i = 0; i <
n_blocks(); ++i)
1724template <
typename VectorType>
1731 for (size_type i = 0; i <
n_blocks(); ++i)
1741template <
typename VectorType>
1748 for (size_type i = 0; i <
n_blocks(); ++i)
1757template <
typename VectorType>
1758template <
typename Number>
1761 const std::vector<Number> &values)
1763 Assert(indices.size() == values.size(),
1765 add(indices.size(), indices.data(), values.data());
1770template <
typename VectorType>
1771template <
typename Number>
1776 Assert(indices.size() == values.size(),
1778 const size_type n_indices = indices.size();
1779 for (size_type i = 0; i < n_indices; ++i)
1780 (*
this)(indices[i]) += values(i);
1785template <
typename VectorType>
1786template <
typename Number>
1789 const size_type *indices,
1790 const Number *values)
1792 for (size_type i = 0; i < n_indices; ++i)
1793 (*
this)(indices[i]) += values[i];
1798template <
typename VectorType>
1804 for (size_type i = 0; i <
n_blocks(); ++i)
1812template <
typename VectorType>
1822 for (size_type i = 0; i <
n_blocks(); ++i)
1830template <
typename VectorType>
1846 for (size_type i = 0; i <
n_blocks(); ++i)
1854template <
typename VectorType>
1864 for (size_type i = 0; i <
n_blocks(); ++i)
1872template <
typename VectorType>
1884 for (size_type i = 0; i <
n_blocks(); ++i)
1892template <
typename VectorType>
1909 for (size_type i = 0; i <
n_blocks(); ++i)
1917template <
typename VectorType>
1939 for (size_type i = 0; i <
n_blocks(); ++i)
1948template <
typename VectorType>
1949template <
class BlockVector2>
1955 for (size_type i = 0; i <
n_blocks(); ++i)
1961template <
typename VectorType>
1971template <
typename VectorType>
1972template <
class BlockVector2>
1981 for (size_type i = 0; i <
n_blocks(); ++i)
1987template <
typename VectorType>
1991 for (size_type i = 0; i <
n_blocks(); ++i)
1992 block(i).update_ghost_values();
1997template <
typename VectorType>
2003 for (size_type i = 0; i <
n_blocks(); ++i)
2010template <
typename VectorType>
2016 for (size_type i = 0; i <
n_blocks(); ++i)
2023template <
typename VectorType>
2024template <
typename VectorType2>
2030 for (size_type i = 0; i <
n_blocks(); ++i)
2038template <
typename VectorType>
2044 size_type index_v = 0;
2046 for (size_type i = 0; i <
block(b).size(); ++i, ++index_v)
2047 block(b)(i) = v(index_v);
2054template <
typename VectorType>
2055template <
typename VectorType2>
2062 for (size_type i = 0; i <
n_blocks(); ++i)
2071template <
typename VectorType>
2077 for (size_type i = 0; i <
n_blocks(); ++i)
2085template <
typename VectorType>
2092 const value_type inverse_factor =
value_type(1.) / factor;
2094 for (size_type i = 0; i <
n_blocks(); ++i)
2101template <
typename VectorType>
2105 const std::pair<unsigned int, size_type> local_index =
2107 return components[local_index.first](local_index.second);
2112template <
typename VectorType>
2116 const std::pair<unsigned int, size_type> local_index =
2118 return components[local_index.first](local_index.second);
2123template <
typename VectorType>
2132template <
typename VectorType>
2141template <
typename VectorType>
2142template <
typename OtherNumber>
2145 const std::vector<size_type> &indices,
2146 std::vector<OtherNumber> &values)
const
2148 for (size_type i = 0; i < indices.size(); ++i)
2149 values[i] =
operator()(indices[i]);
2154template <
typename VectorType>
2161 for (
unsigned int i = 0; i < indices.
size(); ++i)
2164 entries[i] = (*this)[indices[i]];
2170template <
typename VectorType>
2171template <
typename ForwardIterator,
typename OutputIterator>
2174 ForwardIterator indices_begin,
2175 const ForwardIterator indices_end,
2176 OutputIterator values_begin)
const
2178 while (indices_begin != indices_end)
size_type block_size(const unsigned int i) const
size_type total_size() const
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
unsigned int size() const
size_type local_to_global(const unsigned int block, const size_type index) const
size_type block_start(const unsigned int i) const
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
BlockVectorBase & operator+=(const BlockVectorBase &V)
BlockVectorBase()=default
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< OtherNumber > &values) const
value_type mean_value() const
virtual size_type size() const override
void add(const std::vector< size_type > &indices, const std::vector< Number > &values)
::internal::BlockVectorIterators::Iterator< BlockVectorBase, false > iterator
BlockVectorBase & operator/=(const value_type factor)
unsigned int n_blocks() const
real_type norm_sqr() const
::internal::BlockVectorIterators::Iterator< BlockVectorBase, true > const_iterator
const value_type * const_pointer
typename BlockType::const_reference const_reference
void update_ghost_values() const
std::size_t memory_consumption() const
real_type l1_norm() const
real_type linfty_norm() const
value_type add_and_dot(const value_type a, const BlockVectorBase &V, const BlockVectorBase &W)
typename BlockType::value_type value_type
value_type operator*(const BlockVectorBase &V) const
void compress(VectorOperation::values operation)
BlockVectorBase & operator=(const value_type s)
BlockVectorBase & operator*=(const value_type factor)
void sadd(const value_type s, const BlockVectorBase &V)
std::vector< VectorType > components
bool in_local_range(const size_type global_index) const
BlockVectorBase(BlockVectorBase &&) noexcept=default
BlockVectorBase(const BlockVectorBase &)=default
typename BlockType::reference reference
value_type operator[](const size_type i) const
BlockIndices block_indices
typename BlockType::real_type real_type
void scale(const BlockVector2 &v)
BlockType & block(const unsigned int i)
real_type l2_norm() const
value_type operator()(const size_type i) const
IndexSet locally_owned_elements() const
bool is_non_negative() const
BlockVectorBase & operator-=(const BlockVectorBase &V)
bool operator==(const BlockVectorBase< VectorType2 > &v) const
std::size_t locally_owned_size() const
const BlockIndices & get_block_indices() const
void equ(const value_type a, const BlockVector2 &V)
bool operator==(const Iterator< BlockVectorType, OtherConstness > &i) const
dereference_type operator[](const difference_type d) const
std::conditional_t< Constness, value_type, typename BlockVectorType::BlockType::reference > dereference_type
Iterator & operator=(const Iterator &c)
bool operator<(const Iterator< BlockVectorType, OtherConstness > &i) const
Iterator(BlockVector &parent, const size_type global_index, const size_type current_block, const size_type index_within_block, const size_type next_break_forward, const size_type next_break_backward)
Iterator(const Iterator< BlockVectorType, !Constness > &c)
Iterator(const Iterator &c)
dereference_type operator*() const
Iterator operator-(const difference_type &d) const
difference_type operator-(const Iterator< BlockVectorType, OtherConstness > &i) const
std::random_access_iterator_tag iterator_category
Iterator(BlockVector &parent, const size_type global_index)
size_type next_break_forward
bool operator>=(const Iterator< BlockVectorType, OtherConstness > &i) const
unsigned int current_block
Iterator & operator-=(const difference_type &d)
size_type next_break_backward
std::conditional_t< Constness, const typename BlockVectorType::value_type, typename BlockVectorType::value_type > value_type
std::ptrdiff_t difference_type
size_type index_within_block
bool operator<=(const Iterator< BlockVectorType, OtherConstness > &i) const
typename BlockVectorType::reference reference
Iterator operator+(const difference_type &d) const
Iterator & operator+=(const difference_type &d)
bool operator!=(const Iterator< BlockVectorType, OtherConstness > &i) const
bool operator>(const Iterator< BlockVectorType, OtherConstness > &i) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcPointerToDifferentVectors()
#define AssertIsFinite(number)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcDifferentBlockIndices()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
T sum(const T &t, const MPI_Comm mpi_communicator)
constexpr bool has_n_blocks
decltype(std::declval< const T >().block(0)) has_block_t
constexpr bool is_supported_operation
decltype(std::declval< const T >().n_blocks()) has_n_blocks_t
constexpr bool is_block_vector
const types::global_dof_index invalid_size_type
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index