16 #ifndef dealii__block_matrix_base_h 17 #define dealii__block_matrix_base_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/table.h> 22 #include <deal.II/base/thread_management.h> 23 #include <deal.II/base/utilities.h> 24 #include <deal.II/base/smartpointer.h> 25 #include <deal.II/base/memory_consumption.h> 26 #include <deal.II/lac/block_indices.h> 27 #include <deal.II/lac/exceptions.h> 28 #include <deal.II/lac/full_matrix.h> 29 #include <deal.II/lac/matrix_iterator.h> 30 #include <deal.II/lac/vector.h> 34 DEAL_II_NAMESPACE_OPEN
56 template <
class BlockMatrixType>
108 template <
class BlockMatrixType,
bool Constness>
115 template <
class BlockMatrixType>
187 bool operator == (
const Accessor &a)
const;
190 friend class Accessor<BlockMatrixType, true>;
197 template <
class BlockMatrixType>
268 bool operator == (
const Accessor &a)
const;
274 friend class ::MatrixIterator;
340 template <
typename MatrixType>
395 template <
class BlockMatrixType>
397 copy_from (
const BlockMatrixType &source);
403 block (
const unsigned int row,
404 const unsigned int column);
412 block (
const unsigned int row,
413 const unsigned int column)
const;
419 size_type
m ()
const;
425 size_type
n ()
const;
445 void set (
const size_type i,
464 template <
typename number>
465 void set (
const std::vector<size_type> &indices,
467 const bool elide_zero_values =
false);
474 template <
typename number>
475 void set (
const std::vector<size_type> &row_indices,
476 const std::vector<size_type> &col_indices,
478 const bool elide_zero_values =
false);
490 template <
typename number>
491 void set (
const size_type row,
492 const std::vector<size_type> &col_indices,
493 const std::vector<number> &values,
494 const bool elide_zero_values =
false);
505 template <
typename number>
506 void set (
const size_type row,
507 const size_type n_cols,
508 const size_type *col_indices,
509 const number *values,
510 const bool elide_zero_values =
false);
517 void add (
const size_type i,
535 template <
typename number>
536 void add (
const std::vector<size_type> &indices,
538 const bool elide_zero_values =
true);
545 template <
typename number>
546 void add (
const std::vector<size_type> &row_indices,
547 const std::vector<size_type> &col_indices,
549 const bool elide_zero_values =
true);
560 template <
typename number>
561 void add (
const size_type row,
562 const std::vector<size_type> &col_indices,
563 const std::vector<number> &values,
564 const bool elide_zero_values =
true);
575 template <
typename number>
576 void add (
const size_type row,
577 const size_type n_cols,
578 const size_type *col_indices,
579 const number *values,
580 const bool elide_zero_values =
true,
581 const bool col_indices_are_sorted =
false);
604 const size_type j)
const;
615 const size_type j)
const;
653 template <
class BlockVectorType>
655 const BlockVectorType &src)
const;
662 template <
class BlockVectorType>
664 const BlockVectorType &src)
const;
678 template <
class BlockVectorType>
685 template <
class BlockVectorType>
688 const BlockVectorType &v)
const;
693 template <
class BlockVectorType>
695 const BlockVectorType &x,
696 const BlockVectorType &b)
const;
704 void print (std::ostream &out,
705 const bool alternative_output =
false)
const;
773 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" 774 << arg3 <<
',' << arg4 <<
"] have differing row numbers.");
780 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" 781 << arg3 <<
',' << arg4 <<
"] have differing column numbers.");
839 template <
class BlockVectorType>
841 const BlockVectorType &src)
const;
853 template <
class BlockVectorType,
856 const VectorType &src)
const;
868 template <
class BlockVectorType,
871 const BlockVectorType &src)
const;
883 template <
class VectorType>
885 const VectorType &src)
const;
898 template <
class BlockVectorType>
900 const BlockVectorType &src)
const;
912 template <
class BlockVectorType,
915 const VectorType &src)
const;
927 template <
class BlockVectorType,
930 const BlockVectorType &src)
const;
942 template <
class VectorType>
944 const VectorType &src)
const;
1027 template <
typename,
bool>
1043 template <
class BlockMatrixType>
1052 template <
class BlockMatrixType>
1055 AccessorBase<BlockMatrixType>::block_row()
const 1064 template <
class BlockMatrixType>
1067 AccessorBase<BlockMatrixType>::block_column()
const 1076 template <
class BlockMatrixType>
1078 Accessor<BlockMatrixType, true>::Accessor (
1079 const BlockMatrixType *matrix,
1080 const size_type row,
1081 const size_type col)
1084 base_iterator(matrix->block(0,0).begin())
1091 if (row < matrix->
m())
1093 const std::pair<unsigned int,size_type> indices
1094 = matrix->row_block_indices.global_to_local(row);
1098 for (
unsigned int bc=0; bc<matrix->n_block_cols(); ++bc)
1101 = matrix->block(indices.first, bc).begin(indices.second);
1102 if (base_iterator !=
1103 matrix->block(indices.first, bc).end(indices.second))
1105 this->row_block = indices.first;
1106 this->col_block = bc;
1117 *
this = Accessor (matrix, row+1, 0);
1141 template <
class BlockMatrixType>
1143 Accessor<BlockMatrixType, true>::Accessor (
const Accessor<BlockMatrixType, false> &other)
1145 matrix(other.matrix),
1146 base_iterator(other.base_iterator)
1148 this->row_block = other.row_block;
1149 this->col_block = other.col_block;
1153 template <
class BlockMatrixType>
1155 typename Accessor<BlockMatrixType, true>::size_type
1156 Accessor<BlockMatrixType, true>::row()
const 1161 return (matrix->row_block_indices.local_to_global(this->row_block, 0) +
1162 base_iterator->row());
1166 template <
class BlockMatrixType>
1168 typename Accessor<BlockMatrixType, true>::size_type
1169 Accessor<BlockMatrixType, true>::column()
const 1174 return (matrix->column_block_indices.local_to_global(this->col_block,0) +
1175 base_iterator->column());
1179 template <
class BlockMatrixType>
1181 typename Accessor<BlockMatrixType, true>::value_type
1182 Accessor<BlockMatrixType, true>::value ()
const 1189 return base_iterator->value();
1194 template <
class BlockMatrixType>
1197 Accessor<BlockMatrixType, true>::advance ()
1205 size_type local_row = base_iterator->row();
1216 while (base_iterator ==
1217 matrix->block(this->row_block, this->col_block).end(local_row))
1222 if (this->col_block < matrix->n_block_cols()-1)
1226 = matrix->block(this->row_block, this->col_block).begin(local_row);
1232 this->col_block = 0;
1241 if (local_row == matrix->block(this->row_block, this->col_block).m())
1245 if (this->row_block == matrix->n_block_rows())
1254 = matrix->block(this->row_block, this->col_block).begin(local_row);
1260 template <
class BlockMatrixType>
1263 Accessor<BlockMatrixType, true>::operator == (
const Accessor &a)
const 1265 if (matrix != a.matrix)
1268 if (this->row_block == a.row_block
1269 && this->col_block == a.col_block)
1278 (base_iterator == a.base_iterator));
1286 template <
class BlockMatrixType>
1288 Accessor<BlockMatrixType, false>::Accessor (
1289 BlockMatrixType *matrix,
1290 const size_type row,
1291 const size_type col)
1294 base_iterator(matrix->block(0,0).begin())
1300 if (row < matrix->
m())
1302 const std::pair<unsigned int,size_type> indices
1303 = matrix->row_block_indices.global_to_local(row);
1307 for (size_type bc=0; bc<matrix->n_block_cols(); ++bc)
1310 = matrix->block(indices.first, bc).begin(indices.second);
1311 if (base_iterator !=
1312 matrix->block(indices.first, bc).end(indices.second))
1314 this->row_block = indices.first;
1315 this->col_block = bc;
1326 *
this = Accessor (matrix, row+1, 0);
1338 template <
class BlockMatrixType>
1340 typename Accessor<BlockMatrixType, false>::size_type
1341 Accessor<BlockMatrixType, false>::row()
const 1346 return (matrix->row_block_indices.local_to_global(this->row_block, 0) +
1347 base_iterator->row());
1351 template <
class BlockMatrixType>
1353 typename Accessor<BlockMatrixType, false>::size_type
1354 Accessor<BlockMatrixType, false>::column()
const 1359 return (matrix->column_block_indices.local_to_global(this->col_block,0) +
1360 base_iterator->column());
1364 template <
class BlockMatrixType>
1366 typename Accessor<BlockMatrixType, false>::value_type
1367 Accessor<BlockMatrixType, false>::value ()
const 1374 return base_iterator->value();
1379 template <
class BlockMatrixType>
1382 Accessor<BlockMatrixType, false>::set_value (
typename Accessor<BlockMatrixType, false>::value_type newval)
const 1389 base_iterator->value() = newval;
1394 template <
class BlockMatrixType>
1397 Accessor<BlockMatrixType, false>::advance ()
1405 size_type local_row = base_iterator->row();
1416 while (base_iterator ==
1417 matrix->block(this->row_block, this->col_block).end(local_row))
1422 if (this->col_block < matrix->n_block_cols()-1)
1426 = matrix->block(this->row_block, this->col_block).begin(local_row);
1432 this->col_block = 0;
1441 if (local_row == matrix->block(this->row_block, this->col_block).m())
1445 if (this->row_block == matrix->n_block_rows())
1454 = matrix->block(this->row_block, this->col_block).begin(local_row);
1461 template <
class BlockMatrixType>
1464 Accessor<BlockMatrixType, false>::operator == (
const Accessor &a)
const 1466 if (matrix != a.matrix)
1469 if (this->row_block == a.row_block
1470 && this->col_block == a.col_block)
1479 (base_iterator == a.base_iterator));
1489 template <
typename MatrixType>
1494 template <
typename MatrixType>
1507 template <
class MatrixType>
1508 template <
class BlockMatrixType>
1512 copy_from (
const BlockMatrixType &source)
1514 for (
unsigned int r=0; r<n_block_rows(); ++r)
1515 for (
unsigned int c=0; c<n_block_cols(); ++c)
1516 block(r,c).
copy_from (source.block(r,c));
1522 template <
class MatrixType>
1533 sizeof(temporary_data.mutex);
1535 for (
unsigned int r=0; r<n_block_rows(); ++r)
1536 for (
unsigned int c=0; c<n_block_cols(); ++c)
1538 MatrixType *p = this->sub_objects[r][c];
1547 template <
class MatrixType>
1552 for (
unsigned int r=0; r<n_block_rows(); ++r)
1553 for (
unsigned int c=0; c<n_block_cols(); ++c)
1555 MatrixType *p = this->sub_objects[r][c];
1556 this->sub_objects[r][c] = 0;
1559 sub_objects.reinit (0,0);
1562 row_block_indices = column_block_indices =
BlockIndices ();
1567 template <
class MatrixType>
1571 const unsigned int column)
1573 Assert (row<n_block_rows(),
1575 Assert (column<n_block_cols(),
1578 return *sub_objects[row][column];
1583 template <
class MatrixType>
1587 const unsigned int column)
const 1589 Assert (row<n_block_rows(),
1591 Assert (column<n_block_cols(),
1594 return *sub_objects[row][column];
1598 template <
class MatrixType>
1600 typename BlockMatrixBase<MatrixType>::size_type
1603 return row_block_indices.total_size();
1608 template <
class MatrixType>
1610 typename BlockMatrixBase<MatrixType>::size_type
1613 return column_block_indices.total_size();
1618 template <
class MatrixType>
1623 return column_block_indices.size();
1628 template <
class MatrixType>
1633 return row_block_indices.size();
1641 template <
class MatrixType>
1646 const value_type value)
1648 prepare_set_operation();
1652 const std::pair<unsigned int,size_type>
1653 row_index = row_block_indices.global_to_local (i),
1654 col_index = column_block_indices.global_to_local (j);
1655 block(row_index.first,col_index.first).set (row_index.second,
1662 template <
class MatrixType>
1663 template <
typename number>
1667 const std::vector<size_type> &col_indices,
1669 const bool elide_zero_values)
1671 Assert (row_indices.size() == values.
m(),
1673 Assert (col_indices.size() == values.
n(),
1676 for (size_type i=0; i<row_indices.size(); ++i)
1677 set (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1683 template <
class MatrixType>
1684 template <
typename number>
1689 const bool elide_zero_values)
1691 Assert (indices.size() == values.
m(),
1695 for (size_type i=0; i<indices.size(); ++i)
1696 set (indices[i], indices.size(), &indices[0], &values(i,0),
1702 template <
class MatrixType>
1703 template <
typename number>
1707 const std::vector<size_type> &col_indices,
1708 const std::vector<number> &values,
1709 const bool elide_zero_values)
1711 Assert (col_indices.size() == values.size(),
1714 set (row, col_indices.size(), &col_indices[0], &values[0],
1723 template <
class MatrixType>
1724 template <
typename number>
1728 const size_type n_cols,
1729 const size_type *col_indices,
1730 const number *values,
1731 const bool elide_zero_values)
1733 prepare_set_operation();
1740 if (temporary_data.column_indices.size() < this->n_block_cols())
1742 temporary_data.column_indices.resize (this->n_block_cols());
1743 temporary_data.column_values.resize (this->n_block_cols());
1744 temporary_data.counter_within_block.resize (this->n_block_cols());
1757 if (temporary_data.column_indices[0].size() < n_cols)
1759 for (
unsigned int i=0; i<this->n_block_cols(); ++i)
1761 temporary_data.column_indices[i].resize(n_cols);
1762 temporary_data.column_values[i].resize(n_cols);
1768 for (
unsigned int i=0; i<this->n_block_cols(); ++i)
1769 temporary_data.counter_within_block[i] = 0;
1780 for (size_type j=0; j<n_cols; ++j)
1782 number value = values[j];
1784 if (value == number() && elide_zero_values ==
true)
1787 const std::pair<unsigned int, size_type>
1788 col_index = this->column_block_indices.global_to_local(col_indices[j]);
1790 const size_type local_index = temporary_data.counter_within_block[col_index.first]++;
1792 temporary_data.column_indices[col_index.first][local_index] = col_index.second;
1793 temporary_data.column_values[col_index.first][local_index] = value;
1800 for (
unsigned int i=0; i<this->n_block_cols(); ++i)
1801 length += temporary_data.counter_within_block[i];
1810 const std::pair<unsigned int,size_type>
1811 row_index = this->row_block_indices.global_to_local (row);
1812 for (
unsigned int block_col=0; block_col<n_block_cols(); ++block_col)
1814 if (temporary_data.counter_within_block[block_col] == 0)
1817 block(row_index.first, block_col).set
1819 temporary_data.counter_within_block[block_col],
1820 &temporary_data.column_indices[block_col][0],
1821 &temporary_data.column_values[block_col][0],
1828 template <
class MatrixType>
1833 const value_type value)
1838 prepare_add_operation();
1843 typedef typename MatrixType::Traits MatrixTraits;
1844 if ((MatrixTraits::zero_addition_can_be_elided ==
true)
1846 (value == value_type()))
1849 const std::pair<unsigned int,size_type>
1850 row_index = row_block_indices.global_to_local (i),
1851 col_index = column_block_indices.global_to_local (j);
1852 block(row_index.first,col_index.first).add (row_index.second,
1859 template <
class MatrixType>
1860 template <
typename number>
1864 const std::vector<size_type> &col_indices,
1866 const bool elide_zero_values)
1868 Assert (row_indices.size() == values.
m(),
1870 Assert (col_indices.size() == values.
n(),
1873 for (size_type i=0; i<row_indices.size(); ++i)
1874 add (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1880 template <
class MatrixType>
1881 template <
typename number>
1886 const bool elide_zero_values)
1888 Assert (indices.size() == values.
m(),
1892 for (size_type i=0; i<indices.size(); ++i)
1893 add (indices[i], indices.size(), &indices[0], &values(i,0),
1899 template <
class MatrixType>
1900 template <
typename number>
1904 const std::vector<size_type> &col_indices,
1905 const std::vector<number> &values,
1906 const bool elide_zero_values)
1908 Assert (col_indices.size() == values.size(),
1911 add (row, col_indices.size(), &col_indices[0], &values[0],
1920 template <
class MatrixType>
1921 template <
typename number>
1925 const size_type n_cols,
1926 const size_type *col_indices,
1927 const number *values,
1928 const bool elide_zero_values,
1929 const bool col_indices_are_sorted)
1931 prepare_add_operation();
1936 if (col_indices_are_sorted ==
true)
1942 for (size_type i=1; i<n_cols; ++i)
1943 if (col_indices[i] <= before)
1945 "indices appear to not be sorted."))
1947 before = col_indices[i];
1949 const std::pair<unsigned int,size_type>
1950 row_index = this->row_block_indices.global_to_local (row);
1952 if (this->n_block_cols() > 1)
1956 this->column_block_indices.block_start(1));
1958 const size_type n_zero_block_indices = first_block - col_indices;
1959 block(row_index.first, 0).add (row_index.second,
1960 n_zero_block_indices,
1964 col_indices_are_sorted);
1966 if (n_zero_block_indices < n_cols)
1967 this->add(row, n_cols - n_zero_block_indices, first_block,
1968 values + n_zero_block_indices, elide_zero_values,
1973 block(row_index.first, 0). add (row_index.second,
1978 col_indices_are_sorted);
1987 if (temporary_data.column_indices.size() < this->n_block_cols())
1989 temporary_data.column_indices.resize (this->n_block_cols());
1990 temporary_data.column_values.resize (this->n_block_cols());
1991 temporary_data.counter_within_block.resize (this->n_block_cols());
2004 if (temporary_data.column_indices[0].size() < n_cols)
2006 for (
unsigned int i=0; i<this->n_block_cols(); ++i)
2008 temporary_data.column_indices[i].resize(n_cols);
2009 temporary_data.column_values[i].resize(n_cols);
2015 for (
unsigned int i=0; i<this->n_block_cols(); ++i)
2016 temporary_data.counter_within_block[i] = 0;
2027 for (size_type j=0; j<n_cols; ++j)
2029 number value = values[j];
2031 if (value == number() && elide_zero_values ==
true)
2034 const std::pair<unsigned int, size_type>
2035 col_index = this->column_block_indices.global_to_local(col_indices[j]);
2037 const size_type local_index = temporary_data.counter_within_block[col_index.first]++;
2039 temporary_data.column_indices[col_index.first][local_index] = col_index.second;
2040 temporary_data.column_values[col_index.first][local_index] = value;
2047 for (
unsigned int i=0; i<this->n_block_cols(); ++i)
2048 length += temporary_data.counter_within_block[i];
2057 const std::pair<unsigned int,size_type>
2058 row_index = this->row_block_indices.global_to_local (row);
2059 for (
unsigned int block_col=0; block_col<n_block_cols(); ++block_col)
2061 if (temporary_data.counter_within_block[block_col] == 0)
2064 block(row_index.first, block_col).add
2066 temporary_data.counter_within_block[block_col],
2067 &temporary_data.column_indices[block_col][0],
2068 &temporary_data.column_values[block_col][0],
2070 col_indices_are_sorted);
2076 template <
class MatrixType>
2084 prepare_add_operation();
2089 typedef typename MatrixType::Traits MatrixTraits;
2090 if ((MatrixTraits::zero_addition_can_be_elided ==
true)
2095 for (
unsigned int row=0; row<n_block_rows(); ++row)
2096 for (
unsigned int col=0; col<n_block_cols(); ++col)
2099 block(row, col).add(factor, matrix.block(row,col));
2104 template <
class MatrixType>
2108 const size_type j)
const 2110 const std::pair<unsigned int,size_type>
2111 row_index = row_block_indices.global_to_local (i),
2112 col_index = column_block_indices.global_to_local (j);
2113 return block(row_index.first,col_index.first) (row_index.second,
2119 template <
class MatrixType>
2123 const size_type j)
const 2125 const std::pair<unsigned int,size_type>
2126 row_index = row_block_indices.global_to_local (i),
2127 col_index = column_block_indices.global_to_local (j);
2128 return block(row_index.first,col_index.first).el (row_index.second,
2134 template <
class MatrixType>
2139 Assert (n_block_rows() == n_block_cols(),
2142 const std::pair<unsigned int,size_type>
2143 index = row_block_indices.global_to_local (i);
2144 return block(index.first,index.first).diag_element(index.second);
2149 template <
class MatrixType>
2154 for (
unsigned int r=0; r<n_block_rows(); ++r)
2155 for (
unsigned int c=0; c<n_block_cols(); ++c)
2156 block(r,c).compress (operation);
2161 template <
class MatrixType>
2169 for (
unsigned int r=0; r<n_block_rows(); ++r)
2170 for (
unsigned int c=0; c<n_block_cols(); ++c)
2171 block(r,c) *= factor;
2178 template <
class MatrixType>
2187 const value_type factor_inv = 1. / factor;
2189 for (
unsigned int r=0; r<n_block_rows(); ++r)
2190 for (
unsigned int c=0; c<n_block_cols(); ++c)
2191 block(r,c) *= factor_inv;
2198 template <
class MatrixType>
2202 return this->row_block_indices;
2207 template <
class MatrixType>
2211 return this->column_block_indices;
2216 template <
class MatrixType>
2217 template <
class BlockVectorType>
2221 const BlockVectorType &src)
const 2223 Assert (dst.n_blocks() == n_block_rows(),
2225 Assert (src.n_blocks() == n_block_cols(),
2228 for (size_type row=0; row<n_block_rows(); ++row)
2230 block(row,0).vmult (dst.block(row),
2232 for (size_type col=1; col<n_block_cols(); ++col)
2233 block(row,col).vmult_add (dst.block(row),
2240 template <
class MatrixType>
2241 template <
class BlockVectorType,
2246 const BlockVectorType &src)
const 2248 Assert (n_block_rows() == 1,
2250 Assert (src.n_blocks() == n_block_cols(),
2253 block(0,0).vmult (dst, src.block(0));
2254 for (size_type col=1; col<n_block_cols(); ++col)
2255 block(0,col).vmult_add (dst, src.block(col));
2260 template <
class MatrixType>
2261 template <
class BlockVectorType,
2266 const VectorType &src)
const 2268 Assert (dst.n_blocks() == n_block_rows(),
2270 Assert (1 == n_block_cols(),
2273 for (size_type row=0; row<n_block_rows(); ++row)
2274 block(row,0).vmult (dst.block(row),
2280 template <
class MatrixType>
2281 template <
class VectorType>
2285 const VectorType &src)
const 2287 Assert (1 == n_block_rows(),
2289 Assert (1 == n_block_cols(),
2292 block(0,0).vmult (dst, src);
2297 template <
class MatrixType>
2298 template <
class BlockVectorType>
2301 const BlockVectorType &src)
const 2303 Assert (dst.n_blocks() == n_block_rows(),
2305 Assert (src.n_blocks() == n_block_cols(),
2308 for (
unsigned int row=0; row<n_block_rows(); ++row)
2309 for (
unsigned int col=0; col<n_block_cols(); ++col)
2310 block(row,col).vmult_add (dst.block(row),
2317 template <
class MatrixType>
2318 template <
class BlockVectorType>
2322 const BlockVectorType &src)
const 2324 Assert (dst.n_blocks() == n_block_cols(),
2326 Assert (src.n_blocks() == n_block_rows(),
2331 for (
unsigned int row=0; row<n_block_rows(); ++row)
2333 for (
unsigned int col=0; col<n_block_cols(); ++col)
2334 block(row,col).Tvmult_add (dst.block(col),
2341 template <
class MatrixType>
2342 template <
class BlockVectorType,
2347 const VectorType &src)
const 2349 Assert (dst.n_blocks() == n_block_cols(),
2351 Assert (1 == n_block_rows(),
2356 for (
unsigned int col=0; col<n_block_cols(); ++col)
2357 block(0,col).Tvmult_add (dst.block(col), src);
2362 template <
class MatrixType>
2363 template <
class BlockVectorType,
2368 const BlockVectorType &src)
const 2370 Assert (1 == n_block_cols(),
2372 Assert (src.n_blocks() == n_block_rows(),
2375 block(0,0).Tvmult (dst, src.block(0));
2377 for (size_type row=1; row<n_block_rows(); ++row)
2378 block(row,0).Tvmult_add (dst, src.block(row));
2383 template <
class MatrixType>
2384 template <
class VectorType>
2388 const VectorType &src)
const 2390 Assert (1 == n_block_cols(),
2392 Assert (1 == n_block_rows(),
2395 block(0,0).Tvmult (dst, src);
2400 template <
class MatrixType>
2401 template <
class BlockVectorType>
2404 const BlockVectorType &src)
const 2406 Assert (dst.n_blocks() == n_block_cols(),
2408 Assert (src.n_blocks() == n_block_rows(),
2411 for (
unsigned int row=0; row<n_block_rows(); ++row)
2412 for (
unsigned int col=0; col<n_block_cols(); ++col)
2413 block(row,col).Tvmult_add (dst.block(col),
2419 template <
class MatrixType>
2420 template <
class BlockVectorType>
2425 Assert (v.n_blocks() == n_block_rows(),
2428 value_type norm_sqr = 0;
2429 for (
unsigned int row=0; row<n_block_rows(); ++row)
2430 for (
unsigned int col=0; col<n_block_cols(); ++col)
2432 norm_sqr += block(row,col).matrix_norm_square (v.block(row));
2434 norm_sqr += block(row,col).matrix_scalar_product (v.block(row),
2441 template <
class MatrixType>
2442 template <
class BlockVectorType>
2446 const BlockVectorType &v)
const 2448 Assert (u.n_blocks() == n_block_rows(),
2450 Assert (v.n_blocks() == n_block_cols(),
2453 value_type result = 0;
2454 for (
unsigned int row=0; row<n_block_rows(); ++row)
2455 for (
unsigned int col=0; col<n_block_cols(); ++col)
2456 result += block(row,col).matrix_scalar_product (u.block(row),
2463 template <
class MatrixType>
2464 template <
class BlockVectorType>
2468 const BlockVectorType &x,
2469 const BlockVectorType &b)
const 2471 Assert (dst.n_blocks() == n_block_rows(),
2473 Assert (
b.n_blocks() == n_block_rows(),
2475 Assert (x.n_blocks() == n_block_cols(),
2491 for (
unsigned int row=0; row<n_block_rows(); ++row)
2493 block(row,0).residual (dst.block(row),
2497 for (size_type i=0; i<dst.block(row).size(); ++i)
2498 dst.block(row)(i) = -dst.block(row)(i);
2500 for (
unsigned int col=1; col<n_block_cols(); ++col)
2501 block(row,col).vmult_add (dst.block(row),
2504 for (size_type i=0; i<dst.block(row).size(); ++i)
2505 dst.block(row)(i) = -dst.block(row)(i);
2509 for (size_type row=0; row<n_block_rows(); ++row)
2510 res += dst.block(row).norm_sqr ();
2511 return std::sqrt(res);
2516 template <
class MatrixType>
2520 const bool alternative_output)
const 2522 for (
unsigned int row=0; row<n_block_rows(); ++row)
2523 for (
unsigned int col=0; col<n_block_cols(); ++col)
2525 if (!alternative_output)
2526 out <<
"Block (" << row <<
", " << col <<
")" << std::endl;
2528 block(row, col).print(out, alternative_output);
2534 template <
class MatrixType>
2539 return const_iterator(
this, 0);
2544 template <
class MatrixType>
2549 return const_iterator(
this, m());
2554 template <
class MatrixType>
2560 return const_iterator(
this, r);
2565 template <
class MatrixType>
2571 return const_iterator(
this, r+1);
2576 template <
class MatrixType>
2581 return iterator(
this, 0);
2586 template <
class MatrixType>
2591 return iterator(
this, m());
2596 template <
class MatrixType>
2602 return iterator(
this, r);
2607 template <
class MatrixType>
2613 return iterator(
this, r+1);
2618 template <
class MatrixType>
2622 std::vector<size_type> row_sizes (this->n_block_rows());
2623 std::vector<size_type> col_sizes (this->n_block_cols());
2627 for (
unsigned int r=0; r<this->n_block_rows(); ++r)
2628 row_sizes[r] = sub_objects[r][0]->m();
2632 for (
unsigned int c=1; c<this->n_block_cols(); ++c)
2633 for (
unsigned int r=0; r<this->n_block_rows(); ++r)
2634 Assert (row_sizes[r] == sub_objects[r][c]->m(),
2635 ExcIncompatibleRowNumbers (r,0,r,c));
2639 this->row_block_indices.reinit (row_sizes);
2643 for (
unsigned int c=0; c<this->n_block_cols(); ++c)
2644 col_sizes[c] = sub_objects[0][c]->n();
2645 for (
unsigned int r=1; r<this->n_block_rows(); ++r)
2646 for (
unsigned int c=0; c<this->n_block_cols(); ++c)
2647 Assert (col_sizes[c] == sub_objects[r][c]->n(),
2648 ExcIncompatibleRowNumbers (0,c,r,c));
2652 this->column_block_indices.reinit (col_sizes);
2657 template <
class MatrixType>
2661 for (
unsigned int row=0; row<n_block_rows(); ++row)
2662 for (
unsigned int col=0; col<n_block_cols(); ++col)
2663 block(row, col).prepare_add();
2668 template <
class MatrixType>
2672 for (
unsigned int row=0; row<n_block_rows(); ++row)
2673 for (
unsigned int col=0; col<n_block_cols(); ++col)
2674 block(row, col).prepare_set();
2680 DEAL_II_NAMESPACE_CLOSE
2682 #endif // dealii__block_matrix_base_h BlockMatrixType MatrixType
Iterator lower_bound(Iterator first, Iterator last, const T &val)
types::global_dof_index size_type
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 &)
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
BlockMatrixType::value_type value_type
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
const BlockMatrixType MatrixType
TemporaryData temporary_data
static ::ExceptionBase & ExcDivideByZero()
types::global_dof_index size_type
std::size_t memory_consumption() const
BlockMatrixType::value_type value_type
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
BlockMatrixBase & operator*=(const value_type factor)
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
types::global_dof_index size_type
BlockMatrixBase & operator/=(const value_type factor)
unsigned int n_block_cols() const
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
types::global_dof_index size_type
void prepare_add_operation()
value_type el(const size_type i, const size_type j) const
const BlockIndices & get_column_indices() const
BlockType::value_type value_type
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()
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
unsigned int block_row() const
BlockMatrixType::value_type value_type
void vmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
BlockType & block(const unsigned int row, const unsigned int column)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Table< 2, SmartPointer< BlockType, BlockMatrixBase< MatrixType > > > sub_objects
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)
void compress(::VectorOperation::values operation)
static ::ExceptionBase & ExcIncompatibleRowNumbers(int arg1, int arg2, int arg3, int arg4)
void prepare_set_operation()
static ::ExceptionBase & ExcInternalError()