16 #ifndef dealii_block_matrix_base_h 17 #define dealii_block_matrix_base_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/memory_consumption.h> 23 #include <deal.II/base/smartpointer.h> 24 #include <deal.II/base/table.h> 25 #include <deal.II/base/thread_management.h> 26 #include <deal.II/base/utilities.h> 28 #include <deal.II/lac/block_indices.h> 29 #include <deal.II/lac/exceptions.h> 30 #include <deal.II/lac/full_matrix.h> 31 #include <deal.II/lac/matrix_iterator.h> 32 #include <deal.II/lac/vector.h> 33 #include <deal.II/lac/vector_operation.h> 37 DEAL_II_NAMESPACE_OPEN
60 template <
class BlockMatrixType>
114 template <
class BlockMatrixType,
bool Constness>
121 template <
class BlockMatrixType>
195 operator==(
const Accessor &a)
const;
199 friend class Accessor<BlockMatrixType, true>;
206 template <
class BlockMatrixType>
281 operator==(
const Accessor &a)
const;
287 friend class ::MatrixIterator;
353 template <
typename MatrixType>
406 template <
class BlockMatrixType>
408 copy_from(
const BlockMatrixType &source);
414 block(
const unsigned int row,
const unsigned int column);
422 block(
const unsigned int row,
const unsigned int column)
const;
459 set(
const size_type i,
const size_type j,
const value_type value);
476 template <
typename number>
478 set(
const std::vector<size_type> &indices,
480 const bool elide_zero_values =
false);
487 template <
typename number>
489 set(
const std::vector<size_type> &row_indices,
490 const std::vector<size_type> &col_indices,
492 const bool elide_zero_values =
false);
504 template <
typename number>
506 set(
const size_type row,
507 const std::vector<size_type> &col_indices,
508 const std::vector<number> & values,
509 const bool elide_zero_values =
false);
520 template <
typename number>
522 set(
const size_type row,
523 const size_type n_cols,
524 const size_type *col_indices,
525 const number * values,
526 const bool elide_zero_values =
false);
534 add(
const size_type i,
const size_type j,
const value_type value);
550 template <
typename number>
552 add(
const std::vector<size_type> &indices,
554 const bool elide_zero_values =
true);
561 template <
typename number>
563 add(
const std::vector<size_type> &row_indices,
564 const std::vector<size_type> &col_indices,
566 const bool elide_zero_values =
true);
577 template <
typename number>
579 add(
const size_type row,
580 const std::vector<size_type> &col_indices,
581 const std::vector<number> & values,
582 const bool elide_zero_values =
true);
593 template <
typename number>
595 add(
const size_type row,
596 const size_type n_cols,
597 const size_type *col_indices,
598 const number * values,
599 const bool elide_zero_values =
true,
600 const bool col_indices_are_sorted =
false);
623 operator()(
const size_type i,
const size_type j)
const;
634 el(
const size_type i,
const size_type j)
const;
676 template <
class BlockVectorType>
678 vmult_add(BlockVectorType &dst,
const BlockVectorType &src)
const;
685 template <
class BlockVectorType>
687 Tvmult_add(BlockVectorType &dst,
const BlockVectorType &src)
const;
701 template <
class BlockVectorType>
708 template <
class BlockVectorType>
711 const BlockVectorType &v)
const;
716 template <
class BlockVectorType>
719 const BlockVectorType &x,
720 const BlockVectorType &b)
const;
729 print(std::ostream &out,
const bool alternative_output =
false)
const;
747 begin(
const size_type r);
753 end(
const size_type r);
770 begin(
const size_type r)
const;
776 end(
const size_type r)
const;
811 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" << arg3
812 <<
',' << arg4 <<
"] have differing row numbers.");
821 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" << arg3
822 <<
',' << arg4 <<
"] have differing column numbers.");
882 template <
class BlockVectorType>
896 template <
class BlockVectorType,
class VectorType>
910 template <
class BlockVectorType,
class VectorType>
924 template <
class VectorType>
939 template <
class BlockVectorType>
953 template <
class BlockVectorType,
class VectorType>
967 template <
class BlockVectorType,
class VectorType>
981 template <
class VectorType>
1067 template <
typename,
bool>
1083 template <
class BlockMatrixType>
1090 template <
class BlockMatrixType>
1092 AccessorBase<BlockMatrixType>::block_row()
const 1100 template <
class BlockMatrixType>
1102 AccessorBase<BlockMatrixType>::block_column()
const 1110 template <
class BlockMatrixType>
1111 inline Accessor<BlockMatrixType, true>::Accessor(
1112 const BlockMatrixType *matrix,
1113 const size_type row,
1114 const size_type col)
1116 , base_iterator(
matrix->block(0, 0).begin())
1123 if (row < matrix->
m())
1125 const std::pair<unsigned int, size_type> indices =
1126 matrix->row_block_indices.global_to_local(row);
1130 for (
unsigned int bc = 0; bc <
matrix->n_block_cols(); ++bc)
1133 matrix->block(indices.first, bc).begin(indices.second);
1134 if (base_iterator !=
1135 matrix->block(indices.first, bc).end(indices.second))
1137 this->row_block = indices.first;
1138 this->col_block = bc;
1149 *
this = Accessor(matrix, row + 1, 0);
1174 template <
class BlockMatrixType>
1175 inline Accessor<BlockMatrixType, true>::Accessor(
1176 const Accessor<BlockMatrixType, false> &other)
1178 , base_iterator(other.base_iterator)
1180 this->row_block = other.row_block;
1181 this->col_block = other.col_block;
1185 template <
class BlockMatrixType>
1186 inline typename Accessor<BlockMatrixType, true>::size_type
1187 Accessor<BlockMatrixType, true>::row()
const 1192 return (
matrix->row_block_indices.local_to_global(this->row_block, 0) +
1193 base_iterator->row());
1197 template <
class BlockMatrixType>
1198 inline typename Accessor<BlockMatrixType, true>::size_type
1199 Accessor<BlockMatrixType, true>::column()
const 1204 return (
matrix->column_block_indices.local_to_global(this->col_block, 0) +
1205 base_iterator->column());
1209 template <
class BlockMatrixType>
1210 inline typename Accessor<BlockMatrixType, true>::value_type
1211 Accessor<BlockMatrixType, true>::value()
const 1218 return base_iterator->value();
1223 template <
class BlockMatrixType>
1225 Accessor<BlockMatrixType, true>::advance()
1233 size_type local_row = base_iterator->row();
1244 while (base_iterator ==
1245 matrix->block(this->row_block, this->col_block).end(local_row))
1250 if (this->col_block < matrix->n_block_cols() - 1)
1254 matrix->block(this->row_block, this->col_block).begin(local_row);
1260 this->col_block = 0;
1270 matrix->block(this->row_block, this->col_block).m())
1274 if (this->row_block ==
matrix->n_block_rows())
1283 matrix->block(this->row_block, this->col_block).begin(local_row);
1289 template <
class BlockMatrixType>
1291 Accessor<BlockMatrixType, true>::operator==(
const Accessor &a)
const 1293 if (matrix != a.matrix)
1296 if (this->row_block == a.row_block && this->col_block == a.col_block)
1303 (base_iterator == a.base_iterator));
1311 template <
class BlockMatrixType>
1312 inline Accessor<BlockMatrixType, false>::Accessor(BlockMatrixType *matrix,
1313 const size_type row,
1314 const size_type col)
1316 , base_iterator(
matrix->block(0, 0).begin())
1322 if (row < matrix->
m())
1324 const std::pair<unsigned int, size_type> indices =
1325 matrix->row_block_indices.global_to_local(row);
1329 for (size_type bc = 0; bc <
matrix->n_block_cols(); ++bc)
1332 matrix->block(indices.first, bc).begin(indices.second);
1333 if (base_iterator !=
1334 matrix->block(indices.first, bc).end(indices.second))
1336 this->row_block = indices.first;
1337 this->col_block = bc;
1348 *
this = Accessor(matrix, row + 1, 0);
1360 template <
class BlockMatrixType>
1361 inline typename Accessor<BlockMatrixType, false>::size_type
1362 Accessor<BlockMatrixType, false>::row()
const 1366 return (
matrix->row_block_indices.local_to_global(this->row_block, 0) +
1367 base_iterator->row());
1371 template <
class BlockMatrixType>
1372 inline typename Accessor<BlockMatrixType, false>::size_type
1373 Accessor<BlockMatrixType, false>::column()
const 1377 return (
matrix->column_block_indices.local_to_global(this->col_block, 0) +
1378 base_iterator->column());
1382 template <
class BlockMatrixType>
1383 inline typename Accessor<BlockMatrixType, false>::value_type
1384 Accessor<BlockMatrixType, false>::value()
const 1389 return base_iterator->value();
1394 template <
class BlockMatrixType>
1396 Accessor<BlockMatrixType, false>::set_value(
1397 typename Accessor<BlockMatrixType, false>::value_type newval)
const 1402 base_iterator->value() = newval;
1407 template <
class BlockMatrixType>
1409 Accessor<BlockMatrixType, false>::advance()
1415 size_type local_row = base_iterator->row();
1426 while (base_iterator ==
1427 matrix->block(this->row_block, this->col_block).end(local_row))
1432 if (this->col_block < matrix->n_block_cols() - 1)
1436 matrix->block(this->row_block, this->col_block).begin(local_row);
1442 this->col_block = 0;
1452 matrix->block(this->row_block, this->col_block).m())
1456 if (this->row_block ==
matrix->n_block_rows())
1465 matrix->block(this->row_block, this->col_block).begin(local_row);
1472 template <
class BlockMatrixType>
1474 Accessor<BlockMatrixType, false>::operator==(
const Accessor &a)
const 1476 if (matrix != a.matrix)
1479 if (this->row_block == a.row_block && this->col_block == a.col_block)
1486 (base_iterator == a.base_iterator));
1495 template <
typename MatrixType>
1507 template <
class MatrixType>
1508 template <
class BlockMatrixType>
1512 for (
unsigned int r = 0; r < n_block_rows(); ++r)
1513 for (
unsigned int c = 0; c < n_block_cols(); ++c)
1514 block(r, c).
copy_from(source.block(r, c));
1520 template <
class MatrixType>
1531 sizeof(temporary_data.mutex);
1533 for (
unsigned int r = 0; r < n_block_rows(); ++r)
1534 for (
unsigned int c = 0; c < n_block_cols(); ++c)
1536 MatrixType *p = this->sub_objects[r][c];
1545 template <
class MatrixType>
1549 for (
unsigned int r = 0; r < n_block_rows(); ++r)
1550 for (
unsigned int c = 0; c < n_block_cols(); ++c)
1552 MatrixType *p = this->sub_objects[r][c];
1553 this->sub_objects[r][c] =
nullptr;
1556 sub_objects.reinit(0, 0);
1559 row_block_indices = column_block_indices =
BlockIndices();
1564 template <
class MatrixType>
1567 const unsigned int column)
1572 return *sub_objects[row][column];
1577 template <
class MatrixType>
1580 const unsigned int column)
const 1585 return *sub_objects[row][column];
1589 template <
class MatrixType>
1590 inline typename BlockMatrixBase<MatrixType>::size_type
1593 return row_block_indices.total_size();
1598 template <
class MatrixType>
1599 inline typename BlockMatrixBase<MatrixType>::size_type
1602 return column_block_indices.total_size();
1607 template <
class MatrixType>
1611 return column_block_indices.size();
1616 template <
class MatrixType>
1620 return row_block_indices.size();
1628 template <
class MatrixType>
1632 const value_type value)
1634 prepare_set_operation();
1638 const std::pair<unsigned int, size_type>
1639 row_index = row_block_indices.global_to_local(i),
1640 col_index = column_block_indices.global_to_local(j);
1641 block(row_index.first, col_index.first)
1642 .set(row_index.second, col_index.second, value);
1647 template <
class MatrixType>
1648 template <
typename number>
1651 const std::vector<size_type> &col_indices,
1653 const bool elide_zero_values)
1655 Assert(row_indices.size() == values.
m(),
1657 Assert(col_indices.size() == values.
n(),
1660 for (size_type i = 0; i < row_indices.size(); ++i)
1670 template <
class MatrixType>
1671 template <
typename number>
1675 const bool elide_zero_values)
1677 Assert(indices.size() == values.
m(),
1681 for (size_type i = 0; i < indices.size(); ++i)
1691 template <
class MatrixType>
1692 template <
typename number>
1695 const std::vector<size_type> &col_indices,
1696 const std::vector<number> & values,
1697 const bool elide_zero_values)
1699 Assert(col_indices.size() == values.size(),
1714 template <
class MatrixType>
1715 template <
typename number>
1718 const size_type n_cols,
1719 const size_type *col_indices,
1720 const number * values,
1721 const bool elide_zero_values)
1723 prepare_set_operation();
1727 std::lock_guard<std::mutex> lock(temporary_data.mutex);
1730 if (temporary_data.column_indices.size() < this->n_block_cols())
1732 temporary_data.column_indices.resize(this->n_block_cols());
1733 temporary_data.column_values.resize(this->n_block_cols());
1734 temporary_data.counter_within_block.resize(this->n_block_cols());
1747 if (temporary_data.column_indices[0].size() < n_cols)
1749 for (
unsigned int i = 0; i < this->n_block_cols(); ++i)
1751 temporary_data.column_indices[i].resize(n_cols);
1752 temporary_data.column_values[i].resize(n_cols);
1758 for (
unsigned int i = 0; i < this->n_block_cols(); ++i)
1759 temporary_data.counter_within_block[i] = 0;
1770 for (size_type j = 0; j < n_cols; ++j)
1772 number value = values[j];
1774 if (value == number() && elide_zero_values ==
true)
1777 const std::pair<unsigned int, size_type> col_index =
1778 this->column_block_indices.global_to_local(col_indices[j]);
1780 const size_type local_index =
1781 temporary_data.counter_within_block[col_index.first]++;
1783 temporary_data.column_indices[col_index.first][local_index] =
1785 temporary_data.column_values[col_index.first][local_index] = value;
1791 size_type length = 0;
1792 for (
unsigned int i = 0; i < this->n_block_cols(); ++i)
1793 length += temporary_data.counter_within_block[i];
1802 const std::pair<unsigned int, size_type> row_index =
1803 this->row_block_indices.global_to_local(row);
1804 for (
unsigned int block_col = 0; block_col < n_block_cols(); ++block_col)
1806 if (temporary_data.counter_within_block[block_col] == 0)
1809 block(row_index.first, block_col)
1810 .set(row_index.second,
1811 temporary_data.counter_within_block[block_col],
1812 temporary_data.column_indices[block_col].data(),
1813 temporary_data.column_values[block_col].data(),
1820 template <
class MatrixType>
1824 const value_type value)
1828 prepare_add_operation();
1833 using MatrixTraits =
typename MatrixType::Traits;
1834 if ((MatrixTraits::zero_addition_can_be_elided ==
true) &&
1835 (value == value_type()))
1838 const std::pair<unsigned int, size_type>
1839 row_index = row_block_indices.global_to_local(i),
1840 col_index = column_block_indices.global_to_local(j);
1841 block(row_index.first, col_index.first)
1842 .add(row_index.second, col_index.second, value);
1847 template <
class MatrixType>
1848 template <
typename number>
1851 const std::vector<size_type> &col_indices,
1853 const bool elide_zero_values)
1855 Assert(row_indices.size() == values.
m(),
1857 Assert(col_indices.size() == values.
n(),
1860 for (size_type i = 0; i < row_indices.size(); ++i)
1870 template <
class MatrixType>
1871 template <
typename number>
1875 const bool elide_zero_values)
1877 Assert(indices.size() == values.
m(),
1881 for (size_type i = 0; i < indices.size(); ++i)
1891 template <
class MatrixType>
1892 template <
typename number>
1895 const std::vector<size_type> &col_indices,
1896 const std::vector<number> & values,
1897 const bool elide_zero_values)
1899 Assert(col_indices.size() == values.size(),
1914 template <
class MatrixType>
1915 template <
typename number>
1918 const size_type n_cols,
1919 const size_type *col_indices,
1920 const number * values,
1921 const bool elide_zero_values,
1922 const bool col_indices_are_sorted)
1924 prepare_add_operation();
1929 if (col_indices_are_sorted ==
true)
1934 size_type before = col_indices[0];
1935 for (size_type i = 1; i < n_cols; ++i)
1936 if (col_indices[i] <= before)
1938 ExcMessage(
"Flag col_indices_are_sorted is set, but " 1939 "indices appear to not be sorted.")) else before =
1942 const std::pair<unsigned int, size_type> row_index =
1943 this->row_block_indices.global_to_local(row);
1945 if (this->n_block_cols() > 1)
1947 const size_type *first_block =
1949 col_indices + n_cols,
1950 this->column_block_indices.block_start(1));
1952 const size_type n_zero_block_indices = first_block - col_indices;
1953 block(row_index.first, 0)
1954 .add(row_index.second,
1955 n_zero_block_indices,
1959 col_indices_are_sorted);
1961 if (n_zero_block_indices < n_cols)
1963 n_cols - n_zero_block_indices,
1965 values + n_zero_block_indices,
1971 block(row_index.first, 0)
1972 .add(row_index.second,
1977 col_indices_are_sorted);
1984 std::lock_guard<std::mutex> lock(temporary_data.mutex);
1986 if (temporary_data.column_indices.size() < this->n_block_cols())
1988 temporary_data.column_indices.resize(this->n_block_cols());
1989 temporary_data.column_values.resize(this->n_block_cols());
1990 temporary_data.counter_within_block.resize(this->n_block_cols());
2003 if (temporary_data.column_indices[0].size() < n_cols)
2005 for (
unsigned int i = 0; i < this->n_block_cols(); ++i)
2007 temporary_data.column_indices[i].resize(n_cols);
2008 temporary_data.column_values[i].resize(n_cols);
2014 for (
unsigned int i = 0; i < this->n_block_cols(); ++i)
2015 temporary_data.counter_within_block[i] = 0;
2026 for (size_type j = 0; j < n_cols; ++j)
2028 number value = values[j];
2030 if (value == number() && elide_zero_values ==
true)
2033 const std::pair<unsigned int, size_type> col_index =
2034 this->column_block_indices.global_to_local(col_indices[j]);
2036 const size_type local_index =
2037 temporary_data.counter_within_block[col_index.first]++;
2039 temporary_data.column_indices[col_index.first][local_index] =
2041 temporary_data.column_values[col_index.first][local_index] = value;
2047 size_type length = 0;
2048 for (
unsigned int i = 0; i < this->n_block_cols(); ++i)
2049 length += temporary_data.counter_within_block[i];
2058 const std::pair<unsigned int, size_type> row_index =
2059 this->row_block_indices.global_to_local(row);
2060 for (
unsigned int block_col = 0; block_col < n_block_cols(); ++block_col)
2062 if (temporary_data.counter_within_block[block_col] == 0)
2065 block(row_index.first, block_col)
2066 .add(row_index.second,
2067 temporary_data.counter_within_block[block_col],
2068 temporary_data.column_indices[block_col].data(),
2069 temporary_data.column_values[block_col].data(),
2071 col_indices_are_sorted);
2077 template <
class MatrixType>
2084 prepare_add_operation();
2089 using MatrixTraits =
typename MatrixType::Traits;
2090 if ((MatrixTraits::zero_addition_can_be_elided ==
true) && (factor == 0))
2093 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2094 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2097 block(row, col).add(factor,
matrix.block(row, col));
2102 template <
class MatrixType>
2105 const size_type j)
const 2107 const std::pair<unsigned int, size_type>
2108 row_index = row_block_indices.global_to_local(i),
2109 col_index = column_block_indices.global_to_local(j);
2110 return block(row_index.first, col_index.first)(row_index.second,
2116 template <
class MatrixType>
2120 const std::pair<unsigned int, size_type>
2121 row_index = row_block_indices.global_to_local(i),
2122 col_index = column_block_indices.global_to_local(j);
2123 return block(row_index.first, col_index.first)
2124 .el(row_index.second, col_index.second);
2129 template <
class MatrixType>
2135 const std::pair<unsigned int, size_type> index =
2136 row_block_indices.global_to_local(i);
2137 return block(index.first, index.first).diag_element(index.second);
2142 template <
class MatrixType>
2147 for (
unsigned int r = 0; r < n_block_rows(); ++r)
2148 for (
unsigned int c = 0; c < n_block_cols(); ++c)
2149 block(r, c).compress(operation);
2154 template <
class MatrixType>
2161 for (
unsigned int r = 0; r < n_block_rows(); ++r)
2162 for (
unsigned int c = 0; c < n_block_cols(); ++c)
2163 block(r, c) *= factor;
2170 template <
class MatrixType>
2178 const value_type factor_inv = 1. / factor;
2180 for (
unsigned int r = 0; r < n_block_rows(); ++r)
2181 for (
unsigned int c = 0; c < n_block_cols(); ++c)
2182 block(r, c) *= factor_inv;
2189 template <
class MatrixType>
2193 return this->row_block_indices;
2198 template <
class MatrixType>
2202 return this->column_block_indices;
2207 template <
class MatrixType>
2208 template <
class BlockVectorType>
2211 const BlockVectorType &src)
const 2213 Assert(dst.n_blocks() == n_block_rows(),
2215 Assert(src.n_blocks() == n_block_cols(),
2218 for (size_type row = 0; row < n_block_rows(); ++row)
2220 block(row, 0).vmult(dst.block(row), src.block(0));
2221 for (size_type col = 1; col < n_block_cols(); ++col)
2222 block(row, col).vmult_add(dst.block(row), src.block(col));
2228 template <
class MatrixType>
2229 template <
class BlockVectorType,
class VectorType>
2233 const BlockVectorType &src)
const 2236 Assert(src.n_blocks() == n_block_cols(),
2239 block(0, 0).vmult(dst, src.block(0));
2240 for (size_type col = 1; col < n_block_cols(); ++col)
2241 block(0, col).vmult_add(dst, src.block(col));
2246 template <
class MatrixType>
2247 template <
class BlockVectorType,
class VectorType>
2250 const VectorType &src)
const 2252 Assert(dst.n_blocks() == n_block_rows(),
2256 for (size_type row = 0; row < n_block_rows(); ++row)
2257 block(row, 0).vmult(dst.block(row), src);
2262 template <
class MatrixType>
2263 template <
class VectorType>
2267 const VectorType &src)
const 2272 block(0, 0).vmult(dst, src);
2277 template <
class MatrixType>
2278 template <
class BlockVectorType>
2281 const BlockVectorType &src)
const 2283 Assert(dst.n_blocks() == n_block_rows(),
2285 Assert(src.n_blocks() == n_block_cols(),
2288 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2289 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2290 block(row, col).vmult_add(dst.block(row), src.block(col));
2295 template <
class MatrixType>
2296 template <
class BlockVectorType>
2299 BlockVectorType & dst,
2300 const BlockVectorType &src)
const 2302 Assert(dst.n_blocks() == n_block_cols(),
2304 Assert(src.n_blocks() == n_block_rows(),
2309 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2311 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2312 block(row, col).Tvmult_add(dst.block(col), src.block(row));
2318 template <
class MatrixType>
2319 template <
class BlockVectorType,
class VectorType>
2322 const VectorType &src)
const 2324 Assert(dst.n_blocks() == n_block_cols(),
2330 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2331 block(0, col).Tvmult_add(dst.block(col), src);
2336 template <
class MatrixType>
2337 template <
class BlockVectorType,
class VectorType>
2341 const BlockVectorType &src)
const 2344 Assert(src.n_blocks() == n_block_rows(),
2347 block(0, 0).Tvmult(dst, src.block(0));
2349 for (size_type row = 1; row < n_block_rows(); ++row)
2350 block(row, 0).Tvmult_add(dst, src.block(row));
2355 template <
class MatrixType>
2356 template <
class VectorType>
2360 const VectorType &src)
const 2365 block(0, 0).Tvmult(dst, src);
2370 template <
class MatrixType>
2371 template <
class BlockVectorType>
2374 const BlockVectorType &src)
const 2376 Assert(dst.n_blocks() == n_block_cols(),
2378 Assert(src.n_blocks() == n_block_rows(),
2381 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2382 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2383 block(row, col).Tvmult_add(dst.block(col), src.block(row));
2388 template <
class MatrixType>
2389 template <
class BlockVectorType>
2394 Assert(v.n_blocks() == n_block_rows(),
2397 value_type norm_sqr = 0;
2398 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2399 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2401 norm_sqr += block(row, col).matrix_norm_square(v.block(row));
2404 block(row, col).matrix_scalar_product(v.block(row), v.block(col));
2410 template <
class MatrixType>
2411 template <
class BlockVectorType>
2414 const BlockVectorType &u,
2415 const BlockVectorType &v)
const 2417 Assert(u.n_blocks() == n_block_rows(),
2419 Assert(v.n_blocks() == n_block_cols(),
2422 value_type result = 0;
2423 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2424 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2426 block(row, col).matrix_scalar_product(u.block(row), v.block(col));
2432 template <
class MatrixType>
2433 template <
class BlockVectorType>
2436 const BlockVectorType &x,
2437 const BlockVectorType &b)
const 2439 Assert(dst.n_blocks() == n_block_rows(),
2441 Assert(
b.n_blocks() == n_block_rows(),
2443 Assert(x.n_blocks() == n_block_cols(),
2459 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2461 block(row, 0).residual(dst.block(row), x.block(0),
b.block(row));
2463 for (size_type i = 0; i < dst.block(row).size(); ++i)
2464 dst.block(row)(i) = -dst.block(row)(i);
2466 for (
unsigned int col = 1; col < n_block_cols(); ++col)
2467 block(row, col).vmult_add(dst.block(row), x.block(col));
2469 for (size_type i = 0; i < dst.block(row).size(); ++i)
2470 dst.block(row)(i) = -dst.block(row)(i);
2474 for (size_type row = 0; row < n_block_rows(); ++row)
2475 res += dst.block(row).norm_sqr();
2476 return std::sqrt(res);
2481 template <
class MatrixType>
2484 const bool alternative_output)
const 2486 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2487 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2489 if (!alternative_output)
2490 out <<
"Block (" << row <<
", " << col <<
")" << std::endl;
2492 block(row, col).print(out, alternative_output);
2498 template <
class MatrixType>
2502 return const_iterator(
this, 0);
2507 template <
class MatrixType>
2511 return const_iterator(
this, m());
2516 template <
class MatrixType>
2521 return const_iterator(
this, r);
2526 template <
class MatrixType>
2531 return const_iterator(
this, r + 1);
2536 template <
class MatrixType>
2540 return iterator(
this, 0);
2545 template <
class MatrixType>
2549 return iterator(
this, m());
2554 template <
class MatrixType>
2559 return iterator(
this, r);
2564 template <
class MatrixType>
2569 return iterator(
this, r + 1);
2574 template <
class MatrixType>
2578 std::vector<size_type> row_sizes(this->n_block_rows());
2579 std::vector<size_type> col_sizes(this->n_block_cols());
2583 for (
unsigned int r = 0; r < this->n_block_rows(); ++r)
2584 row_sizes[r] = sub_objects[r][0]->m();
2588 for (
unsigned int c = 1; c < this->n_block_cols(); ++c)
2589 for (
unsigned int r = 0; r < this->n_block_rows(); ++r)
2590 Assert(row_sizes[r] == sub_objects[r][c]->m(),
2591 ExcIncompatibleRowNumbers(r, 0, r, c));
2595 this->row_block_indices.reinit(row_sizes);
2599 for (
unsigned int c = 0; c < this->n_block_cols(); ++c)
2600 col_sizes[c] = sub_objects[0][c]->n();
2601 for (
unsigned int r = 1; r < this->n_block_rows(); ++r)
2602 for (
unsigned int c = 0; c < this->n_block_cols(); ++c)
2603 Assert(col_sizes[c] == sub_objects[r][c]->n(),
2604 ExcIncompatibleRowNumbers(0, c, r, c));
2608 this->column_block_indices.reinit(col_sizes);
2613 template <
class MatrixType>
2617 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2618 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2619 block(row, col).prepare_add();
2624 template <
class MatrixType>
2628 for (
unsigned int row = 0; row < n_block_rows(); ++row)
2629 for (
unsigned int col = 0; col < n_block_cols(); ++col)
2630 block(row, col).prepare_set();
2636 DEAL_II_NAMESPACE_CLOSE
2638 #endif // dealii_block_matrix_base_h Iterator lower_bound(Iterator first, Iterator last, const T &val)
const types::global_dof_index invalid_size_type
static const unsigned int invalid_unsigned_int
void vmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
BlockMatrixBase & copy_from(const BlockMatrixType &source)
void Tvmult_add(BlockVectorType &dst, const BlockVectorType &src) const
static ::ExceptionBase & ExcIncompatibleColNumbers(int arg1, int arg2, int arg3, int arg4)
TemporaryData & operator=(const TemporaryData &)
Contents is actually a matrix.
std::vector< size_type > counter_within_block
std::vector< std::vector< size_type > > column_indices
void set(const size_type i, const size_type j, const value_type value)
value_type operator()(const size_type i, const size_type j) const
const BlockMatrixType MatrixType
void vmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
value_type diag_element(const size_type i) const
static ::ExceptionBase & ExcNotInitialized()
void Tvmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
void add(const size_type i, const size_type j, const value_type value)
value_type matrix_norm_square(const BlockVectorType &v) const
void Tvmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void Tvmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
TemporaryData temporary_data
static ::ExceptionBase & ExcDivideByZero()
std::size_t memory_consumption() const
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
BlockMatrixBase()=default
BlockMatrixBase & operator*=(const value_type factor)
typename BlockType::value_type value_type
const BlockMatrixType * matrix
void print(std::ostream &out, const bool alternative_output=false) const
void vmult_add(BlockVectorType &dst, const BlockVectorType &src) const
void Tvmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
const BlockIndices & get_row_indices() const
std::vector< std::vector< value_type > > column_values
BlockMatrixBase & operator/=(const value_type factor)
Table< 2, SmartPointer< BlockType, BlockMatrixBase< MatrixType > > > sub_objects
unsigned int n_block_cols() const
static ::ExceptionBase & ExcMessage(std::string arg1)
~BlockMatrixBase() override
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void prepare_add_operation()
value_type el(const size_type i, const size_type j) const
typename BlockMatrixType::value_type value_type
const BlockIndices & get_column_indices() const
BlockMatrixType::BlockType::iterator base_iterator
BlockMatrixType::BlockType::const_iterator base_iterator
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
value_type residual(BlockVectorType &dst, const BlockVectorType &x, const BlockVectorType &b) const
unsigned int block_column() const
static ::ExceptionBase & ExcIteratorPastEnd()
static ::ExceptionBase & ExcNotQuadratic()
unsigned int global_dof_index
unsigned int block_row() const
types::global_dof_index size_type
void vmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
BlockMatrixType MatrixType
BlockType & block(const unsigned int row, const unsigned int column)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
static ::ExceptionBase & ExcNotImplemented()
BlockIndices row_block_indices
unsigned int n_block_rows() const
value_type matrix_scalar_product(const BlockVectorType &u, const BlockVectorType &v) const
#define AssertIsFinite(number)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void compress(::VectorOperation::values operation)
static ::ExceptionBase & ExcIncompatibleRowNumbers(int arg1, int arg2, int arg3, int arg4)
void prepare_set_operation()
static ::ExceptionBase & ExcInternalError()