Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/block_matrix_base.h>
Classes | |
struct | TemporaryData |
Public Types | |
using | BlockType = MatrixType |
using | value_type = typename BlockType::value_type |
Public Member Functions | |
BlockMatrixBase ()=default | |
~BlockMatrixBase () override | |
template<class BlockMatrixType > | |
BlockMatrixBase & | copy_from (const BlockMatrixType &source) |
BlockType & | block (const unsigned int row, const unsigned int column) |
const BlockType & | block (const unsigned int row, const unsigned int column) const |
size_type | m () const |
size_type | n () const |
unsigned int | n_block_rows () const |
unsigned int | n_block_cols () const |
void | set (const size_type i, const size_type j, const value_type value) |
template<typename number > | |
void | set (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=false) |
template<typename number > | |
void | set (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=false) |
template<typename number > | |
void | set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=false) |
template<typename number > | |
void | set (const size_type row, const size_type n_cols, const size_type *col_indices, const number *values, const bool elide_zero_values=false) |
void | add (const size_type i, const size_type j, const value_type value) |
template<typename number > | |
void | add (const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true) |
template<typename number > | |
void | add (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true) |
template<typename number > | |
void | add (const size_type row, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=true) |
template<typename number > | |
void | add (const size_type row, const size_type n_cols, const size_type *col_indices, const number *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false) |
void | add (const value_type factor, const BlockMatrixBase< MatrixType > &matrix) |
value_type | operator() (const size_type i, const size_type j) const |
value_type | el (const size_type i, const size_type j) const |
value_type | diag_element (const size_type i) const |
void | compress (::VectorOperation::values operation) |
BlockMatrixBase & | operator*= (const value_type factor) |
BlockMatrixBase & | operator/= (const value_type factor) |
template<class BlockVectorType > | |
void | vmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
template<class BlockVectorType > | |
void | Tvmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
template<class BlockVectorType > | |
value_type | matrix_norm_square (const BlockVectorType &v) const |
template<class BlockVectorType > | |
value_type | matrix_scalar_product (const BlockVectorType &u, const BlockVectorType &v) const |
template<class BlockVectorType > | |
value_type | residual (BlockVectorType &dst, const BlockVectorType &x, const BlockVectorType &b) const |
void | print (std::ostream &out, const bool alternative_output=false) const |
iterator | begin () |
iterator | end () |
iterator | begin (const size_type r) |
iterator | end (const size_type r) |
const_iterator | begin () const |
const_iterator | end () const |
const_iterator | begin (const size_type r) const |
const_iterator | end (const size_type r) const |
const BlockIndices & | get_row_indices () const |
const BlockIndices & | get_column_indices () const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcIncompatibleRowNumbers (int arg1, int arg2, int arg3, int arg4) |
static ::ExceptionBase & | ExcIncompatibleColNumbers (int arg1, int arg2, int arg3, int arg4) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions | |
void | clear () |
void | collect_sizes () |
template<class BlockVectorType > | |
void | vmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const |
template<class BlockVectorType , class VectorType > | |
void | vmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const |
template<class BlockVectorType , class VectorType > | |
void | vmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const |
template<class VectorType > | |
void | vmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const |
template<class BlockVectorType > | |
void | Tvmult_block_block (BlockVectorType &dst, const BlockVectorType &src) const |
template<class BlockVectorType , class VectorType > | |
void | Tvmult_block_nonblock (BlockVectorType &dst, const VectorType &src) const |
template<class BlockVectorType , class VectorType > | |
void | Tvmult_nonblock_block (VectorType &dst, const BlockVectorType &src) const |
template<class VectorType > | |
void | Tvmult_nonblock_nonblock (VectorType &dst, const VectorType &src) const |
void | prepare_add_operation () |
void | prepare_set_operation () |
Protected Attributes | |
BlockIndices | row_block_indices |
Table< 2, SmartPointer< BlockType, BlockMatrixBase< MatrixType > > > | sub_objects |
Private Attributes | |
TemporaryData | temporary_data |
Friends | |
template<typename , bool > | |
class | BlockMatrixIterators::Accessor |
template<typename > | |
class | MatrixIterator |
Blocked matrix class. The behaviour of objects of this type is almost as for the usual matrix objects, with most of the functions being implemented in both classes. The main difference is that the matrix represented by this object is composed of an array of matrices (e.g. of type SparseMatrix<number>) and all accesses to the elements of this object are relayed to accesses of the base matrices. The actual type of the individual blocks of this matrix is the type of the template argument, and can, for example be the usual SparseMatrix or PETScWrappers::SparseMatrix.
In addition to the usual matrix access and linear algebra functions, there are functions block() which allow access to the different blocks of the matrix. This may, for example, be of help when you want to implement Schur complement methods, or block preconditioners, where each block belongs to a specific component of the equation you are presently discretizing.
Note that the numbers of blocks and rows are implicitly determined by the sparsity pattern objects used.
Objects of this type are frequently used when a system of differential equations has solutions with variables that fall into different classes. For example, solutions of the Stokes or Navier-Stokes equations have dim
velocity components and one pressure component. In this case, it may make sense to consider the linear system of equations as a system of 2x2 blocks, and one can construct preconditioners or solvers based on this 2x2 block structure. This class can help you in these cases, as it allows to view the matrix alternatively as one big matrix, or as a number of individual blocks.
Since this class simply forwards its calls to the subobjects (if necessary after adjusting indices denoting which subobject is meant), this class is completely independent of the actual type of the subobject. The functions that set up block matrices and destroy them, however, have to be implemented in derived classes. These functions also have to fill the data members provided by this base class, as they are only used passively in this class.
Most of the functions take a vector or block vector argument. These functions can, in general, only successfully be compiled if the individual blocks of this matrix implement the respective functions operating on the vector type in question. For example, if you have a block sparse matrix over deal.II SparseMatrix objects, then you will likely not be able to form the matrix-vector multiplication with a block vector over PETScWrappers::SparseMatrix objects. If you attempt anyway, you will likely get a number of compiler errors.
<float> and <double>
; others can be generated in application programs (see the section on Template instantiations in the manual).Definition at line 1865 of file affine_constraints.h.
using BlockMatrixBase< MatrixType >::BlockType = MatrixType |
Typedef the type of the underlying matrix.
Definition at line 360 of file block_matrix_base.h.
using BlockMatrixBase< MatrixType >::value_type = typename BlockType::value_type |
Type of matrix entries. These are analogous to alias in the standard library containers.
Definition at line 366 of file block_matrix_base.h.
|
default |
Default constructor.
|
override |
Destructor.
BlockMatrixBase& BlockMatrixBase< MatrixType >::copy_from | ( | const BlockMatrixType & | source | ) |
Copy the matrix given as argument into the current object.
Copying matrices is an expensive operation that we do not want to happen by accident through compiler generated code for operator=
. (This would happen, for example, if one accidentally declared a function argument of the current type by value rather than by reference.) The functionality of copying matrices is implemented in this member function instead. All copy operations of objects of this type therefore require an explicit function call.
The source matrix may be a matrix of arbitrary type, as long as its data type is convertible to the data type of this matrix.
The function returns a reference to this
.
BlockType& BlockMatrixBase< MatrixType >::block | ( | const unsigned int | row, |
const unsigned int | column | ||
) |
Access the block with the given coordinates.
const BlockType& BlockMatrixBase< MatrixType >::block | ( | const unsigned int | row, |
const unsigned int | column | ||
) | const |
Access the block with the given coordinates. Version for constant objects.
size_type BlockMatrixBase< MatrixType >::m | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
size_type BlockMatrixBase< MatrixType >::n | ( | ) | const |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
unsigned int BlockMatrixBase< MatrixType >::n_block_rows | ( | ) | const |
Return the number of blocks in a column. Returns zero if no sparsity pattern is presently associated to this matrix.
unsigned int BlockMatrixBase< MatrixType >::n_block_cols | ( | ) | const |
Return the number of blocks in a row. Returns zero if no sparsity pattern is presently associated to this matrix.
void BlockMatrixBase< MatrixType >::set | ( | const size_type | i, |
const size_type | j, | ||
const value_type | value | ||
) |
Set the element (i,j)
to value
. Throws an error if the entry does not exist or if value
is not a finite number. Still, it is allowed to store zero values in non-existent fields.
void BlockMatrixBase< MatrixType >::set | ( | const std::vector< size_type > & | indices, |
const FullMatrix< number > & | full_matrix, | ||
const bool | elide_zero_values = false |
||
) |
Set all elements given in a FullMatrix into the sparse matrix locations given by indices
. In other words, this function writes the elements in full_matrix
into the calling matrix, using the local-to-global indexing specified by indices
for both the rows and the columns of the matrix. This function assumes a quadratic sparse matrix and a quadratic full_matrix, the usual situation in FE calculations.
The optional parameter elide_zero_values
can be used to specify whether zero values should be set anyway or they should be filtered away (and not change the previous content in the respective element if it exists). The default value is false
, i.e., even zero values are treated.
void BlockMatrixBase< MatrixType >::set | ( | const std::vector< size_type > & | row_indices, |
const std::vector< size_type > & | col_indices, | ||
const FullMatrix< number > & | full_matrix, | ||
const bool | elide_zero_values = false |
||
) |
Same function as before, but now including the possibility to use rectangular full_matrices and different local-to-global indexing on rows and columns, respectively.
void BlockMatrixBase< MatrixType >::set | ( | const size_type | row, |
const std::vector< size_type > & | col_indices, | ||
const std::vector< number > & | values, | ||
const bool | elide_zero_values = false |
||
) |
Set several elements in the specified row of the matrix with column indices as given by col_indices
to the respective value.
The optional parameter elide_zero_values
can be used to specify whether zero values should be set anyway or they should be filtered away (and not change the previous content in the respective element if it exists). The default value is false
, i.e., even zero values are treated.
void BlockMatrixBase< MatrixType >::set | ( | const size_type | row, |
const size_type | n_cols, | ||
const size_type * | col_indices, | ||
const number * | values, | ||
const bool | elide_zero_values = false |
||
) |
Set several elements to values given by values
in a given row in columns given by col_indices into the sparse matrix.
The optional parameter elide_zero_values
can be used to specify whether zero values should be inserted anyway or they should be filtered away. The default value is false
, i.e., even zero values are inserted/replaced.
void BlockMatrixBase< MatrixType >::add | ( | const size_type | i, |
const size_type | j, | ||
const value_type | value | ||
) |
Add value
to the element (i,j). Throws an error if the entry does not exist or if value
is not a finite number. Still, it is allowed to store zero values in non-existent fields.
void BlockMatrixBase< MatrixType >::add | ( | const std::vector< size_type > & | indices, |
const FullMatrix< number > & | full_matrix, | ||
const bool | elide_zero_values = true |
||
) |
Add all elements given in a FullMatrix<double> into sparse matrix locations given by indices
. In other words, this function adds the elements in full_matrix
to the respective entries in calling matrix, using the local-to-global indexing specified by indices
for both the rows and the columns of the matrix. This function assumes a quadratic sparse matrix and a quadratic full_matrix, the usual situation in FE calculations.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
void BlockMatrixBase< MatrixType >::add | ( | const std::vector< size_type > & | row_indices, |
const std::vector< size_type > & | col_indices, | ||
const FullMatrix< number > & | full_matrix, | ||
const bool | elide_zero_values = true |
||
) |
Same function as before, but now including the possibility to use rectangular full_matrices and different local-to-global indexing on rows and columns, respectively.
void BlockMatrixBase< MatrixType >::add | ( | const size_type | row, |
const std::vector< size_type > & | col_indices, | ||
const std::vector< number > & | values, | ||
const bool | elide_zero_values = true |
||
) |
Set several elements in the specified row of the matrix with column indices as given by col_indices
to the respective value.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
void BlockMatrixBase< MatrixType >::add | ( | const size_type | row, |
const size_type | n_cols, | ||
const size_type * | col_indices, | ||
const number * | values, | ||
const bool | elide_zero_values = true , |
||
const bool | col_indices_are_sorted = false |
||
) |
Add an array of values given by values
in the given global matrix row at columns specified by col_indices in the sparse matrix.
The optional parameter elide_zero_values
can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true
, i.e., zero values won't be added into the matrix.
void BlockMatrixBase< MatrixType >::add | ( | const value_type | factor, |
const BlockMatrixBase< MatrixType > & | matrix | ||
) |
Add matrix
scaled by factor
to this matrix, i.e. the matrix factor*matrix
is added to this
. If the sparsity pattern of the calling matrix does not contain all the elements in the sparsity pattern of the input matrix, this function will throw an exception.
Depending on MatrixType, however, additional restrictions might arise. Some sparse matrix formats require matrix
to be based on the same sparsity pattern as the calling matrix.
value_type BlockMatrixBase< MatrixType >::operator() | ( | const size_type | i, |
const size_type | j | ||
) | const |
Return the value of the entry (i,j). This may be an expensive operation and you should always take care where to call this function. In order to avoid abuse, this function throws an exception if the wanted element does not exist in the matrix.
value_type BlockMatrixBase< MatrixType >::el | ( | const size_type | i, |
const size_type | j | ||
) | const |
This function is mostly like operator()() in that it returns the value of the matrix entry (i,j)
. The only difference is that if this entry does not exist in the sparsity pattern, then instead of raising an exception, zero is returned. While this may be convenient in some cases, note that it is simple to write algorithms that are slow compared to an optimal solution, since the sparsity of the matrix is not used.
value_type BlockMatrixBase< MatrixType >::diag_element | ( | const size_type | i | ) | const |
Return the main diagonal element in the ith row. This function throws an error if the matrix is not quadratic and also if the diagonal blocks of the matrix are not quadratic.
This function is considerably faster than the operator()(), since for quadratic matrices, the diagonal entry may be the first to be stored in each row and access therefore does not involve searching for the right column number.
void BlockMatrixBase< MatrixType >::compress | ( | ::VectorOperation::values | operation | ) |
Call the compress() function on all the subblocks of the matrix.
See Compressing distributed objects for more information.
BlockMatrixBase& BlockMatrixBase< MatrixType >::operator*= | ( | const value_type | factor | ) |
Multiply the entire matrix by a fixed factor.
BlockMatrixBase& BlockMatrixBase< MatrixType >::operator/= | ( | const value_type | factor | ) |
Divide the entire matrix by a fixed factor.
void BlockMatrixBase< MatrixType >::vmult_add | ( | BlockVectorType & | dst, |
const BlockVectorType & | src | ||
) | const |
Adding Matrix-vector multiplication. Add \(M*src\) on \(dst\) with \(M\) being this matrix.
void BlockMatrixBase< MatrixType >::Tvmult_add | ( | BlockVectorType & | dst, |
const BlockVectorType & | src | ||
) | const |
Adding Matrix-vector multiplication. Add MTsrc to dst with M being this matrix. This function does the same as vmult_add() but takes the transposed matrix.
value_type BlockMatrixBase< MatrixType >::matrix_norm_square | ( | const BlockVectorType & | v | ) | const |
Return the norm of the vector v with respect to the norm induced by this matrix, i.e. vTMv). This is useful, e.g. in the finite element context, where the LT-norm of a function equals the matrix norm with respect to the mass matrix of the vector representing the nodal values of the finite element function. Note that even though the function's name might suggest something different, for historic reasons not the norm but its square is returned, as defined above by the scalar product.
Obviously, the matrix needs to be square for this operation.
value_type BlockMatrixBase< MatrixType >::matrix_scalar_product | ( | const BlockVectorType & | u, |
const BlockVectorType & | v | ||
) | const |
Compute the matrix scalar product \(\left(u,Mv\right)\).
value_type BlockMatrixBase< MatrixType >::residual | ( | BlockVectorType & | dst, |
const BlockVectorType & | x, | ||
const BlockVectorType & | b | ||
) | const |
Compute the residual r=b-Ax. Write the residual into dst
.
void BlockMatrixBase< MatrixType >::print | ( | std::ostream & | out, |
const bool | alternative_output = false |
||
) | const |
Print the matrix to the given stream, using the format (line,col) value
, i.e. one nonzero entry of the matrix per line. The optional flag outputs the sparsity pattern in a different style according to the underlying sparse matrix type.
iterator BlockMatrixBase< MatrixType >::begin | ( | ) |
Iterator starting at the first entry.
iterator BlockMatrixBase< MatrixType >::end | ( | ) |
Final iterator.
iterator BlockMatrixBase< MatrixType >::begin | ( | const size_type | r | ) |
Iterator starting at the first entry of row r
.
iterator BlockMatrixBase< MatrixType >::end | ( | const size_type | r | ) |
Final iterator of row r
.
const_iterator BlockMatrixBase< MatrixType >::begin | ( | ) | const |
Iterator starting at the first entry.
const_iterator BlockMatrixBase< MatrixType >::end | ( | ) | const |
Final iterator.
const_iterator BlockMatrixBase< MatrixType >::begin | ( | const size_type | r | ) | const |
Iterator starting at the first entry of row r
.
const_iterator BlockMatrixBase< MatrixType >::end | ( | const size_type | r | ) | const |
Final iterator of row r
.
const BlockIndices& BlockMatrixBase< MatrixType >::get_row_indices | ( | ) | const |
Return a reference to the underlying BlockIndices data of the rows.
const BlockIndices& BlockMatrixBase< MatrixType >::get_column_indices | ( | ) | const |
Return a reference to the underlying BlockIndices data of the columns.
std::size_t BlockMatrixBase< MatrixType >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object. Note that only the memory reserved on the current processor is returned in case this is called in an MPI-based program.
|
protected |
Release all memory and return to a state just like after having called the default constructor. It also forgets the sparsity pattern it was previously tied to.
This calls clear for all sub-matrices and then resets this object to have no blocks at all.
This function is protected since it may be necessary to release additional structures. A derived class can make it public again, if it is sufficient.
|
protected |
This function collects the sizes of the sub-objects and stores them in internal arrays, in order to be able to relay global indices into the matrix to indices into the subobjects. You must call this function each time after you have changed the size of the sub-objects.
Derived classes should call this function whenever the size of the sub- objects has changed and the X_block_indices
arrays need to be updated.
Note that this function is not public since not all derived classes need to export its interface. For example, for the usual deal.II SparseMatrix class, the sizes are implicitly determined whenever reinit() is called, and individual blocks cannot be resized. For that class, this function therefore does not have to be public. On the other hand, for the PETSc classes, there is no associated sparsity pattern object that determines the block sizes, and for these the function needs to be publicly available. These classes therefore export this function.
|
protected |
Matrix-vector multiplication: let \(dst = M*src\) with \(M\) being this matrix.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block column.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block row.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Matrix-vector multiplication: let \(dst = M^T*src\) with \(M\) being this matrix. This function does the same as vmult() but takes the transposed matrix.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block row.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block column.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Matrix-vector multiplication. Just like the previous function, but only applicable if the matrix has only one block.
Due to problems with deriving template arguments between the block and non-block versions of the vmult/Tvmult functions, the actual functions are implemented in derived classes, with implementations forwarding the calls to the implementations provided here under a unique name for which template arguments can be derived by the compiler.
|
protected |
Some matrix types, in particular PETSc, need to synchronize set and add operations. This has to be done for all matrices in the BlockMatrix. This routine prepares adding of elements by notifying all blocks. Called by all internal routines before adding elements.
|
protected |
Notifies all blocks to let them prepare for setting elements, see prepare_add_operation().
|
friend |
Make the iterator class a friend. We have to work around a compiler bug here again.
Definition at line 1068 of file block_matrix_base.h.
|
protected |
Index arrays for rows and columns.
Definition at line 843 of file block_matrix_base.h.
|
protected |
Array of sub-matrices.
Definition at line 849 of file block_matrix_base.h.
|
private |
A set of scratch arrays that can be used by the add() and set() functions that take pointers to data to pre-sort indices before use. Access from multiple threads is synchronized via the mutex variable that is part of the structure.
Definition at line 1061 of file block_matrix_base.h.