namespace solver¶
// namespaces namespace xf::solver::internal namespace xf::solver::internalgetrf
geqrf¶
geqrf overload (1)¶
#include "MatrixDecomposition/geqrf.hpp"
template < typename T, int NRMAX, int NCMAX, int NCU > int geqrf ( int m, int n, T* A, int lda, T* tau )
This function computes QR decomposition of matrix \(A\)
where \(A\) is a dense matrix of size \(m \times n\) , \(Q\) is a \(m \times n\) matrix with orthonormal columns, and \(R\) is an upper triangular matrix.
The maximum matrix size supported in FPGA is templated by NRMAX and NCMAX.
Parameters:
T | data type (support float and double) |
NRMAX | maximum number of rows of input matrix |
NCMAX | maximum number of columns of input matrix |
NCU | number of computation unit, NCU = n will be n times faster than NCU = 1, and cost n times resource. |
m | number of rows of matrix A |
n | number of cols of matrix A |
A | input matrix of size \(m \times lda\) , and overwritten by the output triangular R matrix and min(m,n) elementary reflectors |
lda | leading dimension of matrix A |
tau | scalar factors for elementary reflectors |
gesvdj¶
#include "MatrixDecomposition/gesvdj.hpp"
template < typename T, int NMAX, int NCU > void gesvdj ( int m, T* A, int lda, T* S, T* U, int ldu, T* V, int ldv, int& info )
Symmetric Matrix Jacobi based Singular Value Decomposition (GESVDJ) .
where \(A\) is a dense symmetric matrix of size \(m \times m\) , \(U\) and \(V\) are \(m \times m\) matrix with orthonormal columns, and \(\Sigma\) is diagonal matrix.
The maximum matrix size supported in FPGA is templated by NMAX.
Parameters:
T | data type (support float and double). |
NMAX | maximum number of rows/columns of input matrix |
NCU | number of computation unit |
m | number of rows/cols of matrix A |
A | input matrix of size \(m \times m\) |
S | decomposed diagonal singular matrix of size \(m \times m\) |
U | left U matrix of SVD |
V | right V matrix of SVD |
lda | leading dimension of matrix A |
ldu | leading dimension of matrix U |
ldv | leading dimension of matrix V |
info | output info (unused) |
gesvj¶
#include "MatrixDecomposition/gesvj.hpp"
template < typename T, int NRMAX, int NCMAX, int MCU, int NCU > void gesvj ( int m, int n, T* A, T* U, T* S, T* V )
This function implements singular value decomposition of matrix A using one-sided Jacobi algorihtm.
where \(A\) is a dense matrix of size \(m \times n\) , \(U\) is \(m \times m\) matrix with orthonormal columns, \(V\) is \(n \times n\) matrix with orthonormal columns, and \(\Sigma\) is diagonal matrix.
The maximum matrix size supported in FPGA is templated by NCMAX, NRMAX.
Parameters:
T | |
: | the data type of gesvj |
NRMAX | maximum number of rows of input matrix |
NCMAX | maximum number of columns of input matrix |
MCU | number of computation unit of M |
NCU | number of computation unit of N |
m | number of rows of matrix A |
n | number of cols of matrix A |
A | input matrix of size \(m \times n\) |
S | decomposed diagonal singular matrix of size n |
U | left U matrix of SVD of size \(m \times m\) |
V | right V matrix of SVD \(n \times n\) |
getrf¶
#include "MatrixDecomposition/getrf.hpp"
template < typename T, int NMAX, int NCU > void getrf ( int n, T* A, int lda, int* ipiv, int& info )
This function computes the LU decomposition (with partial pivoting) of matrix \(A\)
where \(P\) is a permutation matrix, \(A\) is a dense matrix of size \(n \times n\) , \(L\) is a lower triangular matrix with unit diagonal, and \(U\) is an upper triangular matrix. This function does not implement pivoting.
The maximum matrix size supported in FPGA is templated by NMAX.
Parameters:
T | data type (support float and double) |
NMAX | maximum number of rows/columns of input matrix |
NCU | number of computation unit |
n | number of rows/cols of matrix A |
A | input matrix, and overwritten by the output upper and lower triangular matrix |
lda | leading dimention of input matrix A |
pivot | indices, row i of matrix A is stored in row[i] |
info | output info (unused) |
getrf_nopivot¶
#include "MatrixDecomposition/getrf_nopivot.hpp"
template < typename T, int NMAX, int NCU > void getrf_nopivot ( int n, T* A, int lda, int& info )
This function computes the LU decomposition (without pivoting) of matrix \(A\)
where \(A\) is a dense matrix of size \(n \times n\) , \(L\) is a lower triangular matrix with unit diagonal, and \(U\) is an upper triangular matrix. This function does not implement pivoting.
The maximum matrix size supported in FPGA is templated by NMAX.
Parameters:
T | data type (support float and double) |
NMAX | maximum number of rows/cols of input matrix |
NCU | number of computation unit |
n | number of rows/cols of matrix A |
A | input matrix |
lda | leading dimention of input matrix A |
info | output info (unused) |
potrf¶
#include "MatrixDecomposition/potrf.hpp"
template < typename T, int NMAX, int NCU > void potrf ( int m, T* A, int lda, int& info )
This function computes the Cholesky decomposition of matrix \(A\)
where \(A\) is a dense symmetric positive-definite matrix of size \(m \times m\) , \(L\) is a lower triangular matrix, and \({L}^T\) is the transposed matrix of \(L\) .
The maximum matrix size supported in FPGA is templated by NMAX.
Parameters:
T | data type (support float and double) |
NMAX | maximum number of rows/columns of input matrix |
NCU | number of computation unit |
m | number of rows/cols of matrix A |
A | input matrix of size \(m \times m\) |
lda | leading dimention of input matrix A |
info | output info (unused) |