Important changes to repositories hosted on mbed.com
Mbed hosted mercurial repositories are deprecated and are due to be permanently deleted in July 2026.
To keep a copy of this software download the repository Zip archive or clone locally using Mercurial.
It is also possible to export all your personal repositories from the account settings page.
Dependents: Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more
Eigen::internal Namespace Reference
Applies the clock wise 2D rotation j to the set of 2D vectors of cordinates x and y:
More...
Namespaces | |
namespace | anonymous_namespace{Macros.h} |
Data Structures | |
class | BandMatrix |
Represents a rectangular matrix with a banded storage. More... | |
class | TridiagonalMatrix |
Represents a tridiagonal matrix with a compact banded storage. More... | |
struct | dense_xpr_base_dispatcher_for_doxygen< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > > |
This class is just a workaround for Doxygen and it does not not actually exist. More... | |
struct | dense_xpr_base_dispatcher_for_doxygen< Array< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > > |
This class is just a workaround for Doxygen and it does not not actually exist. More... | |
struct | HessenbergDecompositionMatrixHReturnType |
More... | |
struct | FullPivHouseholderQRMatrixQReturnType |
Expression type for return value of FullPivHouseholderQR::matrixQ() More... | |
Functions | |
template<typename LhsScalar , typename RhsScalar , int KcFactor, typename SizeType > | |
void | computeProductBlockingSizes (SizeType &k, SizeType &m, SizeType &n) |
Computes the blocking parameters for a m x k times k x n matrix product. | |
template<typename MatrixType , typename DiagonalType , typename SubDiagonalType > | |
void | tridiagonalization_inplace (MatrixType &mat, DiagonalType &diag, SubDiagonalType &subdiag, bool extractQ) |
Performs a full tridiagonalization in place. |
Detailed Description
Applies the clock wise 2D rotation j to the set of 2D vectors of cordinates x and y:
Function Documentation
void Eigen::internal::computeProductBlockingSizes | ( | SizeType & | k, |
SizeType & | m, | ||
SizeType & | n | ||
) |
Computes the blocking parameters for a m x k times k x n matrix product.
- Parameters:
-
[in,out] k Input: the third dimension of the product. Output: the blocking size along the same dimension. [in,out] m Input: the number of rows of the left hand side. Output: the blocking size along the same dimension. [in,out] n Input: the number of columns of the right hand side. Output: the blocking size along the same dimension.
Given a m x k times k x n matrix product of scalar types LhsScalar
and RhsScalar
, this function computes the blocking size parameters along the respective dimensions for matrix products and related algorithms. The blocking sizes depends on various parameters:
- the L1 and L2 cache sizes,
- the register level blocking sizes defined by gebp_traits,
- the number of scalars that fit into a packet (when vectorization is enabled).
- See also:
- setCpuCacheSizes
Definition at line 73 of file GeneralBlockPanelKernel.h.
void Eigen::internal::tridiagonalization_inplace | ( | MatrixType & | mat, |
DiagonalType & | diag, | ||
SubDiagonalType & | subdiag, | ||
bool | extractQ | ||
) |
Performs a full tridiagonalization in place.
- Parameters:
-
[in,out] mat On input, the selfadjoint matrix whose tridiagonal decomposition is to be computed. Only the lower triangular part referenced. The rest is left unchanged. On output, the orthogonal matrix Q in the decomposition if extractQ
is true.[out] diag The diagonal of the tridiagonal matrix T in the decomposition. [out] subdiag The subdiagonal of the tridiagonal matrix T in the decomposition. [in] extractQ If true, the orthogonal matrix Q in the decomposition is computed and stored in mat
.
Computes the tridiagonal decomposition of the selfadjoint matrix mat
in place such that where
is unitary and
a real symmetric tridiagonal matrix.
The tridiagonal matrix T is passed to the output parameters diag
and subdiag
. If extractQ
is true, then the orthogonal matrix Q is passed to mat
. Otherwise the lower part of the matrix mat
is destroyed.
The vectors diag
and subdiag
are not resized. The function assumes that they are already of the correct size. The length of the vector diag
should equal the number of rows in mat
, and the length of the vector subdiag
should be one left.
This implementation contains an optimized path for 3-by-3 matrices which is especially useful for plane fitting.
- Note:
- Currently, it requires two temporary vectors to hold the intermediate Householder coefficients, and to reconstruct the matrix Q from the Householder reflectors.
Example (this uses the same matrix as the example in Tridiagonalization::Tridiagonalization(const MatrixType&)):
Output:
- See also:
- class Tridiagonalization
Definition at line 427 of file Tridiagonalization.h.
Generated on Tue Jul 12 2022 17:47:06 by
