Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/full_matrix.h>
Public Types | |
using | size_type = std::size_t |
using | value_type = number |
using | iterator = typename Table< 2, number >::iterator |
using | const_iterator = typename Table< 2, number >::const_iterator |
using | real_type = typename numbers::NumberTraits< number >::real_type |
Public Types inherited from TableBase< N, T > | |
using | size_type = typename AlignedVector< T >::size_type |
Public Member Functions | |
Constructors and initialization. See also the base class Table. | |
FullMatrix (const size_type n=0) | |
FullMatrix (const size_type rows, const size_type cols) | |
FullMatrix (const size_type rows, const size_type cols, const number *entries) | |
FullMatrix (const IdentityMatrix &id) | |
Copying into and out of other matrices | |
template<typename number2 > | |
FullMatrix< number > & | operator= (const FullMatrix< number2 > &) |
FullMatrix< number > & | operator= (const number d) |
FullMatrix< number > & | operator= (const IdentityMatrix &id) |
template<typename number2 > | |
FullMatrix< number > & | operator= (const LAPACKFullMatrix< number2 > &) |
template<typename MatrixType > | |
void | copy_from (const MatrixType &) |
template<typename MatrixType > | |
void | copy_transposed (const MatrixType &) |
template<int dim> | |
void | copy_from (const Tensor< 2, dim > &T, const unsigned int src_r_i=0, const unsigned int src_r_j=dim - 1, const unsigned int src_c_i=0, const unsigned int src_c_j=dim - 1, const size_type dst_r=0, const size_type dst_c=0) |
template<int dim> | |
void | copy_to (Tensor< 2, dim > &T, const size_type src_r_i=0, const size_type src_r_j=dim - 1, const size_type src_c_i=0, const size_type src_c_j=dim - 1, const unsigned int dst_r=0, const unsigned int dst_c=0) const |
template<typename MatrixType , typename index_type > | |
void | extract_submatrix_from (const MatrixType &matrix, const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set) |
template<typename MatrixType , typename index_type > | |
void | scatter_matrix_to (const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set, MatrixType &matrix) const |
template<typename number2 > | |
void | fill (const FullMatrix< number2 > &src, const size_type dst_offset_i=0, const size_type dst_offset_j=0, const size_type src_offset_i=0, const size_type src_offset_j=0) |
template<typename number2 > | |
void | fill (const number2 *) |
template<typename number2 > | |
void | fill_permutation (const FullMatrix< number2 > &src, const std::vector< size_type > &p_rows, const std::vector< size_type > &p_cols) |
void | set (const size_type i, const size_type j, const number value) |
Non-modifying operators | |
bool | operator== (const FullMatrix< number > &) const |
size_type | m () const |
size_type | n () const |
bool | all_zero () const |
template<typename number2 > | |
number2 | matrix_norm_square (const Vector< number2 > &v) const |
template<typename number2 > | |
number2 | matrix_scalar_product (const Vector< number2 > &u, const Vector< number2 > &v) const |
real_type | l1_norm () const |
real_type | linfty_norm () const |
real_type | frobenius_norm () const |
real_type | relative_symmetry_norm2 () const |
number | determinant () const |
number | trace () const |
template<class StreamType > | |
void | print (StreamType &s, const unsigned int width=5, const unsigned int precision=2) const |
void | print_formatted (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1., const double threshold=0.) const |
std::size_t | memory_consumption () const |
Iterator functions | |
iterator | begin (const size_type r) |
iterator | end (const size_type r) |
const_iterator | begin (const size_type r) const |
const_iterator | end (const size_type r) const |
Modifying operators | |
FullMatrix & | operator*= (const number factor) |
FullMatrix & | operator/= (const number factor) |
template<typename number2 > | |
void | add (const number a, const FullMatrix< number2 > &A) |
template<typename number2 > | |
void | add (const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B) |
template<typename number2 > | |
void | add (const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B, const number c, const FullMatrix< number2 > &C) |
template<typename number2 > | |
void | add (const FullMatrix< number2 > &src, const number factor, const size_type dst_offset_i=0, const size_type dst_offset_j=0, const size_type src_offset_i=0, const size_type src_offset_j=0) |
template<typename number2 > | |
void | Tadd (const number s, const FullMatrix< number2 > &B) |
template<typename number2 > | |
void | Tadd (const FullMatrix< number2 > &src, const number factor, const size_type dst_offset_i=0, const size_type dst_offset_j=0, const size_type src_offset_i=0, const size_type src_offset_j=0) |
void | add (const size_type row, const size_type column, const number value) |
template<typename number2 , typename index_type > | |
void | add (const size_type row, const size_type n_cols, const index_type *col_indices, const number2 *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false) |
void | add_row (const size_type i, const number s, const size_type j) |
void | add_row (const size_type i, const number s, const size_type j, const number t, const size_type k) |
void | add_col (const size_type i, const number s, const size_type j) |
void | add_col (const size_type i, const number s, const size_type j, const number t, const size_type k) |
void | swap_row (const size_type i, const size_type j) |
void | swap_col (const size_type i, const size_type j) |
void | diagadd (const number s) |
template<typename number2 > | |
void | equ (const number a, const FullMatrix< number2 > &A) |
template<typename number2 > | |
void | equ (const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B) |
template<typename number2 > | |
void | equ (const number a, const FullMatrix< number2 > &A, const number b, const FullMatrix< number2 > &B, const number c, const FullMatrix< number2 > &C) |
void | symmetrize () |
void | gauss_jordan () |
template<typename number2 > | |
void | invert (const FullMatrix< number2 > &M) |
template<typename number2 > | |
void | cholesky (const FullMatrix< number2 > &A) |
template<typename number2 > | |
void | outer_product (const Vector< number2 > &V, const Vector< number2 > &W) |
template<typename number2 > | |
void | left_invert (const FullMatrix< number2 > &M) |
template<typename number2 > | |
void | right_invert (const FullMatrix< number2 > &M) |
Multiplications | |
template<typename number2 > | |
void | mmult (FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const |
template<typename number2 > | |
void | Tmmult (FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const |
template<typename number2 > | |
void | mTmult (FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const |
template<typename number2 > | |
void | TmTmult (FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const |
void | triple_product (const FullMatrix< number > &A, const FullMatrix< number > &B, const FullMatrix< number > &D, const bool transpose_B=false, const bool transpose_D=false, const number scaling=number(1.)) |
template<typename number2 > | |
void | vmult (Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const |
template<typename number2 > | |
void | vmult_add (Vector< number2 > &w, const Vector< number2 > &v) const |
template<typename number2 > | |
void | Tvmult (Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const |
template<typename number2 > | |
void | Tvmult_add (Vector< number2 > &w, const Vector< number2 > &v) const |
template<typename somenumber > | |
void | precondition_Jacobi (Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const |
template<typename number2 , typename number3 > | |
number | residual (Vector< number2 > &dst, const Vector< number2 > &x, const Vector< number3 > &b) const |
template<typename number2 > | |
void | forward (Vector< number2 > &dst, const Vector< number2 > &src) const |
template<typename number2 > | |
void | backward (Vector< number2 > &dst, const Vector< number2 > &src) const |
Public Member Functions inherited from TableBase< N, T > | |
TableBase ()=default | |
TableBase (const TableIndices< N > &sizes) | |
template<typename InputIterator > | |
TableBase (const TableIndices< N > &sizes, InputIterator entries, const bool C_style_indexing=true) | |
TableBase (const TableBase< N, T > &src) | |
template<typename T2 > | |
TableBase (const TableBase< N, T2 > &src) | |
TableBase (TableBase< N, T > &&src) noexcept | |
~TableBase () override=default | |
TableBase< N, T > & | operator= (const TableBase< N, T > &src) |
template<typename T2 > | |
TableBase< N, T > & | operator= (const TableBase< N, T2 > &src) |
TableBase< N, T > & | operator= (TableBase< N, T > &&src) noexcept |
bool | operator== (const TableBase< N, T > &T2) const |
void | reset_values () |
void | reinit (const TableIndices< N > &new_size, const bool omit_default_initialization=false) |
size_type | size (const unsigned int i) const |
const TableIndices< N > & | size () const |
size_type | n_elements () const |
bool | empty () const |
template<typename InputIterator > | |
void | fill (InputIterator entries, const bool C_style_indexing=true) |
void | fill (const T &value) |
AlignedVector< T >::reference | operator() (const TableIndices< N > &indices) |
AlignedVector< T >::const_reference | operator() (const TableIndices< N > &indices) const |
void | swap (TableBase< N, T > &v) |
std::size_t | memory_consumption () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
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 & | ExcEmptyMatrix () |
static ::ExceptionBase & | ExcNotRegular (number arg1) |
static ::ExceptionBase & | ExcInvalidDestination (size_type arg1, size_type arg2, size_type arg3) |
static ::ExceptionBase & | ExcSourceEqualsDestination () |
static ::ExceptionBase & | ExcMatrixNotPositiveDefinite () |
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) |
Additional Inherited Members | |
Protected Member Functions inherited from TableBase< N, T > | |
size_type | position (const TableIndices< N > &indices) const |
AlignedVector< T >::reference | el (const TableIndices< N > &indices) |
AlignedVector< T >::const_reference | el (const TableIndices< N > &indices) const |
Protected Attributes inherited from TableBase< N, T > | |
AlignedVector< T > | values |
TableIndices< N > | table_size |
Implementation of a classical rectangular scheme of numbers. The data type of the entries is provided in the template argument number
. The interface is quite fat and in fact has grown every time a new feature was needed. So, a lot of functions are provided.
Internal calculations are usually done with the accuracy of the vector argument to functions. If there is no argument with a number type, the matrix number type is used.
<float>, <double>, <std::complex<float>>, <std::complex<double>>
. Others can be generated in application programs, see Template instantiations for details.Definition at line 33 of file dof_accessor.h.
using FullMatrix< number >::size_type = std::size_t |
A type of used to index into this container.
Definition at line 82 of file full_matrix.h.
using FullMatrix< number >::value_type = number |
Type of matrix entries. This alias is analogous to value_type
in the standard library containers.
Definition at line 88 of file full_matrix.h.
using FullMatrix< number >::iterator = typename Table<2, number>::iterator |
Use the base class mutable iterator type.
Definition at line 93 of file full_matrix.h.
using FullMatrix< number >::const_iterator = typename Table<2, number>::const_iterator |
Use the base class constant iterator type.
Definition at line 98 of file full_matrix.h.
using FullMatrix< number >::real_type = typename numbers::NumberTraits<number>::real_type |
Declare a type that has holds real-valued numbers with the same precision as the template argument to this class. If the template argument of this class is a real data type, then real_type equals the template argument. If the template argument is a std::complex type then real_type equals the type underlying the complex numbers.
This alias is used to represent the return type of norms.
Definition at line 119 of file full_matrix.h.
|
explicit |
Constructor. Initialize the matrix as a square matrix with dimension n
.
In order to avoid the implicit conversion of integers and other types to a matrix, this constructor is declared explicit
.
By default, no memory is allocated.
FullMatrix< number >::FullMatrix | ( | const size_type | rows, |
const size_type | cols | ||
) |
Constructor. Initialize the matrix as a rectangular matrix.
FullMatrix< number >::FullMatrix | ( | const size_type | rows, |
const size_type | cols, | ||
const number * | entries | ||
) |
Constructor initializing from an array of numbers. The array is arranged line by line. No range checking is performed.
FullMatrix< number >::FullMatrix | ( | const IdentityMatrix & | id | ) |
Construct a full matrix that equals the identity matrix of the size of the argument. Using this constructor, one can easily create an identity matrix of size n
by saying
FullMatrix<number>& FullMatrix< number >::operator= | ( | const FullMatrix< number2 > & | ) |
Variable assignment operator.
FullMatrix<number>& FullMatrix< number >::operator= | ( | const number | d | ) |
This operator assigns a scalar to a matrix. To avoid confusion with the semantics of this function, zero is the only value allowed for d
, allowing you to clear a matrix in an intuitive way.
FullMatrix<number>& FullMatrix< number >::operator= | ( | const IdentityMatrix & | id | ) |
Copy operator to create a full matrix that equals the identity matrix of the size of the argument. This way, one can easily create an identity matrix of size n
by saying
FullMatrix<number>& FullMatrix< number >::operator= | ( | const LAPACKFullMatrix< number2 > & | ) |
Assignment operator for a LapackFullMatrix. The calling matrix must be of the same size as the LAPACK matrix.
void FullMatrix< number >::copy_from | ( | const MatrixType & | ) |
Assignment from different matrix classes. This assignment operator uses iterators of the typename MatrixType. Therefore, sparse matrices are possible sources.
void FullMatrix< number >::copy_transposed | ( | const MatrixType & | ) |
Transposing assignment from different matrix classes. This assignment operator uses iterators of the typename MatrixType. Therefore, sparse matrices are possible sources.
void FullMatrix< number >::copy_from | ( | const Tensor< 2, dim > & | T, |
const unsigned int | src_r_i = 0 , |
||
const unsigned int | src_r_j = dim - 1 , |
||
const unsigned int | src_c_i = 0 , |
||
const unsigned int | src_c_j = dim - 1 , |
||
const size_type | dst_r = 0 , |
||
const size_type | dst_c = 0 |
||
) |
Fill matrix with elements extracted from a tensor, taking rows included between r_i
and r_j
and columns between c_i
and c_j
. The resulting matrix is then inserted in the destination matrix at position (dst_r, dst_c)
Checks on the indices are made.
void FullMatrix< number >::copy_to | ( | Tensor< 2, dim > & | T, |
const size_type | src_r_i = 0 , |
||
const size_type | src_r_j = dim - 1 , |
||
const size_type | src_c_i = 0 , |
||
const size_type | src_c_j = dim - 1 , |
||
const unsigned int | dst_r = 0 , |
||
const unsigned int | dst_c = 0 |
||
) | const |
Insert a submatrix (also rectangular) into a tensor, putting its upper left element at the specified position (dst_r, dst_c)
and the other elements consequently. Default values are chosen so that no parameter needs to be specified if the size of the tensor and that of the matrix coincide.
void FullMatrix< number >::extract_submatrix_from | ( | const MatrixType & | matrix, |
const std::vector< index_type > & | row_index_set, | ||
const std::vector< index_type > & | column_index_set | ||
) |
Copy a subset of the rows and columns of another matrix into the current object.
matrix | The matrix from which a subset is to be taken from. |
row_index_set | The set of rows of matrix from which to extract. |
column_index_set | The set of columns of matrix from which to extract. |
row_index_set
and column_index_set
shall be equal to the number of rows and columns in the current object. In other words, the current object is not resized for this operation. void FullMatrix< number >::scatter_matrix_to | ( | const std::vector< index_type > & | row_index_set, |
const std::vector< index_type > & | column_index_set, | ||
MatrixType & | matrix | ||
) | const |
Copy the elements of the current matrix object into a specified set of rows and columns of another matrix. Thus, this is a scatter operation.
row_index_set | The rows of matrix into which to write. |
column_index_set | The columns of matrix into which to write. |
matrix | The matrix within which certain elements are to be replaced. |
row_index_set
and column_index_set
shall be equal to the number of rows and columns in the current object. In other words, the current object is not resized for this operation. void FullMatrix< number >::fill | ( | const FullMatrix< number2 > & | src, |
const size_type | dst_offset_i = 0 , |
||
const size_type | dst_offset_j = 0 , |
||
const size_type | src_offset_i = 0 , |
||
const size_type | src_offset_j = 0 |
||
) |
Fill rectangular block.
A rectangular block of the matrix src
is copied into this
. The upper left corner of the block being copied is (src_offset_i,src_offset_j)
. The upper left corner of the copied block is (dst_offset_i,dst_offset_j)
. The size of the rectangular block being copied is the maximum size possible, determined either by the size of this
or src
.
void FullMatrix< number >::fill | ( | const number2 * | ) |
Make function of base class available.
void FullMatrix< number >::fill_permutation | ( | const FullMatrix< number2 > & | src, |
const std::vector< size_type > & | p_rows, | ||
const std::vector< size_type > & | p_cols | ||
) |
Fill with permutation of another matrix.
The matrix src
is copied into the target. The two permutation p_r
and p_c
operate in a way, such that result(i,j) = src(p_r[i], p_c[j])
.
The vectors may also be a selection from a larger set of integers, if the matrix src
is bigger. It is also possible to duplicate rows or columns by this method.
void FullMatrix< number >::set | ( | const size_type | i, |
const size_type | j, | ||
const number | value | ||
) |
Set a particular entry of the matrix to a value. Thus, calling A.set(1,2,3.141);
is entirely equivalent to the operation A(1,2) = 3.141;
. This function exists for compatibility with the various sparse matrix objects.
i | The row index of the element to be set. |
j | The columns index of the element to be set. |
value | The value to be written into the element. |
bool FullMatrix< number >::operator== | ( | const FullMatrix< number > & | ) | const |
Comparison operator. Be careful with this thing, it may eat up huge amounts of computing time! It is most commonly used for internal consistency checks of programs.
size_type FullMatrix< number >::m | ( | ) | const |
Number of rows of this matrix. Note that the matrix is of dimension m x n.
size_type FullMatrix< number >::n | ( | ) | const |
Number of columns of this matrix. Note that the matrix is of dimension m x n.
bool FullMatrix< number >::all_zero | ( | ) | const |
Return whether the matrix contains only elements with value zero. This function is mainly for internal consistency checks and should seldom be used when not in debug mode since it uses quite some time.
number2 FullMatrix< number >::matrix_norm_square | ( | const Vector< number2 > & | v | ) | const |
Return the square of the norm of the vector v
induced by this matrix, i.e. (v,Mv). This is useful, e.g. in the finite element context, where the L2 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.
Obviously, the matrix needs to be quadratic for this operation, and for the result to actually be a norm it also needs to be either real symmetric or complex hermitian.
The underlying template types of both this matrix and the given vector should either both be real or complex-valued, but not mixed, for this function to make sense.
number2 FullMatrix< number >::matrix_scalar_product | ( | const Vector< number2 > & | u, |
const Vector< number2 > & | v | ||
) | const |
Build the matrix scalar product uT M v
. This function is mostly useful when building the cellwise scalar product of two functions in the finite element context.
The underlying template types of both this matrix and the given vector should either both be real or complex-valued, but not mixed, for this function to make sense.
real_type FullMatrix< number >::l1_norm | ( | ) | const |
Return the l1-norm of the matrix, where \(||M||_1 = \max_j \sum_i |M_{ij}|\) (maximum of the sums over columns).
real_type FullMatrix< number >::linfty_norm | ( | ) | const |
Return the \(l_\infty\)-norm of the matrix, where \(||M||_\infty = \max_i \sum_j |M_{ij}|\) (maximum of the sums over rows).
real_type FullMatrix< number >::frobenius_norm | ( | ) | const |
Compute the Frobenius norm of the matrix. Return value is the root of the square sum of all matrix entries.
real_type FullMatrix< number >::relative_symmetry_norm2 | ( | ) | const |
Compute the relative norm of the skew-symmetric part. The return value is the Frobenius norm of the skew-symmetric part of the matrix divided by that of the matrix.
Main purpose of this function is to check, if a matrix is symmetric within a certain accuracy, or not.
number FullMatrix< number >::determinant | ( | ) | const |
Compute the determinant of a matrix. This is only implemented for one, two, and three dimensions, since for higher dimensions the numerical work explodes. Obviously, the matrix needs to be quadratic for this function.
number FullMatrix< number >::trace | ( | ) | const |
Return the trace of the matrix, i.e. the sum of the diagonal values (which happens to also equal the sum of the eigenvalues of a matrix). Obviously, the matrix needs to be quadratic for this function.
void FullMatrix< number >::print | ( | StreamType & | s, |
const unsigned int | width = 5 , |
||
const unsigned int | precision = 2 |
||
) | const |
Output of the matrix in user-defined format given by the specified precision and width. This function saves width and precision of the stream before setting these given values for output, and restores the previous values after output.
void FullMatrix< number >::print_formatted | ( | std::ostream & | out, |
const unsigned int | precision = 3 , |
||
const bool | scientific = true , |
||
const unsigned int | width = 0 , |
||
const char * | zero_string = " " , |
||
const double | denominator = 1. , |
||
const double | threshold = 0. |
||
) | const |
Print the matrix and allow formatting of entries.
The parameters allow for a flexible setting of the output format:
precision
denotes the number of trailing digits.scientific
is used to determine the number format, where scientific
= false
means fixed point notation.width
denotes the with of each column. A zero entry for width
makes the function compute a width, but it may be changed to a positive value, if output is crude.zero_string
specifies a string printed for zero entries.denominator
Multiply the whole matrix by this common denominator to get nicer numbers.threshold
: all entries with absolute value smaller than this are considered zero. std::size_t FullMatrix< number >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
iterator FullMatrix< number >::begin | ( | const size_type | r | ) |
Mutable iterator starting at the first entry of row r
.
iterator FullMatrix< number >::end | ( | const size_type | r | ) |
One past the end mutable iterator of row r
.
const_iterator FullMatrix< number >::begin | ( | const size_type | r | ) | const |
Constant iterator starting at the first entry of row r
.
const_iterator FullMatrix< number >::end | ( | const size_type | r | ) | const |
One past the end constant iterator of row r
.
FullMatrix& FullMatrix< number >::operator*= | ( | const number | factor | ) |
Scale the entire matrix by a fixed factor.
FullMatrix& FullMatrix< number >::operator/= | ( | const number | factor | ) |
Scale the entire matrix by the inverse of the given factor.
void FullMatrix< number >::add | ( | const number | a, |
const FullMatrix< number2 > & | A | ||
) |
Simple addition of a scaled matrix, i.e. *this += a*A
.
The matrix A
may be a full matrix over an arbitrary underlying scalar type, as long as its data type is convertible to the data type of this matrix.
void FullMatrix< number >::add | ( | const number | a, |
const FullMatrix< number2 > & | A, | ||
const number | b, | ||
const FullMatrix< number2 > & | B | ||
) |
Multiple addition of scaled matrices, i.e. *this += a*A + b*B
.
The matrices A
and B
may be a full matrix over an arbitrary underlying scalar type, as long as its data type is convertible to the data type of this matrix.
void FullMatrix< number >::add | ( | const number | a, |
const FullMatrix< number2 > & | A, | ||
const number | b, | ||
const FullMatrix< number2 > & | B, | ||
const number | c, | ||
const FullMatrix< number2 > & | C | ||
) |
Multiple addition of scaled matrices, i.e. *this += a*A + b*B + c*C
.
The matrices A
, B
and C
may be a full matrix over an arbitrary underlying scalar type, as long as its data type is convertible to the data type of this matrix.
void FullMatrix< number >::add | ( | const FullMatrix< number2 > & | src, |
const number | factor, | ||
const size_type | dst_offset_i = 0 , |
||
const size_type | dst_offset_j = 0 , |
||
const size_type | src_offset_i = 0 , |
||
const size_type | src_offset_j = 0 |
||
) |
Add rectangular block.
A rectangular block of the matrix src
is added to this
. The upper left corner of the block being copied is (src_offset_i,src_offset_j)
. The upper left corner of the copied block is (dst_offset_i,dst_offset_j)
. The size of the rectangular block being copied is the maximum size possible, determined either by the size of this
or src
and the given offsets.
void FullMatrix< number >::Tadd | ( | const number | s, |
const FullMatrix< number2 > & | B | ||
) |
Weighted addition of the transpose of B
to this
.
A += s BT
void FullMatrix< number >::Tadd | ( | const FullMatrix< number2 > & | src, |
const number | factor, | ||
const size_type | dst_offset_i = 0 , |
||
const size_type | dst_offset_j = 0 , |
||
const size_type | src_offset_i = 0 , |
||
const size_type | src_offset_j = 0 |
||
) |
Add transpose of a rectangular block.
A rectangular block of the matrix src
is transposed and addedadded to this
. The upper left corner of the block being copied is (src_offset_i,src_offset_j)
in the coordinates of the non-transposed matrix. The upper left corner of the copied block is (dst_offset_i,dst_offset_j)
. The size of the rectangular block being copied is the maximum size possible, determined either by the size of this
or src
.
void FullMatrix< number >::add | ( | const size_type | row, |
const size_type | column, | ||
const number | value | ||
) |
Add a single element at the given position.
void FullMatrix< number >::add | ( | const size_type | row, |
const size_type | n_cols, | ||
const index_type * | col_indices, | ||
const number2 * | 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 full matrix. This function is present for compatibility with the various sparse matrices in deal.II. In particular, the two boolean fields elide_zero_values
and col_indices_are_sorted
do not impact the performance of this routine, as opposed to the sparse matrix case and are indeed ignored in the implementation.
void FullMatrix< number >::add_row | ( | const size_type | i, |
const number | s, | ||
const size_type | j | ||
) |
A(i,1...n) += s*A(j,1...n). Simple addition of rows of this
void FullMatrix< number >::add_row | ( | const size_type | i, |
const number | s, | ||
const size_type | j, | ||
const number | t, | ||
const size_type | k | ||
) |
A(i,1...n) += s*A(j,1...n) + t*A(k,1...n). Multiple addition of rows of this.
void FullMatrix< number >::add_col | ( | const size_type | i, |
const number | s, | ||
const size_type | j | ||
) |
A(1...n,i) += s*A(1...n,j). Simple addition of columns of this.
void FullMatrix< number >::add_col | ( | const size_type | i, |
const number | s, | ||
const size_type | j, | ||
const number | t, | ||
const size_type | k | ||
) |
A(1...n,i) += s*A(1...n,j) + t*A(1...n,k). Multiple addition of columns of this.
void FullMatrix< number >::swap_row | ( | const size_type | i, |
const size_type | j | ||
) |
Swap A(i,1...n) <-> A(j,1...n). Swap rows i and j of this
void FullMatrix< number >::swap_col | ( | const size_type | i, |
const size_type | j | ||
) |
Swap A(1...n,i) <-> A(1...n,j). Swap columns i and j of this
void FullMatrix< number >::diagadd | ( | const number | s | ) |
Add constant to diagonal elements of this, i.e. add a multiple of the identity matrix.
void FullMatrix< number >::equ | ( | const number | a, |
const FullMatrix< number2 > & | A | ||
) |
Assignment *this = a*A
.
void FullMatrix< number >::equ | ( | const number | a, |
const FullMatrix< number2 > & | A, | ||
const number | b, | ||
const FullMatrix< number2 > & | B | ||
) |
Assignment *this = a*A + b*B
.
void FullMatrix< number >::equ | ( | const number | a, |
const FullMatrix< number2 > & | A, | ||
const number | b, | ||
const FullMatrix< number2 > & | B, | ||
const number | c, | ||
const FullMatrix< number2 > & | C | ||
) |
Assignment *this = a*A + b*B + c*C
.
void FullMatrix< number >::symmetrize | ( | ) |
Symmetrize the matrix by forming the mean value between the existing matrix and its transpose, A = 1/2(A+AT).
Obviously the matrix must be quadratic for this operation.
void FullMatrix< number >::gauss_jordan | ( | ) |
A=Inverse(A). A must be a square matrix. Inversion of this matrix by Gauss-Jordan algorithm with partial pivoting. This process is well- behaved for positive definite matrices, but be aware of round-off errors in the indefinite case.
In case deal.II was configured with LAPACK, the functions Xgetrf and Xgetri build an LU factorization and invert the matrix upon that factorization, providing best performance up to matrices with a few hundreds rows and columns.
The numerical effort to invert an n x n
matrix is of the order n**3
.
template void FullMatrix< number >::invert< long double > | ( | const FullMatrix< number2 > & | M | ) |
Assign the inverse of the given matrix to *this
. This function is hardcoded for quadratic matrices of dimension one to four. However, since the amount of code needed grows quickly, the method gauss_jordan() is invoked implicitly if the dimension is larger.
void FullMatrix< number >::cholesky | ( | const FullMatrix< number2 > & | A | ) |
Assign the Cholesky decomposition \(A=:L L^T\) of the given matrix \(A\) to *this
, where \(L\) is lower triangular matrix. The given matrix must be symmetric positive definite.
ExcMatrixNotPositiveDefinite will be thrown in the case that the matrix is not positive definite.
void FullMatrix< number >::outer_product | ( | const Vector< number2 > & | V, |
const Vector< number2 > & | W | ||
) |
*this(i,j)
= \(V(i) W(j)\) where \(V,W\) are vectors of the same length.
void FullMatrix< number >::left_invert | ( | const FullMatrix< number2 > & | M | ) |
Assign the left_inverse of the given matrix to *this
. The calculation being performed is (AT*A)-1 *AT.
void FullMatrix< number >::right_invert | ( | const FullMatrix< number2 > & | M | ) |
Assign the right_inverse of the given matrix to *this
. The calculation being performed is AT*(A*AT) -1.
template void FullMatrix< number >::mmult< long double > | ( | FullMatrix< number2 > & | C, |
const FullMatrix< number2 > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += A*B
if (!adding) C = A*B
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm if the product of the three matrix dimensions is larger than 300 and BLAS was detected during configuration. Using BLAS usually results in considerable performance gains.
template void FullMatrix< number >::Tmmult< long double > | ( | FullMatrix< number2 > & | C, |
const FullMatrix< number2 > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication using transpose of this
.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += AT*B
if (!adding) C = AT*B
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm if the product of the three matrix dimensions is larger than 300 and BLAS was detected during configuration. Using BLAS usually results in considerable performance gains.
template void FullMatrix< number >::mTmult< long double > | ( | FullMatrix< number2 > & | C, |
const FullMatrix< number2 > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication using transpose of B
.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += A*BT
if (!adding) C = A*BT
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm if the product of the three matrix dimensions is larger than 300 and BLAS was detected during configuration. Using BLAS usually results in considerable performance gains.
template void FullMatrix< number >::TmTmult< long double > | ( | FullMatrix< number2 > & | C, |
const FullMatrix< number2 > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication using transpose of this
and B
.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += AT*BT
if (!adding) C = AT*BT
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm if the product of the three matrix dimensions is larger than 300 and BLAS was detected during configuration. Using BLAS usually results in considerable performance gains.
void FullMatrix< number >::triple_product | ( | const FullMatrix< number > & | A, |
const FullMatrix< number > & | B, | ||
const FullMatrix< number > & | D, | ||
const bool | transpose_B = false , |
||
const bool | transpose_D = false , |
||
const number | scaling = number(1.) |
||
) |
Add to the current matrix the triple product B A D. Optionally, use the transposes of the matrices B and D. The scaling factor scales the whole product, which is helpful when adding a multiple of the triple product to the matrix.
This product was written with the Schur complement BT A-1 D in mind. Note that in this case the argument for A
must be the inverse of the matrix A.
template void FullMatrix< number >::vmult< long double > | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v, | ||
const bool | adding = false |
||
) | const |
Matrix-vector-multiplication.
The optional parameter adding
determines, whether the result is stored in w
or added to w
.
if (adding) w += A*v
if (!adding) w = A*v
Source and destination must not be the same vector.
void FullMatrix< number >::vmult_add | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v | ||
) | const |
Adding Matrix-vector-multiplication. w += A*v
Source and destination must not be the same vector.
template void FullMatrix< number >::Tvmult< long double > | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v, | ||
const bool | adding = false |
||
) | const |
Transpose matrix-vector-multiplication.
The optional parameter adding
determines, whether the result is stored in w
or added to w
.
if (adding) w += AT*v
if (!adding) w = AT*v
Source and destination must not be the same vector.
void FullMatrix< number >::Tvmult_add | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v | ||
) | const |
Adding transpose matrix-vector-multiplication. w += AT*v
Source and destination must not be the same vector.
void FullMatrix< number >::precondition_Jacobi | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src, | ||
const number | omega = 1. |
||
) | const |
Apply the Jacobi preconditioner, which multiplies every element of the src
vector by the inverse of the respective diagonal element and multiplies the result with the damping factor omega
.
number FullMatrix< number >::residual | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | x, | ||
const Vector< number3 > & | b | ||
) | const |
dst=b-A*x. Residual calculation, returns the l2-norm |dst|.
Source x and destination dst must not be the same vector.
void FullMatrix< number >::forward | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | src | ||
) | const |
Forward elimination of lower triangle. Inverts the lower triangle of a rectangular matrix for a given right hand side.
If the matrix has more columns than rows, this function only operates on the left quadratic submatrix. If there are more rows, the upper quadratic part of the matrix is considered.
dst
and src
. void FullMatrix< number >::backward | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | src | ||
) | const |
Backward elimination of upper triangle.
See forward()
dst
and src
.