MAGMA  2.0.2
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_queue_t queue, 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...
 

Detailed Description

Function Documentation

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_queue_t  queue,
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.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 0.
[in,out]dADOUBLE 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 orthogonal matrix Q as a product of elementary reflectors (see Further Details).
[in]lddaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[out]dtauDOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).
dwork(workspace) DOUBLE PRECISION array, dimension (N)
[in]queuemagma_queue_t Queue to execute in.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value

Further Details

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.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 0.
[in,out]dADOUBLE PRECISION array, dimension (LDDA,N) On entry, the m by n matrix A. On exit, the orthogonal 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 orthogonal matrix Q as a product of elementary reflectors (see Further Details).
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]dtauDOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).
[out]dTDOUBLE 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]ddADOUBLE 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]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value

Further Details

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.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. N >= 0.
[in,out]dADOUBLE PRECISION array, dimension (LDDA,N) On entry, the m by n matrix A. On exit, the orthogonal 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 orthogonal matrix Q as a product of elementary reflectors (see Further Details).
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]dtauDOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).
[out]dTDOUBLE 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]ddADOUBLE 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]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value

Further Details

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_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.

Parameters
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]nINTEGER The number of columns of the matrix A. 0 <= N <= min(M, 128).
[in,out]dADOUBLE PRECISION array, dimension (LDDA,N) On entry, the m by n matrix A. On exit, the orthogonal 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 orthogonal matrix Q as a product of elementary reflectors (see Further Details).
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]dtauDOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).
[out]dTDOUBLE 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]ddADOUBLE 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]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value

Further Details

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).