![]() |
MAGMA
2.0.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
magma_int_t | magma_dgeqr2x2_gpu (magma_int_t m, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dtau, magmaDouble_ptr dT, magmaDouble_ptr ddA, magmaDouble_ptr dwork, magma_int_t *info) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R. More... | |
magma_int_t | magma_dgeqr2x3_gpu (magma_int_t m, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dtau, magmaDouble_ptr dT, magmaDouble_ptr ddA, magmaDouble_ptr dwork, magma_int_t *info) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R. More... | |
magma_int_t | magma_dgeqr2x_gpu (magma_int_t m, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dtau, magmaDouble_ptr dT, magmaDouble_ptr ddA, magmaDouble_ptr dwork, magma_int_t *info) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R. More... | |
magma_int_t | magma_dgeqr2_gpu (magma_int_t m, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dtau, magmaDouble_ptr dwork, magma_int_t *info) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R using the non-blocking Householder QR. More... | |
magma_int_t | magma_dgeqr2_batched (magma_int_t m, magma_int_t n, double **dA_array, magma_int_t ldda, double **dtau_array, magma_int_t *info_array, magma_int_t batchCount, magma_queue_t queue) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R. More... | |
void | dgeqrf_copy_upper_batched (magma_int_t n, magma_int_t nb, double **dV_array, magma_int_t lddv, double **dR_array, magma_int_t lddr, magma_int_t batchCount, magma_queue_t queue) |
These are internal routines that might have many assumption. More... | |
magma_int_t | magma_dgeqr2x4_gpu (magma_int_t m, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dtau, magmaDouble_ptr dT, magmaDouble_ptr ddA, magmaDouble_ptr dwork, magma_queue_t queue, magma_int_t *info) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R. More... | |
void dgeqrf_copy_upper_batched | ( | magma_int_t | n, |
magma_int_t | nb, | ||
double ** | dV_array, | ||
magma_int_t | lddv, | ||
double ** | dR_array, | ||
magma_int_t | lddr, | ||
magma_int_t | batchCount, | ||
magma_queue_t | queue | ||
) |
These are internal routines that might have many assumption.
They are used in dgeqrf_batched.cpp
Copy part of the data in dV to dR
[in] | n | INTEGER The order of the matrix . N >= 0. |
[in] | nb | INTEGER Tile size in matrix. nb <= N. |
[in] | dV_array | Array of pointers, dimension (batchCount). Each is a DOUBLE_PRECISION array on the GPU, dimension (LDDA,N). |
[in] | lddv | INTEGER The leading dimension of each array V. LDDV >= max(1,N). |
[in,out] | dR_array | Array of pointers, dimension (batchCount). Each is a DOUBLE_PRECISION array on the GPU, dimension (LDDR,N). |
[in] | lddr | INTEGER The leading dimension of each array R. LDDR >= max(1,N). |
[in] | batchCount | INTEGER The number of matrices to operate on. |
[in] | queue | magma_queue_t Queue to execute in. |
magma_int_t magma_dgeqr2_batched | ( | magma_int_t | m, |
magma_int_t | n, | ||
double ** | dA_array, | ||
magma_int_t | ldda, | ||
double ** | dtau_array, | ||
magma_int_t * | info_array, | ||
magma_int_t | batchCount, | ||
magma_queue_t | queue | ||
) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R.
This version implements the right-looking QR with non-blocking.
[in] | m | INTEGER The number of rows of the matrix A. M >= 0. |
[in] | n | INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | dA_array | Array of pointers, dimension (batchCount). Each is a DOUBLE_PRECISION array on the GPU, dimension (LDDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors (see Further Details). |
[in] | ldda | INTEGER The leading dimension of the array dA. LDDA >= max(1,M). To benefit from coalescent memory accesses LDDA must be divisible by 16. |
[out] | dtau_array | Array of pointers, dimension (batchCount). Each is a DOUBLE_PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | info_array | Array of INTEGERs, dimension (batchCount), for corresponding matrices.
|
[in] | batchCount | INTEGER The number of matrices to operate on. |
[in] | queue | magma_queue_t Queue to execute in. |
The matrix Q is represented as a product of elementary reflectors
Q = H(1) H(2) . . . H(k), where k = min(m,n).
Each H(i) has the form
H(i) = I - tau * v * v'
where tau is a real scalar, and v is a real vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
magma_int_t magma_dgeqr2_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magmaDouble_ptr | dA, | ||
magma_int_t | ldda, | ||
magmaDouble_ptr | dtau, | ||
magmaDouble_ptr | dwork, | ||
magma_int_t * | info | ||
) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R using the non-blocking Householder QR.
[in] | m | INTEGER The number of rows of the matrix A. M >= 0. |
[in] | n | INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | dA | DOUBLE PRECISION array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the elements on and above the diagonal of the array contain the min(m,n) by n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | dtau | DOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
dwork | (workspace) DOUBLE_PRECISION array, dimension (N) | |
[out] | info | INTEGER
|
The matrix Q is represented as a product of elementary reflectors
Q = H(1) H(2) . . . H(k), where k = min(m,n).
Each H(i) has the form
H(i) = I - tau * v * v**H
where tau is a real scalar, and v is a real vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
magma_int_t magma_dgeqr2x2_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magmaDouble_ptr | dA, | ||
magma_int_t | ldda, | ||
magmaDouble_ptr | dtau, | ||
magmaDouble_ptr | dT, | ||
magmaDouble_ptr | ddA, | ||
magmaDouble_ptr | dwork, | ||
magma_int_t * | info | ||
) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R.
This expert routine requires two more arguments than the standard dgeqr2, namely, dT and ddA, explained below. The storage for A is also not as in the LAPACK's dgeqr2 routine (see below).
The first is used to output the triangular n x n factor T of the block reflector used in the factorization. The second holds the diagonal nxn blocks of A, i.e., the diagonal submatrices of R. This routine implements the left looking QR.
[in] | m | INTEGER The number of rows of the matrix A. M >= 0. |
[in] | n | INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | dA | DOUBLE_PRECISION array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the unitary matrix Q as a product of elementary reflectors (see Further Details). the elements on and above the diagonal of the array contain the min(m,n) by n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | dtau | DOUBLE_PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | dT | DOUBLE_PRECISION array, dimension N x N. Stores the triangular N x N factor T of the block reflector used in the factorization. The lower triangular part is 0. |
[out] | ddA | DOUBLE_PRECISION array, dimension N x N. Stores the elements of the upper N x N diagonal block of A. LAPACK stores this array in A. There are 0s below the diagonal. |
dwork | (workspace) DOUBLE_PRECISION array, dimension (3 N) | |
[out] | info | INTEGER
|
The matrix Q is represented as a product of elementary reflectors
Q = H(1) H(2) . . . H(k), where k = min(m,n).
Each H(i) has the form
H(i) = I - tau * v * v'
where tau is a real scalar, and v is a real vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
magma_int_t magma_dgeqr2x3_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magmaDouble_ptr | dA, | ||
magma_int_t | ldda, | ||
magmaDouble_ptr | dtau, | ||
magmaDouble_ptr | dT, | ||
magmaDouble_ptr | ddA, | ||
magmaDouble_ptr | dwork, | ||
magma_int_t * | info | ||
) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R.
This expert routine requires two more arguments than the standard dgeqr2, namely, dT and ddA, explained below. The storage for A is also not as in the LAPACK's dgeqr2 routine (see below).
The first is used to output the triangular n x n factor T of the block reflector used in the factorization. The second holds the diagonal nxn blocks of A, i.e., the diagonal submatrices of R. This routine implements the left looking QR.
This version adds internal blocking.
[in] | m | INTEGER The number of rows of the matrix A. M >= 0. |
[in] | n | INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | dA | DOUBLE_PRECISION array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the unitary matrix Q as a product of elementary reflectors (see Further Details). the elements on and above the diagonal of the array contain the min(m,n) by n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | dtau | DOUBLE_PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | dT | DOUBLE_PRECISION array, dimension N x N. Stores the triangular N x N factor T of the block reflector used in the factorization. The lower triangular part is 0. |
[out] | ddA | DOUBLE_PRECISION array, dimension N x N. Stores the elements of the upper N x N diagonal block of A. LAPACK stores this array in A. There are 0s below the diagonal. |
dwork | (workspace) DOUBLE_PRECISION array, dimension (3 N) | |
[out] | info | INTEGER
|
The matrix Q is represented as a product of elementary reflectors
Q = H(1) H(2) . . . H(k), where k = min(m,n).
Each H(i) has the form
H(i) = I - tau * v * v'
where tau is a real scalar, and v is a real vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
magma_int_t magma_dgeqr2x4_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magmaDouble_ptr | dA, | ||
magma_int_t | ldda, | ||
magmaDouble_ptr | dtau, | ||
magmaDouble_ptr | dT, | ||
magmaDouble_ptr | ddA, | ||
magmaDouble_ptr | dwork, | ||
magma_queue_t | queue, | ||
magma_int_t * | info | ||
) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R.
This expert routine requires two more arguments than the standard dgeqr2, namely, dT and ddA, explained below. The storage for A is also not as in the LAPACK's dgeqr2 routine (see below).
The first is used to output the triangular n x n factor T of the block reflector used in the factorization. The second holds the diagonal nxn blocks of A, i.e., the diagonal submatrices of R. This routine implements the left looking QR.
This version adds internal blocking.
[in] | m | INTEGER The number of rows of the matrix A. M >= 0. |
[in] | n | INTEGER The number of columns of the matrix A. N >= 0. |
[in,out] | dA | DOUBLE_PRECISION array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the unitary matrix Q as a product of elementary reflectors (see Further Details). the elements on and above the diagonal of the array contain the min(m,n) by n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | dtau | DOUBLE_PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | dT | DOUBLE_PRECISION array, dimension N x N. Stores the triangular N x N factor T of the block reflector used in the factorization. The lower triangular part is 0. |
[out] | ddA | DOUBLE_PRECISION array, dimension N x N. Stores the elements of the upper N x N diagonal block of A. LAPACK stores this array in A. There are 0s below the diagonal. |
dwork | (workspace) DOUBLE_PRECISION array, dimension (3 N) | |
[out] | info | INTEGER
|
[in] | queue | magma_queue_t Queue to execute in. |
The matrix Q is represented as a product of elementary reflectors
Q = H(1) H(2) . . . H(k), where k = min(m,n).
Each H(i) has the form
H(i) = I - tau * v * v**H
where tau is a real scalar, and v is a real vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).
magma_int_t magma_dgeqr2x_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magmaDouble_ptr | dA, | ||
magma_int_t | ldda, | ||
magmaDouble_ptr | dtau, | ||
magmaDouble_ptr | dT, | ||
magmaDouble_ptr | ddA, | ||
magmaDouble_ptr | dwork, | ||
magma_int_t * | info | ||
) |
DGEQR2 computes a QR factorization of a real m by n matrix A: A = Q * R.
This expert routine requires two more arguments than the standard dgeqr2, namely, dT and ddA, explained below. The storage for A is also not as in the LAPACK's dgeqr2 routine (see below).
The first is used to output the triangular n x n factor T of the block reflector used in the factorization. The second holds the diagonal nxn blocks of A, i.e., the diagonal submatrices of R.
This version implements the right-looking QR. A hard-coded requirement for N is to be <= min(M, 128). For larger N one should use a blocking QR version.
[in] | m | INTEGER The number of rows of the matrix A. M >= 0. |
[in] | n | INTEGER The number of columns of the matrix A. 0 <= N <= min(M, 128). |
[in,out] | dA | DOUBLE_PRECISION array, dimension (LDA,N) On entry, the m by n matrix A. On exit, the unitary matrix Q as a product of elementary reflectors (see Further Details). the elements on and above the diagonal of the array contain the min(m,n) by n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the unitary matrix Q as a product of elementary reflectors (see Further Details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[out] | dtau | DOUBLE_PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details). |
[out] | dT | DOUBLE_PRECISION array, dimension N x N. Stores the triangular N x N factor T of the block reflector used in the factorization. The lower triangular part is 0. |
[out] | ddA | DOUBLE_PRECISION array, dimension N x N. Stores the elements of the upper N x N diagonal block of A. LAPACK stores this array in A. There are 0s below the diagonal. |
dwork | (workspace) DOUBLE_PRECISION array, dimension (N) | |
[out] | info | INTEGER
|
The matrix Q is represented as a product of elementary reflectors
Q = H(1) H(2) . . . H(k), where k = min(m,n).
Each H(i) has the form
H(i) = I - tau * v * v'
where tau is a real scalar, and v is a real vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and tau in TAU(i).