Class DecompositionFactory_DDRM
Contains operations related to creating and evaluating the quality of common matrix decompositions. Except in specialized situations, matrix decompositions should be instantiated from this factory instead of being directly constructed. Low level implementations are more prone to changes and new algorithms will be automatically placed here.
Several functions are also provided to evaluate the quality of a decomposition. This is provided as a way to sanity check a decomposition. Often times a significant error in a decomposition will result in a poor (larger) quality value. Typically a quality value of around 1e-15 means it is within machine precision.
-
Constructor Summary
-
Method Summary
Modifier and TypeMethodDescriptionstatic CholeskyDecomposition_F64<DMatrixRMaj>
chol
(boolean lower) Returns aCholeskyDecomposition_F64
that isn't specialized for any specific matrix size.static CholeskyDecomposition_F64<DMatrixRMaj>
chol
(int matrixSize, boolean lower) Returns aCholeskyDecomposition_F64
that has been optimized for the specified matrix size.cholLDL()
cholLDL
(int matrixSize) Returns aCholeskyDecompositionLDL_DDRM
that has been optimized for the specified matrix size.static <T extends DMatrix>
booleandecomposeSafe
(DecompositionInterface<T> decomp, T M) A simple convenience function that decomposes the matrix but automatically checks the input to make sure it is not being modified.static EigenDecomposition_F64<DMatrixRMaj>
eig
(boolean needVectors) static EigenDecomposition_F64<DMatrixRMaj>
eig
(boolean computeVectors, boolean isSymmetric) static EigenDecomposition_F64<DMatrixRMaj>
eig
(int matrixSize, boolean needVectors) Returns anEigenDecomposition
that has been optimized for the specified matrix size.static EigenDecomposition_F64<DMatrixRMaj>
eig
(int matrixSize, boolean computeVectors, boolean isSymmetric) Returns anEigenDecomposition
which is specialized for symmetric matrices or the general problem.static LUDecomposition_F64<DMatrixRMaj>
lu()
static LUDecomposition_F64<DMatrixRMaj>
lu
(int numRows, int numCol) Returns aLUDecomposition
that has been optimized for the specified matrix size.static QRDecomposition<DMatrixRMaj>
qr()
static QRDecomposition<DMatrixRMaj>
qr
(int numRows, int numCols) Returns aQRDecomposition
that has been optimized for the specified matrix size.static QRPDecomposition_F64<DMatrixRMaj>
qrp()
static QRPDecomposition_F64<DMatrixRMaj>
qrp
(int numRows, int numCols) Returns aQRPDecomposition_F64
that has been optimized for the specified matrix size.static double
quality
(DMatrixRMaj orig, DMatrixRMaj U, DMatrixRMaj W, DMatrixRMaj Vt) static double
quality
(DMatrixRMaj orig, EigenDecomposition_F64<DMatrixRMaj> eig) Computes a metric which measures the the quality of an eigen value decomposition.static double
quality
(DMatrixRMaj orig, SingularValueDecomposition<DMatrixRMaj> svd) Computes a metric which measures the the quality of a singular value decomposition.svd
(boolean needU, boolean needV, boolean compact) Returns aSingularValueDecomposition
that is NOT optimized for any specified matrix size.svd
(int numRows, int numCols, boolean needU, boolean needV, boolean compact) Returns aSingularValueDecomposition
that has been optimized for the specified matrix size.tridiagonal
(int matrixSize) Checks to see if the passed in tridiagonal decomposition is of the appropriate type for the matrix of the provided size.
-
Constructor Details
-
DecompositionFactory_DDRM
public DecompositionFactory_DDRM()
-
-
Method Details
-
chol
Returns a
CholeskyDecomposition_F64
that has been optimized for the specified matrix size.- Parameters:
matrixSize
- Number of rows and columns that the returned decomposition is optimized for.lower
- should a lower or upper triangular matrix be used. If not sure set to true.- Returns:
- A new CholeskyDecomposition.
-
chol
Returns aCholeskyDecomposition_F64
that isn't specialized for any specific matrix size.- Parameters:
lower
- should a lower or upper triangular matrix be used. If not sure set to true.- Returns:
- A new CholeskyDecomposition.
-
cholLDL
Returns a
CholeskyDecompositionLDL_DDRM
that has been optimized for the specified matrix size.- Parameters:
matrixSize
- Number of rows and columns that the returned decomposition is optimized for.- Returns:
- CholeskyLDLDecomposition_F64
-
cholLDL
-
lu
Returns a
LUDecomposition
that has been optimized for the specified matrix size.- Parameters:
numRows
- Shape of the matrix that the code should be targeted towards. Does not need to be exact.numCol
- Shape of the matrix that the code should be targeted towards. Does not need to be exact.- Returns:
- LUDecomposition
-
lu
-
svd
public static SingularValueDecomposition_F64<DMatrixRMaj> svd(int numRows, int numCols, boolean needU, boolean needV, boolean compact) Returns a
SingularValueDecomposition
that has been optimized for the specified matrix size. For improved performance only the portion of the decomposition that the user requests will be computed.- Parameters:
numRows
- Number of rows the returned decomposition is optimized for.numCols
- Number of columns that the returned decomposition is optimized for.needU
- Should it compute the U matrix. If not sure set to true.needV
- Should it compute the V matrix. If not sure set to true.compact
- Should it compute the SVD in compact form. If not sure set to false.- Returns:
- SVD
-
svd
public static SingularValueDecomposition_F64<DMatrixRMaj> svd(boolean needU, boolean needV, boolean compact) Returns aSingularValueDecomposition
that is NOT optimized for any specified matrix size.- Parameters:
needU
- Should it compute the U matrix. If not sure set to true.needV
- Should it compute the V matrix. If not sure set to true.compact
- Should it compute the SVD in compact form. If not sure set to false.- Returns:
- SVD
-
qr
Returns a
QRDecomposition
that has been optimized for the specified matrix size.- Parameters:
numRows
- Number of rows the returned decomposition is optimized for.numCols
- Number of columns that the returned decomposition is optimized for.- Returns:
- QRDecomposition
-
qr
-
qrp
Returns a
QRPDecomposition_F64
that has been optimized for the specified matrix size.- Parameters:
numRows
- Number of rows the returned decomposition is optimized for.numCols
- Number of columns that the returned decomposition is optimized for.- Returns:
- QRPDecomposition_F64
-
qrp
-
eig
Returns an
EigenDecomposition
that has been optimized for the specified matrix size. If the input matrix is symmetric within tolerance then the symmetric algorithm will be used, otherwise a general purpose eigenvalue decomposition is used.- Parameters:
matrixSize
- Number of rows and columns that the returned decomposition is optimized for.needVectors
- Should eigenvectors be computed or not. If not sure set to true.- Returns:
- A new EigenDecomposition
-
eig
-
eig
public static EigenDecomposition_F64<DMatrixRMaj> eig(int matrixSize, boolean computeVectors, boolean isSymmetric) Returns an
EigenDecomposition
which is specialized for symmetric matrices or the general problem.- Parameters:
matrixSize
- Number of rows and columns that the returned decomposition is optimized for.computeVectors
- Should it compute the eigenvectors or just eigenvalues.isSymmetric
- If true then the returned algorithm is specialized only for symmetric matrices, if false then a general purpose algorithm is returned.- Returns:
- EVD for any matrix.
-
eig
-
quality
Computes a metric which measures the the quality of a singular value decomposition. If a value is returned that is close to or smaller than 1e-15 then it is within machine precision.
SVD quality is defined as:
Quality = || A - U W VT|| / || A ||
where A is the original matrix , U W V is the decomposition, and ||A|| is the norm-f of A.- Parameters:
orig
- The original matrix which was decomposed. Not modified.svd
- The decomposition after processing 'orig'. Not modified.- Returns:
- The quality of the decomposition.
-
quality
-
quality
Computes a metric which measures the the quality of an eigen value decomposition. If a value is returned that is close to or smaller than 1e-15 then it is within machine precision.
EVD quality is defined as:
Quality = ||A*V - V*D|| / ||A*V||.- Parameters:
orig
- The original matrix. Not modified.eig
- EVD of the original matrix. Not modified.- Returns:
- The quality of the decomposition.
-
tridiagonal
Checks to see if the passed in tridiagonal decomposition is of the appropriate type for the matrix of the provided size. Returns the same instance or a new instance.- Parameters:
matrixSize
- Number of rows and columns that the returned decomposition is optimized for.
-
decomposeSafe
A simple convenience function that decomposes the matrix but automatically checks the input to make sure it is not being modified.- Type Parameters:
T
- Matrix type.- Parameters:
decomp
- Decomposition which is being wrappedM
- The matrix being decomposed.- Returns:
- If the decomposition was successful or not.
-