![]() |
Open CASCADE Technology Reference Manual 8.0.0
|
Namespaces | |
| namespace | Internal |
Data Structures | |
| struct | CroutResult |
| Result for Crout LDL^T decomposition. Specialized for symmetric matrices. More... | |
| struct | EigenResult |
| Result for eigenvalue decomposition of tridiagonal matrix. More... | |
| struct | LeastSquaresResult |
| Result for least squares problems. More... | |
| struct | LUResult |
| Result for LU decomposition. More... | |
| struct | QRResult |
| Result for QR decomposition using Householder reflections. More... | |
| struct | SVDResult |
| Result for SVD decomposition. More... | |
Enumerations | |
| enum class | LeastSquaresMethod { NormalEquations , QR , SVD } |
| Method for solving least squares problems. More... | |
Functions | |
| LUResult | LU (const math_Matrix &theA, double theMinPivot=1.0e-20) |
| Perform LU decomposition of matrix A with partial pivoting. Decomposes A into L*U where L is lower triangular with unit diagonal and U is upper triangular. The result stores L and U in a combined matrix. | |
| LinearResult | Solve (const math_Matrix &theA, const math_Vector &theB, double theMinPivot=1.0e-20) |
| Solve linear system AX = B using LU decomposition. | |
| LinearMultipleResult | SolveMultiple (const math_Matrix &theA, const math_Matrix &theB, double theMinPivot=1.0e-20) |
| Solve multiple linear systems AX = B where B is a matrix. Each column of B is a separate right-hand side. | |
| LinearResult | Determinant (const math_Matrix &theA, double theMinPivot=1.0e-20) |
| Compute determinant of matrix A. | |
| InverseResult | Invert (const math_Matrix &theA, double theMinPivot=1.0e-20) |
| Compute inverse of matrix A. | |
| CroutResult | Crout (const math_Matrix &theA, double theMinPivot=1.0e-20) |
| Crout decomposition for symmetric matrices: A = L * D * L^T. | |
| LinearResult | SolveCrout (const math_Matrix &theA, const math_Vector &theB, double theMinPivot=1.0e-20) |
| Solve symmetric linear system Ax = b using Crout decomposition. | |
| InverseResult | InvertCrout (const math_Matrix &theA, double theMinPivot=1.0e-20) |
| Compute inverse of symmetric matrix using Crout decomposition. | |
| SVDResult | SVD (const math_Matrix &theA, double theTolerance=1.0e-15) |
| Singular Value Decomposition: A = U * diag(S) * V^T. | |
| LinearResult | SolveSVD (const math_Matrix &theA, const math_Vector &theB, double theTolerance=1.0e-6) |
| Solve linear system Ax = b using SVD decomposition. This is particularly useful for ill-conditioned or singular systems. | |
| InverseResult | PseudoInverse (const math_Matrix &theA, double theTolerance=1.0e-6) |
| Compute pseudo-inverse (Moore-Penrose inverse) of matrix A. A^+ = V * diag(1/w) * U^T where singular values below threshold are set to 0. | |
| double | ConditionNumber (const math_Matrix &theA) |
| Compute condition number of matrix using SVD. Condition number = sigma_max / sigma_min (ratio of largest to smallest singular value). | |
| int | NumericalRank (const math_Matrix &theA, double theTolerance=1.0e-15) |
| Compute numerical rank of matrix using SVD. Rank is the number of singular values above the threshold. | |
| QRResult | QR (const math_Matrix &theA, double theTolerance=1.0e-20) |
| QR decomposition using Householder reflections: A = Q * R. | |
| LinearResult | SolveQR (const math_Matrix &theA, const math_Vector &theB, double theTolerance=1.0e-20) |
| Solve overdetermined system Ax = b using QR decomposition (least squares). | |
| LinearMultipleResult | SolveQRMultiple (const math_Matrix &theA, const math_Matrix &theB, double theTolerance=1.0e-20) |
| Solve multiple right-hand sides using QR decomposition. | |
| EigenResult | Jacobi (const math_Matrix &theA, bool theSortDescending=true) |
| Compute eigenvalues and eigenvectors of a symmetric matrix using the Jacobi iterative method. | |
| EigenResult | EigenValues (const math_Matrix &theA, bool theSortDescending=true) |
| Compute only eigenvalues of a symmetric matrix (faster). | |
| EigenResult | SpectralDecomposition (const math_Matrix &theA) |
| Compute spectral decomposition A = V * D * V^T. | |
| std::optional< math_Matrix > | MatrixPower (const math_Matrix &theA, double thePower) |
| Compute matrix power A^p for symmetric positive semi-definite matrix. | |
| std::optional< math_Matrix > | MatrixSqrt (const math_Matrix &theA) |
| Compute matrix square root of symmetric positive semi-definite matrix. | |
| std::optional< math_Matrix > | MatrixInvSqrt (const math_Matrix &theA) |
| Compute matrix inverse square root of symmetric positive definite matrix. | |
| LeastSquaresResult | LeastSquares (const math_Matrix &theA, const math_Vector &theB, LeastSquaresMethod theMethod=LeastSquaresMethod::QR, double theTolerance=1.0e-15) |
| Solve overdetermined linear least squares: minimize ||Ax - b||_2. | |
| LeastSquaresResult | WeightedLeastSquares (const math_Matrix &theA, const math_Vector &theB, const math_Vector &theW, LeastSquaresMethod theMethod=LeastSquaresMethod::QR, double theTolerance=1.0e-15) |
| Solve weighted least squares: minimize ||W^{1/2}(Ax - b)||_2. | |
| LeastSquaresResult | RegularizedLeastSquares (const math_Matrix &theA, const math_Vector &theB, double theLambda, double theTolerance=1.0e-15) |
| Solve regularized least squares (Tikhonov/Ridge regression): minimize ||Ax - b||_2^2 + lambda*||x||_2^2. | |
| double | OptimalRegularization (const math_Matrix &theA, const math_Vector &theB, double theLambdaMin=1.0e-10, double theLambdaMax=1.0e2, int theNbPoints=20) |
| Compute optimal regularization parameter using Leave-One-Out Cross-Validation. | |
| EigenResult | EigenTridiagonal (const math_Vector &theDiagonal, const math_Vector &theSubdiagonal, int theMaxIterations=30) |
| Eigenvalue decomposition of symmetric tridiagonal matrix using QL algorithm. | |
| math_Vector | GetEigenVector (const EigenResult &theResult, int theIndex) |
| Get a single eigenvector from the result. | |
|
strong |
|
inline |
Compute condition number of matrix using SVD. Condition number = sigma_max / sigma_min (ratio of largest to smallest singular value).
High condition number (> 1e10) indicates ill-conditioned matrix.
| theA | input matrix |
|
inline |
Crout decomposition for symmetric matrices: A = L * D * L^T.
This algorithm decomposes a symmetric matrix A into:
Properties:
| theA | input symmetric matrix (only lower triangle used) |
| theMinPivot | minimum pivot value (smaller treated as singular) |
|
inline |
Compute determinant of matrix A.
| theA | input square matrix |
| theMinPivot | minimum pivot value |
|
inline |
Eigenvalue decomposition of symmetric tridiagonal matrix using QL algorithm.
The QL algorithm with implicit Wilkinson shifts finds all eigenvalues and eigenvectors of a symmetric tridiagonal matrix T = Q * D * Q^T.
Properties:
| theDiagonal | diagonal elements of the tridiagonal matrix |
| theSubdiagonal | subdiagonal elements (one less than diagonal) |
| theMaxIterations | maximum iterations per eigenvalue (default 30) |
|
inline |
Compute only eigenvalues of a symmetric matrix (faster).
Uses the same Jacobi method but may be optimized to not store eigenvectors if not needed.
| theA | input symmetric matrix (n x n) |
| theSortDescending | if true, eigenvalues are sorted in descending order |
|
inline |
Get a single eigenvector from the result.
| theResult | eigenvalue decomposition result |
| theIndex | 1-based index of eigenvector |
|
inline |
Compute inverse of matrix A.
| theA | input square matrix |
| theMinPivot | minimum pivot value |
|
inline |
Compute inverse of symmetric matrix using Crout decomposition.
| theA | input symmetric matrix |
| theMinPivot | minimum pivot value |
|
inline |
Compute eigenvalues and eigenvectors of a symmetric matrix using the Jacobi iterative method.
The Jacobi method applies a sequence of plane rotations (Givens rotations) to diagonalize the symmetric matrix A: A' = R^T * A * R where R is a rotation that zeroes one off-diagonal element.
After convergence, A is diagonal with eigenvalues on the diagonal, and the accumulated rotations form the eigenvector matrix.
Properties:
Complexity: O(n^3) per sweep, typically needs 5-10 sweeps.
| theA | input symmetric matrix (n x n) |
| theSortDescending | if true, eigenvalues are sorted in descending order |
|
inline |
Solve overdetermined linear least squares: minimize ||Ax - b||_2.
Given m x n matrix A (m >= n) and m-vector b, finds n-vector x that minimizes the 2-norm of the residual r = Ax - b.
Methods:
| theA | coefficient matrix (m x n, m >= n) |
| theB | right-hand side vector (length m) |
| theMethod | solution method (default: QR) |
| theTolerance | for rank/singularity detection |
|
inline |
Perform LU decomposition of matrix A with partial pivoting. Decomposes A into L*U where L is lower triangular with unit diagonal and U is upper triangular. The result stores L and U in a combined matrix.
| theA | input square matrix |
| theMinPivot | minimum pivot value (smaller treated as singular) |
|
inline |
Compute matrix inverse square root of symmetric positive definite matrix.
| theA | input symmetric positive definite matrix |
|
inline |
Compute matrix power A^p for symmetric positive semi-definite matrix.
Uses spectral decomposition: A^p = V * D^p * V^T where D^p is the diagonal matrix with eigenvalues raised to power p.
| theA | input symmetric positive semi-definite matrix |
| thePower | exponent (can be fractional, e.g., 0.5 for sqrt) |
|
inline |
Compute matrix square root of symmetric positive semi-definite matrix.
| theA | input symmetric positive semi-definite matrix |
|
inline |
Compute numerical rank of matrix using SVD. Rank is the number of singular values above the threshold.
| theA | input matrix |
| theTolerance | relative tolerance for singular values |
|
inline |
Compute optimal regularization parameter using Leave-One-Out Cross-Validation.
Minimizes the LOO-CV score: sum_i (a_i^T * x_{-i} - b_i)^2 where x_{-i} is the solution with the i-th observation removed.
| theA | coefficient matrix |
| theB | right-hand side vector |
| theLambdaMin | minimum lambda to consider |
| theLambdaMax | maximum lambda to consider |
| theNbPoints | number of lambda values to try |
|
inline |
Compute pseudo-inverse (Moore-Penrose inverse) of matrix A. A^+ = V * diag(1/w) * U^T where singular values below threshold are set to 0.
Properties:
| theA | input matrix (m x n) |
| theTolerance | for singular value threshold |
|
inline |
QR decomposition using Householder reflections: A = Q * R.
Decomposes an m x n matrix A (m >= n) into:
The Householder method applies orthogonal transformations to reduce A to upper triangular form. It is more numerically stable than Gram-Schmidt orthogonalization.
Uses: Least squares problems, orthogonalization, computing determinant sign.
| theA | input matrix A (m x n, m >= n) |
| theTolerance | for rank determination |
|
inline |
Solve regularized least squares (Tikhonov/Ridge regression): minimize ||Ax - b||_2^2 + lambda*||x||_2^2.
Adds regularization to stabilize ill-conditioned problems. The solution is: x = (A^T*A + lambda*I)^{-1} * A^T * b
| theA | coefficient matrix (m x n) |
| theB | right-hand side vector (length m) |
| theLambda | regularization parameter (>= 0) |
| theTolerance | for singularity detection |
|
inline |
Solve linear system AX = B using LU decomposition.
| theA | coefficient matrix (square) |
| theB | right-hand side vector |
| theMinPivot | minimum pivot value |
|
inline |
Solve symmetric linear system Ax = b using Crout decomposition.
Uses precomputed Crout decomposition to solve for x. More efficient than LU for symmetric positive definite matrices.
| theA | coefficient matrix (symmetric) |
| theB | right-hand side vector |
| theMinPivot | minimum pivot value |
|
inline |
Solve multiple linear systems AX = B where B is a matrix. Each column of B is a separate right-hand side.
| theA | coefficient matrix (square) |
| theB | right-hand side matrix |
| theMinPivot | minimum pivot value |
|
inline |
Solve overdetermined system Ax = b using QR decomposition (least squares).
For m x n system with m > n, finds x that minimizes ||Ax - b||_2.
Algorithm:
| theA | coefficient matrix (m x n, m >= n) |
| theB | right-hand side vector (length m) |
| theTolerance | for singularity detection |
|
inline |
Solve multiple right-hand sides using QR decomposition.
| theA | coefficient matrix (m x n, m >= n) |
| theB | right-hand side matrix (m x p) |
| theTolerance | for singularity detection |
|
inline |
Solve linear system Ax = b using SVD decomposition. This is particularly useful for ill-conditioned or singular systems.
For overdetermined systems (m > n), finds the least squares solution. For underdetermined systems (m < n), finds the minimum norm solution.
| theA | coefficient matrix (m x n) |
| theB | right-hand side vector (length m) |
| theTolerance | for singular value threshold |
|
inline |
Compute spectral decomposition A = V * D * V^T.
For symmetric matrix A, decomposes into:
Such that A = V * D * V^T
| theA | input symmetric matrix (n x n) |
|
inline |
Singular Value Decomposition: A = U * diag(S) * V^T.
Decomposes an m x n matrix A into:
Properties:
| theA | input matrix A (m x n) |
| theTolerance | for rank determination (relative to largest singular value) |
|
inline |
Solve weighted least squares: minimize ||W^{1/2}(Ax - b)||_2.
Equivalent to minimizing sum of w_i * (a_i^T * x - b_i)^2 where w_i are the weights.
| theA | coefficient matrix (m x n) |
| theB | right-hand side vector (length m) |
| theW | weight vector (length m, positive values) |
| theMethod | solution method |
| theTolerance | for rank detection |