Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/lapack_full_matrix.h>
Public Types | |
using | size_type = std::make_unsigned_t<types::blas_int> |
using | value_type |
using | reference |
using | const_reference |
using | const_iterator |
using | iterator |
Public Member Functions | |
LAPACKFullMatrix (const size_type size=0) | |
LAPACKFullMatrix (const size_type rows, const size_type cols) | |
LAPACKFullMatrix (const LAPACKFullMatrix &) | |
LAPACKFullMatrix< number > & | operator= (const LAPACKFullMatrix< number > &) |
template<typename number2 > | |
LAPACKFullMatrix< number > & | operator= (const FullMatrix< number2 > &) |
template<typename number2 > | |
LAPACKFullMatrix< number > & | operator= (const SparseMatrix< number2 > &) |
LAPACKFullMatrix< number > & | operator= (const number d) |
LAPACKFullMatrix< number > & | operator*= (const number factor) |
LAPACKFullMatrix< number > & | operator/= (const number factor) |
void | set (const size_type i, const size_type j, const number value) |
void | add (const number a, const LAPACKFullMatrix< number > &B) |
void | rank1_update (const number a, const Vector< number > &v) |
void | apply_givens_rotation (const std::array< number, 3 > &csr, const size_type i, const size_type k, const bool left=true) |
template<typename MatrixType > | |
void | copy_from (const MatrixType &) |
void | reinit (const size_type size) |
void | grow_or_shrink (const size_type size) |
void | remove_row_and_column (const size_type row, const size_type col) |
void | reinit (const size_type rows, const size_type cols) |
void | set_property (const LAPACKSupport::Property property) |
size_type | m () const |
size_type | n () const |
template<typename MatrixType > | |
void | fill (const MatrixType &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, const number factor=1., const bool transpose=false) |
template<typename number2 > | |
void | vmult (Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const |
void | vmult (Vector< number > &w, const Vector< number > &v, const bool adding=false) const |
template<typename number2 > | |
void | vmult_add (Vector< number2 > &w, const Vector< number2 > &v) const |
void | vmult_add (Vector< number > &w, const Vector< number > &v) const |
template<typename number2 > | |
void | Tvmult (Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const |
void | Tvmult (Vector< number > &w, const Vector< number > &v, const bool adding=false) const |
template<typename number2 > | |
void | Tvmult_add (Vector< number2 > &w, const Vector< number2 > &v) const |
void | Tvmult_add (Vector< number > &w, const Vector< number > &v) const |
void | mmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | mmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | Tmmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | Tmmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | Tmmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const Vector< number > &V, const bool adding=false) const |
void | mTmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | mTmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | TmTmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | TmTmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | transpose (LAPACKFullMatrix< number > &B) const |
void | scale_rows (const Vector< number > &V) |
void | compute_lu_factorization () |
void | compute_cholesky_factorization () |
number | reciprocal_condition_number (const number l1_norm) const |
number | reciprocal_condition_number () const |
number | determinant () const |
number | l1_norm () const |
number | linfty_norm () const |
number | frobenius_norm () const |
number | trace () const |
void | invert () |
void | solve (Vector< number > &v, const bool transposed=false) const |
void | solve (LAPACKFullMatrix< number > &B, const bool transposed=false) const |
void | compute_eigenvalues (const bool right_eigenvectors=false, const bool left_eigenvectors=false) |
void | compute_eigenvalues_symmetric (const number lower_bound, const number upper_bound, const number abs_accuracy, Vector< number > &eigenvalues, FullMatrix< number > &eigenvectors) |
void | compute_generalized_eigenvalues_symmetric (LAPACKFullMatrix< number > &B, const number lower_bound, const number upper_bound, const number abs_accuracy, Vector< number > &eigenvalues, std::vector< Vector< number > > &eigenvectors, const types::blas_int itype=1) |
void | compute_generalized_eigenvalues_symmetric (LAPACKFullMatrix< number > &B, std::vector< Vector< number > > &eigenvectors, const types::blas_int itype=1) |
void | compute_svd () |
void | compute_inverse_svd (const double threshold=0.) |
void | compute_inverse_svd_with_kernel (const unsigned int kernel_size) |
std::complex< typename numbers::NumberTraits< number >::real_type > | eigenvalue (const size_type i) const |
FullMatrix< std::complex< typename numbers::NumberTraits< number >::real_type > > | get_right_eigenvectors () const |
FullMatrix< std::complex< typename numbers::NumberTraits< number >::real_type > > | get_left_eigenvectors () const |
number | singular_value (const size_type i) const |
const LAPACKFullMatrix< number > & | get_svd_u () const |
const LAPACKFullMatrix< number > & | get_svd_vt () 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 char *separator=" ") const |
void | reinit (const size_type size1, const size_type size2, const bool omit_default_initialization=false) |
void | reinit (const TableIndices< N > &new_size, const bool omit_default_initialization=false) |
const_reference | operator() (const size_type i, const size_type j) const |
reference | operator() (const size_type i, const size_type j) |
AlignedVector< number >::reference | operator() (const TableIndices< N > &indices) |
AlignedVector< number >::const_reference | operator() (const TableIndices< N > &indices) const |
size_type | n_rows () const |
size_type | n_cols () const |
iterator | begin () |
const_iterator | begin () const |
iterator | end () |
const_iterator | end () const |
bool | operator== (const TableBase< N, number > &T2) const |
void | reset_values () |
void | clear () |
size_type | size (const unsigned int i) const |
const TableIndices< N > & | size () const |
size_type | n_elements () const |
bool | empty () const |
void | fill (InputIterator entries, const bool C_style_indexing=true) |
void | fill (const number &value) |
void | replicate_across_communicator (const MPI_Comm communicator, const unsigned int root_process) |
void | swap (TableBase< N, number > &v) noexcept |
std::size_t | memory_consumption () const |
void | serialize (Archive &ar, const unsigned int version) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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 |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions | |
reference | el (const size_type i, const size_type j) |
const_reference | el (const size_type i, const size_type j) const |
AlignedVector< number >::reference | el (const TableIndices< N > &indices) |
AlignedVector< number >::const_reference | el (const TableIndices< N > &indices) const |
size_type | position (const TableIndices< N > &indices) const |
Protected Attributes | |
AlignedVector< number > | values |
TableIndices< N > | table_size |
Private Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
number | norm (const char type) const |
void | check_no_subscribers () const noexcept |
Private Attributes | |
LAPACKSupport::State | state |
LAPACKSupport::Property | property |
std::vector< number > | work |
std::vector< types::blas_int > | iwork |
std::vector< types::blas_int > | ipiv |
std::vector< number > | inv_work |
std::vector< typename numbers::NumberTraits< number >::real_type > | wr |
std::vector< number > | wi |
std::vector< number > | vl |
std::vector< number > | vr |
std::unique_ptr< LAPACKFullMatrix< number > > | svd_u |
std::unique_ptr< LAPACKFullMatrix< number > > | svd_vt |
Threads::Mutex | mutex |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
A variant of FullMatrix using LAPACK functions wherever possible. In order to do this, the matrix is stored in transposed order. The element access functions hide this fact by reverting the transposition.
Definition at line 58 of file lapack_full_matrix.h.
using LAPACKFullMatrix< number >::size_type = std::make_unsigned_t<types::blas_int> |
Declare type for container size.
Definition at line 64 of file lapack_full_matrix.h.
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
privateinherited |
The data type used in counter_map.
Definition at line 229 of file subscriptor.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 234 of file subscriptor.h.
|
explicit |
Constructor. Initialize the matrix as a square matrix with dimension size
.
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.
Definition at line 246 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number >::LAPACKFullMatrix | ( | const size_type | rows, |
const size_type | cols ) |
Constructor. Initialize the matrix as a rectangular matrix \(\rm{rows} \times \rm{cols}\).
Definition at line 255 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number >::LAPACKFullMatrix | ( | const LAPACKFullMatrix< number > & | M | ) |
Copy constructor. This constructor does a deep copy of the matrix. Therefore, it poses a possible efficiency problem, if for example, function arguments are passed by value rather than by reference. Unfortunately, we can't mark this copy constructor explicit
, since that prevents the use of this class in containers, such as std::vector
. The responsibility to check performance of programs must therefore remain with the user of this class.
Definition at line 264 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const LAPACKFullMatrix< number > & | M | ) |
Assignment operator.
Definition at line 274 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const FullMatrix< number2 > & | M | ) |
Assignment operator from a regular FullMatrix.
Definition at line 382 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const SparseMatrix< number2 > & | M | ) |
Assignment operator from a regular SparseMatrix.
Definition at line 400 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const number | d | ) |
This operator assigns a scalar to a matrix. To avoid confusion with constructors, zero (when cast to the number
type) is the only value allowed for d
.
Definition at line 417 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator*= | ( | const number | factor | ) |
This operator multiplies all entries by a fixed factor
.
Definition at line 433 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator/= | ( | const number | factor | ) |
This operator divides all entries by a fixed factor
.
Definition at line 462 of file lapack_full_matrix.cc.
|
inline |
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 column index of the element to be set. |
value | The value to be written into the element. |
Definition at line 1035 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::add | ( | const number | a, |
const LAPACKFullMatrix< number > & | B ) |
Simple addition of a scaled matrix, i.e. \(\mathbf A \mathrel{+}= a \, \mathbf B\).
Definition at line 493 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::rank1_update | ( | const number | a, |
const Vector< number > & | v ) |
Perform a rank-1 update of a symmetric matrix \( \mathbf A \leftarrow \mathbf A + a \, \mathbf v \mathbf v^T \).
This function also works for Cholesky factorization. In that case, updating ( \(a>0\)) is performed via Givens rotations, whereas downdating ( \(a<0\)) via hyperbolic rotations. Note that the latter case might lead to a negative definite matrix in which case the error will be thrown (because Cholesky factorizations are only valid for symmetric and positive definite matrices).
Definition at line 620 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::apply_givens_rotation | ( | const std::array< number, 3 > & | csr, |
const size_type | i, | ||
const size_type | k, | ||
const bool | left = true ) |
Apply Givens rotation csr
(a triplet of cosine, sine and radius, see Utilities::LinearAlgebra::givens_rotation() for the definition of the rotation matrix \(\mathbf G\)) to this matrix in the plane spanned by the i'th
and k'th
unit vectors. If left
is true
, the rotation is applied from left \(\mathbf A \leftarrow \mathbf G \mathbf A\) and only rows i
and k
are affected. Otherwise, transpose of the rotation matrix is applied from right \(\mathbf A \leftarrow \mathbf A \mathbf G^T\) and only columns i
and k
are affected.
Definition at line 310 of file lapack_full_matrix.cc.
|
inline |
Assignment from different matrix classes, performing the usual conversion to the transposed format expected by LAPACK. This assignment operator uses iterators of the typename MatrixType. Therefore, sparse matrices are possible sources.
Definition at line 1065 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::reinit | ( | const size_type | size | ) |
Regenerate the current matrix by one that has the same properties as if it were created by the constructor of this class with the same argument list as this present function.
Definition at line 286 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::grow_or_shrink | ( | const size_type | size | ) |
Same as above but will preserve the values of matrix upon resizing. The original values of the matrix are kept on increasing the size
\[ \mathbf A \rightarrow \left( \begin{array}{cc} \mathbf A & \mathbf 0 \\ \mathbf 0 & \mathbf 0 \end{array} \right) \]
Whereas if the new size is smaller, the matrix will contain the upper left block of the original one
\[ \left( \begin{array}{cc} \mathbf A_{11} & \mathbf A_{12} \\ \mathbf A_{21} & \mathbf A_{22} \end{array} \right) \rightarrow \mathbf A_{11} \]
Definition at line 296 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::remove_row_and_column | ( | const size_type | row, |
const size_type | col ) |
Remove row row
and column col
from the matrix.
\[ \left( \begin{array}{ccc} \mathbf A_{11} & \mathbf a_{12} & \mathbf A_{13} \\ \mathbf a_{21}^T & a_{22} & \mathbf a_{23}^T \\ \mathbf A_{31} & \mathbf a_{32} & \mathbf A_{33} \end{array} \right) \rightarrow \left( \begin{array}{cc} \mathbf A_{11} & \mathbf A_{13} \\ \mathbf A_{31} & \mathbf A_{33} \end{array} \right) \]
Definition at line 344 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::reinit | ( | const size_type | rows, |
const size_type | cols ) |
Regenerate the current matrix by one that has the same properties as if it were created by the constructor of this class with the same argument list as this present function.
Definition at line 371 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::set_property | ( | const LAPACKSupport::Property | property | ) |
Assign property
to this matrix.
Definition at line 1409 of file lapack_full_matrix.cc.
|
inline |
Return the dimension of the codomain (or range) space.
Definition at line 1046 of file lapack_full_matrix.h.
|
inline |
Return the dimension of the domain space.
Definition at line 1055 of file lapack_full_matrix.h.
|
inline |
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
.
The final two arguments allow to enter a multiple of the source or its transpose.
Definition at line 1089 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::vmult | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v, | ||
const bool | adding = false ) const |
Matrix-vector-multiplication.
Depending on previous transformations recorded in state, the result of this function is one of
The optional parameter adding
determines, whether the result is stored in the vector \(\mathbf w = \mathbf A \cdot \mathbf v\) or added to it \(\mathbf w \mathrel{+}= \mathbf A \cdot \mathbf v\).
number2
only exists for compile-time compatibility with FullMatrix. Only the case number2
= number
is implemented due to limitations in the underlying LAPACK interface. All other variants throw an error upon invocation. Definition at line 1124 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::vmult | ( | Vector< number > & | w, |
const Vector< number > & | v, | ||
const bool | adding = false ) const |
Specialization of above function for compatible Vector::value_type.
Definition at line 660 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::vmult_add | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v ) const |
Adding Matrix-vector-multiplication \(\mathbf w \mathrel{+}= \mathbf A \cdot \mathbf v\).
See the documentation of vmult() for details on the implementation.
Definition at line 1138 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::vmult_add | ( | Vector< number > & | w, |
const Vector< number > & | v ) const |
Specialization of above function for compatible Vector::value_type.
Definition at line 935 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tvmult | ( | 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 the vector \(\mathbf w = \mathbf A^T \cdot \mathbf v\) or added to it \(\mathbf w \mathrel{+}= \mathbf A^T \cdot \mathbf v\).
See the documentation of vmult() for details on the implementation.
Definition at line 1151 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::Tvmult | ( | Vector< number > & | w, |
const Vector< number > & | v, | ||
const bool | adding = false ) const |
Specialization of above function for compatible Vector::value_type.
Definition at line 796 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tvmult_add | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v ) const |
Adding transpose matrix-vector-multiplication \(\mathbf w \mathrel{+}= \mathbf A^T \cdot \mathbf v\).
See the documentation of vmult() for details on the implementation.
Definition at line 1165 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::Tvmult_add | ( | Vector< number > & | w, |
const Vector< number > & | v ) const |
Specialization of above function for compatible Vector::value_type.
Definition at line 945 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false ) const |
Matrix-matrix-multiplication.
The optional parameter adding
determines, whether the result is stored in the matrix \(\mathbf C = \mathbf A \cdot \mathbf B\) or added to it \(\mathbf C \mathrel{+}= \mathbf A \cdot \mathbf B\).
A
and B
have compatible sizes and that C
already has the right size.This
function uses the BLAS function Xgemm.
Definition at line 955 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false ) const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 990 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tmmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false ) const |
Matrix-matrix-multiplication using transpose of this
.
The optional parameter adding
determines, whether the result is stored in the matrix \(\mathbf C = \mathbf A^T \cdot \mathbf B\) or added to it \(\mathbf C \mathrel{+}= \mathbf A^T \cdot \mathbf B\).
A
and B
have compatible sizes and that C
already has the right size.Definition at line 1126 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tmmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false ) const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 1184 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tmmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const Vector< number > & | V, | ||
const bool | adding = false ) const |
Matrix-matrix-multiplication using transpose of this
and a diagonal vector V
.
If the adding=false
then the result is stored in the matrix \(\mathbf C = \mathbf A^T \cdot \rm{diag}(\mathbf V) \cdot \mathbf B\) otherwise it is added \(\mathbf C \mathrel{+}= \mathbf A^T \cdot
\rm{diag}(\mathbf V) \cdot \mathbf B\).
A
, B
and V
have compatible sizes and that C
already has the right size.Definition at line 1026 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mTmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false ) const |
Matrix-matrix-multiplication using transpose of B
.
The optional parameter adding
determines, whether the result is stored in the matrix \(\mathbf C = \mathbf A \cdot \mathbf B^T\) or added to it \(\mathbf C \mathrel{+}= \mathbf A \cdot \mathbf B^T\).
A
and B
have compatible sizes and that C
already has the right size.Definition at line 1220 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mTmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false ) const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 1278 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::TmTmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | 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 the matrix \(\mathbf C = \mathbf A^T \cdot \mathbf B^T\) or added to it \(\mathbf C \mathrel{+}= \mathbf A^T \cdot \mathbf B^T\).
A
and B
have compatible sizes and that C
already has the right size.Definition at line 1314 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::TmTmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false ) const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 1349 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::transpose | ( | LAPACKFullMatrix< number > & | B | ) | const |
Performs out-place transposition. Matrix B
should be appropriately sized.
mkl_?omatcopy
will be used, otherwise transposition is done element by element. Definition at line 1088 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::scale_rows | ( | const Vector< number > & | V | ) |
Scale rows of this matrix by V
. This is equivalent to premultiplication with a diagonal matrix \(\mathbf A\leftarrow {\rm diag}(\mathbf V)\mathbf
A\).
Definition at line 1109 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_lu_factorization | ( | ) |
Compute the LU factorization of the matrix using LAPACK function Xgetrf.
Definition at line 1385 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_cholesky_factorization | ( | ) |
Compute the Cholesky factorization of the matrix using LAPACK function Xpotrf.
Definition at line 1499 of file lapack_full_matrix.cc.
number LAPACKFullMatrix< number >::reciprocal_condition_number | ( | const number | l1_norm | ) | const |
Estimate the reciprocal of the condition number \(1/k(\mathbf A)\) in \(L_1\) norm ( \(1/(||\mathbf A||_1 \, ||\mathbf A^{-1}||_1)\)) of a symmetric positive definite matrix using Cholesky factorization. This function can only be called if the matrix is already factorized.
error = std::numeric_limits<Number>::epsilon * k
. Alternatively one can get the number of accurate digits std::floor(std::log10(k))
.[in] | l1_norm | Is the \(L_1\) norm of the matrix before calling Cholesky factorization. It can be obtained by calling l1_norm(). |
Definition at line 1526 of file lapack_full_matrix.cc.
number LAPACKFullMatrix< number >::reciprocal_condition_number | ( | ) | const |
Estimate the reciprocal of the condition number \(1/k(\mathbf A)\) in \(L_1\) norm for triangular matrices. The matrix has to have the LAPACKSupport::Property set to either LAPACKSupport::Property::upper_triangular or LAPACKSupport::Property::lower_triangular, see set_property().
Definition at line 1559 of file lapack_full_matrix.cc.
number LAPACKFullMatrix< number >::determinant | ( | ) | const |
Compute the determinant of a matrix. As it requires the LU factorization of the matrix, this function can only be called after compute_lu_factorization() has been called.
Definition at line 1872 of file lapack_full_matrix.cc.
number LAPACKFullMatrix< number >::l1_norm | ( | ) | const |
Compute \(L_1\) norm.
Definition at line 1418 of file lapack_full_matrix.cc.
number LAPACKFullMatrix< number >::linfty_norm | ( | ) | const |
Compute \(L_\infty\) norm.
Definition at line 1428 of file lapack_full_matrix.cc.
number LAPACKFullMatrix< number >::frobenius_norm | ( | ) | const |
Compute Frobenius norm
Definition at line 1438 of file lapack_full_matrix.cc.
number LAPACKFullMatrix< number >::trace | ( | ) | const |
Compute trace of the matrix, i.e. the sum of the diagonal values. Obviously, the matrix needs to be quadratic for this function.
Definition at line 1481 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::invert | ( | ) |
Invert the matrix by first computing an LU/Cholesky factorization with the LAPACK function Xgetrf/Xpotrf and then building the actual inverse using Xgetri/Xpotri.
Definition at line 1718 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::solve | ( | Vector< number > & | v, |
const bool | transposed = false ) const |
Solve the linear system with right hand side v
and put the solution back to v
. The matrix should be either triangular or LU/Cholesky factorization should be previously computed.
The flag transposed indicates whether the solution of the transposed system is to be performed.
Definition at line 1761 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::solve | ( | LAPACKFullMatrix< number > & | B, |
const bool | transposed = false ) const |
Same as above but for multiple right hand sides (as many as there are columns in the matrix B
).
Definition at line 1804 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_eigenvalues | ( | const bool | right_eigenvectors = false, |
const bool | left_eigenvectors = false ) |
Compute eigenvalues of the matrix. After this routine has been called, eigenvalues can be retrieved using the eigenvalue() function. The matrix itself will be LAPACKSupport::unusable after this operation.
The optional arguments allow to compute left and right eigenvectors as well.
Note that the function does not return the computed eigenvalues right away since that involves copying data around between the output arrays of the LAPACK functions and any return array. This is often unnecessary since one may not be interested in all eigenvalues at once, but for example only the extreme ones. In that case, it is cheaper to just have this function compute the eigenvalues and have a separate function that returns whatever eigenvalue is requested.
Definition at line 1904 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_eigenvalues_symmetric | ( | const number | lower_bound, |
const number | upper_bound, | ||
const number | abs_accuracy, | ||
Vector< number > & | eigenvalues, | ||
FullMatrix< number > & | eigenvectors ) |
Compute eigenvalues and eigenvectors of a real symmetric matrix. Only eigenvalues in the interval \((\rm{lower\_bound}, \rm{upper\_bound}]\) are computed with the absolute tolerance \(\rm abs\_accuracy\). An approximate eigenvalue is accepted as converged when it is determined to lie in an interval \([a,b]\) of width less than or equal to \(\rm{abs\_accuracy} + eps * \rm{max}(|a|,|b|)\), where \(eps\) is the machine precision. If \(\rm{abs\_accuracy}\) is less than or equal to zero, then \(eps\,|\mathbf{T}|_1\) will be used in its place, where \(|\mathbf{T}|_1\) is the 1-norm of the tridiagonal matrix obtained by reducing \(\mathbf A\) to tridiagonal form. Eigenvalues will be computed most accurately when \(\rm{abs\_accuracy}\) is set to twice the underflow threshold, not zero. After this routine has been called, all eigenvalues in \((\rm{lower\_bound}, \rm{upper\_bound}]\) will be stored in eigenvalues and the corresponding eigenvectors will be stored in the columns of eigenvectors, whose dimension is set accordingly.
Definition at line 2103 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_generalized_eigenvalues_symmetric | ( | LAPACKFullMatrix< number > & | B, |
const number | lower_bound, | ||
const number | upper_bound, | ||
const number | abs_accuracy, | ||
Vector< number > & | eigenvalues, | ||
std::vector< Vector< number > > & | eigenvectors, | ||
const types::blas_int | itype = 1 ) |
Compute generalized eigenvalues and eigenvectors of a real generalized symmetric eigenproblem of the form
where \(\mathbf A\) is this matrix. \(\mathbf A\) and \(\mathbf B\) are assumed to be symmetric, and \(\mathbf B\) has to be positive definite. Only eigenvalues in the interval \((\rm{lower\_bound}, \rm{upper\_bound}]\) are computed with the absolute tolerance \(\rm{abs\_accuracy}\). An approximate eigenvalue is accepted as converged when it is determined to lie in an interval \([a,b]\) of width less than or equal to \(\rm{abs\_accuracy} + eps * \rm{max}( |a|,|b| )\), where \(eps\) is the machine precision. If \(\rm{abs\_accuracy}\) is less than or equal to zero, then \(eps \, |\mathbf{T}|_1\) will be used in its place, where \(|\mathbf{T}|_1\) is the 1-norm of the tridiagonal matrix obtained by reducing \(\mathbf A\) to tridiagonal form. Eigenvalues will be computed most accurately when \(\rm{abs\_accuracy}\) is set to twice the underflow threshold, not zero. After this routine has been called, all eigenvalues in \((\rm{lower\_bound}, \rm{upper\_bound}]\) will be stored in eigenvalues and the corresponding eigenvectors will be stored in eigenvectors, whose dimension is set accordingly.
Definition at line 2235 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_generalized_eigenvalues_symmetric | ( | LAPACKFullMatrix< number > & | B, |
std::vector< Vector< number > > & | eigenvectors, | ||
const types::blas_int | itype = 1 ) |
Same as the other compute_generalized_eigenvalues_symmetric function except that all eigenvalues are computed and the tolerance is set automatically. Note that this function does not return the computed eigenvalues right away since that involves copying data around between the output arrays of the LAPACK functions and any return array. This is often unnecessary since one may not be interested in all eigenvalues at once, but for example only the extreme ones. In that case, it is cheaper to just have this function compute the eigenvalues and have a separate function that returns whatever eigenvalue is requested. Eigenvalues can be retrieved using the eigenvalue() function. The number of computed eigenvectors is equal to eigenvectors.size()
Definition at line 2393 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_svd | ( | ) |
Compute the singular value decomposition of the matrix using LAPACK function Xgesdd.
Requires that the state is LAPACKSupport::matrix, fills the data members wr, svd_u, and svd_vt, and leaves the object in the state LAPACKSupport::svd.
The singular value decomposition factorizes the provided matrix (A) into three parts: U, sigma, and the transpose of V (V^T), such that A = U sigma V^T. Sigma is a MxN matrix which contains the singular values of A on the diagonal while all the other elements are zero. U is a MxM orthogonal matrix containing the left singular vectors corresponding to the singular values of A. V is a NxN orthonal matrix containing the right singular vectors corresponding the singular values of A.
Note that the variable svd_vt contains the transpose of V and can be accessed by get_svd_vt(), while U is accessed with get_svd_u().
Definition at line 1597 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_inverse_svd | ( | const double | threshold = 0. | ) |
Compute the inverse of the matrix by singular value decomposition.
Requires that state is either LAPACKSupport::matrix or LAPACKSupport::svd. In the first case, this function calls compute_svd(). After this function, the object will have the state LAPACKSupport::inverse_svd.
For a singular value decomposition, the inverse is simply computed by replacing all singular values by their reciprocal values. If the matrix does not have maximal rank, singular values 0 are not touched, thus computing the minimal norm right inverse of the matrix.
The parameter threshold
determines, when a singular value should be considered zero. It is the ratio of the smallest to the largest nonzero singular value \(s_{max}\). Thus, the inverses of all singular values less than \(s_{max}/\rm{threshold}\) will be set to zero.
Definition at line 1674 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_inverse_svd_with_kernel | ( | const unsigned int | kernel_size | ) |
Same as above but provide the size of the kernel instead of a threshold, i.e. the kernel_size
smallest eigenvalues.
Definition at line 1697 of file lapack_full_matrix.cc.
|
inline |
Retrieve eigenvalue after compute_eigenvalues() was called. The return type is always a complex number: In case the underlying matrix is real-valued, the type std::complex<number>
is returned, whereas number
is returned for complex numbers.
Definition at line 1201 of file lapack_full_matrix.h.
FullMatrix< std::complex< typename numbers::NumberTraits< number >::real_type > > LAPACKFullMatrix< number >::get_right_eigenvectors | ( | ) | const |
After a call to compute_eigenvalues(), this function returns the \(n\times
n\) matrix of (right) eigenvectors in a decomposition of the form \(A V = V
\Lambda\). Note that this function constructs the associated matrix on the fly, since LAPACK packs complex-conjugate eigenvalue/eigenvector pairs of real-valued matrices into a real-valued return matrix. This call only succeeds in case the respective flag right_eigenvectors
in compute_eigenvalues() has been set to true.
Definition at line 2063 of file lapack_full_matrix.cc.
FullMatrix< std::complex< typename numbers::NumberTraits< number >::real_type > > LAPACKFullMatrix< number >::get_left_eigenvectors | ( | ) | const |
Return the matrix of left eigenvectors after a call to compute_eigenvalues(), following the same principles as get_right_eigenvectors().
Definition at line 2083 of file lapack_full_matrix.cc.
|
inline |
Retrieve singular values after compute_svd() or compute_inverse_svd() was called.
Definition at line 1215 of file lapack_full_matrix.h.
|
inline |
Retrieve the matrix svd_u after compute_svd() or compute_inverse_svd() was called.
Definition at line 1228 of file lapack_full_matrix.h.
|
inline |
Retrieve the matrix svd_vt after compute_svd() or compute_inverse_svd() was called.
Definition at line 1240 of file lapack_full_matrix.h.
void LAPACKFullMatrix< 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 char * | separator = " " ) const |
Print the matrix and allow formatting of entries.
The parameters allow for a flexible setting of the output format:
out | This specifies the stream to write to. |
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. |
separator | specifies a string printed to separate row entries. |
Definition at line 2519 of file lapack_full_matrix.cc.
|
private |
Internal function to compute various norms.
Definition at line 1448 of file lapack_full_matrix.cc.
|
inherited |
Reinitialize the object. This function is mostly here for compatibility with the earlier vector2d
class. Passes down to the base class by converting the arguments to the data type requested by the base class.
|
inherited |
Set the dimensions of this object to the sizes given in the first argument, and allocate the required memory for table entries to accommodate these sizes. If omit_default_initialization
is set to false
, all elements of the table are set to a default constructed object for the element type. Otherwise the memory is left in an uninitialized or otherwise undefined state.
|
inherited |
Direct access to one element of the table by specifying all indices at the same time. Range checks are performed.
This version of the function only allows read access.
|
inherited |
Direct access to one element of the table by specifying all indices at the same time. Range checks are performed.
This version of the function allows read-write access.
|
inherited |
Return a read-write reference to the indicated element.
|
inherited |
Return the value of the indicated element as a read-only reference.
We return the requested value as a constant reference rather than by value since this object may hold data types that may be large, and we don't know here whether copying is expensive or not.
|
inherited |
Number of rows. This function really makes only sense since we have a two-dimensional object here.
|
inherited |
Number of columns. This function really makes only sense since we have a two-dimensional object here.
|
inherited |
Return an iterator pointing to the first entry.
|
inherited |
Return a constant iterator pointing to the first entry.
|
inherited |
Return an iterator pointing to one past the last entry.
|
inherited |
Return a constant iterator pointing to one past the last entry.
|
protectedinherited |
Return a read-write reference to the element (i,j)
.
This function does no bounds checking and is only to be used internally and in functions already checked.
These functions are mainly here for compatibility with a former implementation of these table classes for 2d arrays, then called vector2d
.
|
protectedinherited |
Return the value of the element (i,j)
as a read-only reference.
This function does no bounds checking and is only to be used internally and in functions already checked.
We return the requested value as a constant reference rather than by value since this object may hold data types that may be large, and we don't know here whether copying is expensive or not.
These functions are mainly here for compatibility with a former implementation of these table classes for 2d arrays, then called vector2d
.
|
protectedinherited |
Return a read-write reference to the indicated element.
This function does no bounds checking and is only to be used internally and in functions already checked.
|
protectedinherited |
Return the value of the indicated element as a read-only reference.
This function does no bounds checking and is only to be used internally and in functions already checked.
We return the requested value as a constant reference rather than by value since this object may hold data types that may be large, and we don't know here whether copying is expensive or not.
Test for equality of two tables.
|
inherited |
Set all entries to their default value (i.e. copy them over with default constructed objects). Do not change the size of the table, though.
|
inherited |
Set all dimensions to zero.
Size of the table in direction i
.
|
inherited |
Return the sizes of this object in each direction.
Return the number of elements stored in this object, which is the product of the extensions in each dimension.
Return whether the object is empty, i.e. one of the directions is zero. This is equivalent to n_elements()==0
.
|
inherited |
Fill this table (which is assumed to already have the correct size) from a source given by dereferencing the given forward iterator (which could, for example, be a pointer to the first element of an array, or an inserting std::istream_iterator). The second argument denotes whether the elements pointed to are arranged in a way that corresponds to the last index running fastest or slowest. The default is to use C-style indexing where the last index runs fastest (as opposed to Fortran-style where the first index runs fastest when traversing multidimensional arrays. For example, if you try to fill an object of type Table<2,T>, then calling this function with the default value for the second argument will result in the equivalent of doing
On the other hand, if the second argument to this function is false, then this would result in code of the following form:
Note the switched order in which we fill the table elements by traversing the given set of iterators.
entries | An iterator to a set of elements from which to initialize this table. It is assumed that iterator can be incremented and dereferenced a sufficient number of times to fill this table. |
C_style_indexing | If true, run over elements of the table with the last index changing fastest as we dereference subsequent elements of the input range. If false, change the first index fastest. |
|
inherited |
Fill all table entries with the same value.
|
inherited |
This function replicates the state found on the process indicated by root_process
across all processes of the MPI communicator. The current state found on any of the processes other than root_process
is lost in this process. One can imagine this operation to act like a call to Utilities::MPI::broadcast() from the root process to all other processes, though in practice the function may try to move the data into shared memory regions on each of the machines that host MPI processes and let all MPI processes on this machine then access this shared memory region instead of keeping their own copy. See the general documentation of this class for a code example.
The intent of this function is to quickly exchange large arrays from one process to others, rather than having to compute or create it on all processes. This is specifically the case for data loaded from disk – say, large data tables – that are more easily dealt with by reading once and then distributing across all processes in an MPI universe, than letting each process read the data from disk itself. Specifically, the use of shared memory regions allows for replicating the data only once per multicore machine in the MPI universe, rather than replicating data once for each MPI process. This results in large memory savings if the data is large on today's machines that can easily house several dozen MPI processes per shared memory space.
This function does not imply a model of keeping data on different processes in sync, as LinearAlgebra::distributed::Vector and other vector classes do where there exists a notion of certain elements of the vector owned by each process and possibly ghost elements that are mirrored from its owning process to other processes. Rather, the elements of the current object are simply copied to the other processes, and it is useful to think of this operation as creating a set of const
TableBase objects on all processes that should not be changed any more after the replication operation, as this is the only way to ensure that the vectors remain the same on all processes. This is particularly true because of the use of shared memory regions where any modification of a vector element on one MPI process may also result in a modification of elements visible on other processes, assuming they are located within one shared memory node.
communicator
object, which is generally an expensive operation. Likewise, the generation of shared memory spaces is not a cheap operation. As a consequence, this function primarily makes sense when the goal is to share large read-only data tables among processes; examples are data tables that are loaded at start-up time and then used over the course of the run time of the program. In such cases, the start-up cost of running this function can be amortized over time, and the potential memory savings from not having to store the table on each process may be substantial on machines with large core counts on which many MPI processes run on the same machine.T
is "self-contained", i.e., all of its information is stored in its member variables, and if none of the member variables are pointers to other parts of the memory. This is because if a type T
does have pointers to other parts of memory, then moving T
into a shared memory space does not result in the other processes having access to data that the object points to with its member variable pointers: These continue to live only on one process, and are typically in memory areas not accessible to the other processes. As a consequence, the usual use case for this function is to share arrays of simple objects such as double
s or int
s.root_process
throw away their data, which is not a collective operation. Generally, these restrictions on what can and can not be done hint at the correctness of the comments above: You should treat an AlignedVector on which the current function has been called as const
, on which no further operations can be performed until the destructor is called. Swap the contents of this table and the other table v
. One could do this operation with a temporary variable and copying over the data elements, but this function is significantly more efficient since it only swaps the pointers to the data of the two vectors and therefore does not need to allocate temporary storage and move data around.
This function is analogous to the swap
function of all C++ standard containers. Also, there is a global function swap(u,v)
that simply calls u.swap(v)
, again in analogy to standard functions.
|
inherited |
Determine an estimate for the memory consumption (in bytes) of this object.
Write or read the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
|
protectedinherited |
Return the position of the indicated element within the array of elements stored one after the other. This function does no index checking.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 135 of file subscriptor.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 155 of file subscriptor.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 300 of file subscriptor.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 317 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 203 of file subscriptor.cc.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 52 of file subscriptor.cc.
|
private |
Since LAPACK operations notoriously change the meaning of the matrix entries, we record the current state after the last operation here.
Definition at line 931 of file lapack_full_matrix.h.
|
private |
Additional property of the matrix which may help to select more efficient LAPACK functions.
Definition at line 937 of file lapack_full_matrix.h.
|
mutableprivate |
The working array used for some LAPACK functions.
Definition at line 942 of file lapack_full_matrix.h.
|
mutableprivate |
Integer working array used for some LAPACK functions.
Definition at line 947 of file lapack_full_matrix.h.
|
private |
The vector storing the permutations applied for pivoting in the LU- factorization.
Also used as the scratch array IWORK for LAPACK functions needing it.
Definition at line 955 of file lapack_full_matrix.h.
|
private |
Workspace for calculating the inverse matrix from an LU factorization.
Definition at line 960 of file lapack_full_matrix.h.
|
private |
Real parts of eigenvalues or the singular values. Filled by compute_eigenvalues() or compute_svd().
Definition at line 966 of file lapack_full_matrix.h.
|
private |
Imaginary parts of eigenvalues, or, in the complex scalar case, the eigenvalues themselves. Filled by compute_eigenvalues.
Definition at line 972 of file lapack_full_matrix.h.
|
private |
Space where left eigenvectors can be stored.
Definition at line 977 of file lapack_full_matrix.h.
|
private |
Space where right eigenvectors can be stored.
Definition at line 982 of file lapack_full_matrix.h.
|
private |
The matrix \(\mathbf U\) in the singular value decomposition \(\mathbf U \cdot \mathbf S \cdot \mathbf V^T\).
Definition at line 988 of file lapack_full_matrix.h.
|
private |
The matrix \(\mathbf V^T\) in the singular value decomposition \(\mathbf U \cdot \mathbf S \cdot \mathbf V^T\).
Definition at line 994 of file lapack_full_matrix.h.
|
mutableprivate |
Thread mutex.
Definition at line 999 of file lapack_full_matrix.h.
|
protectedinherited |
|
protectedinherited |
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 218 of file subscriptor.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 224 of file subscriptor.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.
Definition at line 240 of file subscriptor.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 248 of file subscriptor.h.