![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
magma_int_t | magma_clabrd_gpu (magma_int_t m, magma_int_t n, magma_int_t nb, magmaFloatComplex *A, magma_int_t lda, magmaFloatComplex_ptr dA, magma_int_t ldda, float *d, float *e, magmaFloatComplex *tauq, magmaFloatComplex *taup, magmaFloatComplex *X, magma_int_t ldx, magmaFloatComplex_ptr dX, magma_int_t lddx, magmaFloatComplex *Y, magma_int_t ldy, magmaFloatComplex_ptr dY, magma_int_t lddy, magmaFloatComplex *work, magma_int_t lwork, magma_queue_t queue) |
CLABRD reduces the first NB rows and columns of a complex general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A. | |
magma_int_t | magma_dlabrd_gpu (magma_int_t m, magma_int_t n, magma_int_t nb, double *A, magma_int_t lda, magmaDouble_ptr dA, magma_int_t ldda, double *d, double *e, double *tauq, double *taup, double *X, magma_int_t ldx, magmaDouble_ptr dX, magma_int_t lddx, double *Y, magma_int_t ldy, magmaDouble_ptr dY, magma_int_t lddy, double *work, magma_int_t lwork, magma_queue_t queue) |
DLABRD reduces the first NB rows and columns of a real general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A. | |
magma_int_t | magma_slabrd_gpu (magma_int_t m, magma_int_t n, magma_int_t nb, float *A, magma_int_t lda, magmaFloat_ptr dA, magma_int_t ldda, float *d, float *e, float *tauq, float *taup, float *X, magma_int_t ldx, magmaFloat_ptr dX, magma_int_t lddx, float *Y, magma_int_t ldy, magmaFloat_ptr dY, magma_int_t lddy, float *work, magma_int_t lwork, magma_queue_t queue) |
SLABRD reduces the first NB rows and columns of a real general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A. | |
magma_int_t | magma_zlabrd_gpu (magma_int_t m, magma_int_t n, magma_int_t nb, magmaDoubleComplex *A, magma_int_t lda, magmaDoubleComplex_ptr dA, magma_int_t ldda, double *d, double *e, magmaDoubleComplex *tauq, magmaDoubleComplex *taup, magmaDoubleComplex *X, magma_int_t ldx, magmaDoubleComplex_ptr dX, magma_int_t lddx, magmaDoubleComplex *Y, magma_int_t ldy, magmaDoubleComplex_ptr dY, magma_int_t lddy, magmaDoubleComplex *work, magma_int_t lwork, magma_queue_t queue) |
ZLABRD reduces the first NB rows and columns of a complex general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A. | |
magma_int_t magma_clabrd_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magma_int_t | nb, | ||
magmaFloatComplex * | A, | ||
magma_int_t | lda, | ||
magmaFloatComplex_ptr | dA, | ||
magma_int_t | ldda, | ||
float * | d, | ||
float * | e, | ||
magmaFloatComplex * | tauq, | ||
magmaFloatComplex * | taup, | ||
magmaFloatComplex * | X, | ||
magma_int_t | ldx, | ||
magmaFloatComplex_ptr | dX, | ||
magma_int_t | lddx, | ||
magmaFloatComplex * | Y, | ||
magma_int_t | ldy, | ||
magmaFloatComplex_ptr | dY, | ||
magma_int_t | lddy, | ||
magmaFloatComplex * | work, | ||
magma_int_t | lwork, | ||
magma_queue_t | queue ) |
CLABRD reduces the first NB rows and columns of a complex general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A.
If m >= n, A is reduced to upper bidiagonal form; if m < n, to lower bidiagonal form.
This is an auxiliary routine called by CGEBRD.
[in] | m | INTEGER The number of rows in the matrix A. |
[in] | n | INTEGER The number of columns in the matrix A. |
[in] | nb | INTEGER The number of leading rows and columns of A to be reduced. |
[in,out] | A | COMPLEX array, dimension (LDA,N) On entry, the m by n general matrix to be reduced. On exit, the first NB rows and columns of the matrix are overwritten; the rest of the array is unchanged. If m >= n, elements on and below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors; and elements above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. If m < n, elements below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors, and elements on and above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. See Further Details. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[in,out] | dA | COMPLEX array, dimension (LDDA,N) Copy of A on GPU. |
[in] | ldda | INTEGER The leading dimension of the array dA. LDDA >= max(1,M). |
[out] | d | COMPLEX array, dimension (NB) The diagonal elements of the first NB rows and columns of the reduced matrix. D(i) = A(i,i). |
[out] | e | COMPLEX array, dimension (NB) The off-diagonal elements of the first NB rows and columns of the reduced matrix. |
[out] | tauq | COMPLEX array dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix Q. See Further Details. |
[out] | taup | COMPLEX array, dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix P. See Further Details. |
[out] | X | COMPLEX array, dimension (LDX,NB) The m-by-nb matrix X required to update the unreduced part of A. |
[in] | ldx | INTEGER The leading dimension of the array X. LDX >= M. |
[out] | dX | COMPLEX array, dimension (LDDX,NB) Copy of X on GPU. |
[in] | lddx | INTEGER The leading dimension of the array dX. LDDX >= M. |
[out] | Y | COMPLEX array, dimension (LDY,NB) The n-by-nb matrix Y required to update the unreduced part of A. |
[in] | ldy | INTEGER The leading dimension of the array Y. LDY >= N. |
[out] | dY | COMPLEX array, dimension (LDDY,NB) Copy of Y on GPU. |
[in] | lddy | INTEGER The leading dimension of the array dY. LDDY >= N. |
work | COMPLEX array, dimension (LWORK) Workspace. | |
[in] | lwork | INTEGER The dimension of the array WORK. LWORK >= max( M, N ). |
[in] | queue | magma_queue_t Queue to execute in. |
The matrices Q and P are represented as products of elementary reflectors:
Q = H(1) H(2) . . . H(nb) and P = G(1) G(2) . . . G(nb)
Each H(i) and G(i) has the form:
H(i) = I - tauq * v * v' and G(i) = I - taup * u * u'
where tauq and taup are complex scalars, and v and u are complex vectors.
If m >= n, v(1:i-1) = 0, v(i) = 1, and v(i:m) is stored on exit in A(i:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
If m < n, v(1:i) = 0, v(i+1) = 1, and v(i+1:m) is stored on exit in A(i+2:m,i); u(1:i-1) = 0, u(i) = 1, and u(i:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
The elements of the vectors v and u together form the m-by-nb matrix V and the nb-by-n matrix U' which are needed, with X and Y, to apply the transformation to the unreduced part of the matrix, using a block update of the form: A := A - V*Y' - X*U'.
The contents of A on exit are illustrated by the following examples with nb = 2:
m = 6 and n = 5 (m > n): m = 5 and n = 6 (m < n): ( 1 1 u1 u1 u1 ) ( 1 u1 u1 u1 u1 u1 ) ( v1 1 1 u2 u2 ) ( 1 1 u2 u2 u2 u2 ) ( v1 v2 a a a ) ( v1 1 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a )
where a denotes an element of the original matrix which is unchanged, vi denotes an element of the vector defining H(i), and ui an element of the vector defining G(i).
magma_int_t magma_dlabrd_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magma_int_t | nb, | ||
double * | A, | ||
magma_int_t | lda, | ||
magmaDouble_ptr | dA, | ||
magma_int_t | ldda, | ||
double * | d, | ||
double * | e, | ||
double * | tauq, | ||
double * | taup, | ||
double * | X, | ||
magma_int_t | ldx, | ||
magmaDouble_ptr | dX, | ||
magma_int_t | lddx, | ||
double * | Y, | ||
magma_int_t | ldy, | ||
magmaDouble_ptr | dY, | ||
magma_int_t | lddy, | ||
double * | work, | ||
magma_int_t | lwork, | ||
magma_queue_t | queue ) |
DLABRD reduces the first NB rows and columns of a real general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A.
If m >= n, A is reduced to upper bidiagonal form; if m < n, to lower bidiagonal form.
This is an auxiliary routine called by DGEBRD.
[in] | m | INTEGER The number of rows in the matrix A. |
[in] | n | INTEGER The number of columns in the matrix A. |
[in] | nb | INTEGER The number of leading rows and columns of A to be reduced. |
[in,out] | A | DOUBLE PRECISION array, dimension (LDA,N) On entry, the m by n general matrix to be reduced. On exit, the first NB rows and columns of the matrix are overwritten; the rest of the array is unchanged. If m >= n, elements on and below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors; and elements above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. If m < n, elements below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors, and elements on and above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. See Further Details. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[in,out] | dA | DOUBLE PRECISION array, dimension (LDDA,N) Copy of A on GPU. |
[in] | ldda | INTEGER The leading dimension of the array dA. LDDA >= max(1,M). |
[out] | d | DOUBLE PRECISION array, dimension (NB) The diagonal elements of the first NB rows and columns of the reduced matrix. D(i) = A(i,i). |
[out] | e | DOUBLE PRECISION array, dimension (NB) The off-diagonal elements of the first NB rows and columns of the reduced matrix. |
[out] | tauq | DOUBLE PRECISION array dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix Q. See Further Details. |
[out] | taup | DOUBLE PRECISION array, dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix P. See Further Details. |
[out] | X | DOUBLE PRECISION array, dimension (LDX,NB) The m-by-nb matrix X required to update the unreduced part of A. |
[in] | ldx | INTEGER The leading dimension of the array X. LDX >= M. |
[out] | dX | DOUBLE PRECISION array, dimension (LDDX,NB) Copy of X on GPU. |
[in] | lddx | INTEGER The leading dimension of the array dX. LDDX >= M. |
[out] | Y | DOUBLE PRECISION array, dimension (LDY,NB) The n-by-nb matrix Y required to update the unreduced part of A. |
[in] | ldy | INTEGER The leading dimension of the array Y. LDY >= N. |
[out] | dY | DOUBLE PRECISION array, dimension (LDDY,NB) Copy of Y on GPU. |
[in] | lddy | INTEGER The leading dimension of the array dY. LDDY >= N. |
work | DOUBLE PRECISION array, dimension (LWORK) Workspace. | |
[in] | lwork | INTEGER The dimension of the array WORK. LWORK >= max( M, N ). |
[in] | queue | magma_queue_t Queue to execute in. |
The matrices Q and P are represented as products of elementary reflectors:
Q = H(1) H(2) . . . H(nb) and P = G(1) G(2) . . . G(nb)
Each H(i) and G(i) has the form:
H(i) = I - tauq * v * v' and G(i) = I - taup * u * u'
where tauq and taup are real scalars, and v and u are real vectors.
If m >= n, v(1:i-1) = 0, v(i) = 1, and v(i:m) is stored on exit in A(i:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
If m < n, v(1:i) = 0, v(i+1) = 1, and v(i+1:m) is stored on exit in A(i+2:m,i); u(1:i-1) = 0, u(i) = 1, and u(i:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
The elements of the vectors v and u together form the m-by-nb matrix V and the nb-by-n matrix U' which are needed, with X and Y, to apply the transformation to the unreduced part of the matrix, using a block update of the form: A := A - V*Y' - X*U'.
The contents of A on exit are illustrated by the following examples with nb = 2:
m = 6 and n = 5 (m > n): m = 5 and n = 6 (m < n): ( 1 1 u1 u1 u1 ) ( 1 u1 u1 u1 u1 u1 ) ( v1 1 1 u2 u2 ) ( 1 1 u2 u2 u2 u2 ) ( v1 v2 a a a ) ( v1 1 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a )
where a denotes an element of the original matrix which is unchanged, vi denotes an element of the vector defining H(i), and ui an element of the vector defining G(i).
magma_int_t magma_slabrd_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magma_int_t | nb, | ||
float * | A, | ||
magma_int_t | lda, | ||
magmaFloat_ptr | dA, | ||
magma_int_t | ldda, | ||
float * | d, | ||
float * | e, | ||
float * | tauq, | ||
float * | taup, | ||
float * | X, | ||
magma_int_t | ldx, | ||
magmaFloat_ptr | dX, | ||
magma_int_t | lddx, | ||
float * | Y, | ||
magma_int_t | ldy, | ||
magmaFloat_ptr | dY, | ||
magma_int_t | lddy, | ||
float * | work, | ||
magma_int_t | lwork, | ||
magma_queue_t | queue ) |
SLABRD reduces the first NB rows and columns of a real general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A.
If m >= n, A is reduced to upper bidiagonal form; if m < n, to lower bidiagonal form.
This is an auxiliary routine called by SGEBRD.
[in] | m | INTEGER The number of rows in the matrix A. |
[in] | n | INTEGER The number of columns in the matrix A. |
[in] | nb | INTEGER The number of leading rows and columns of A to be reduced. |
[in,out] | A | REAL array, dimension (LDA,N) On entry, the m by n general matrix to be reduced. On exit, the first NB rows and columns of the matrix are overwritten; the rest of the array is unchanged. If m >= n, elements on and below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors; and elements above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. If m < n, elements below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors, and elements on and above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. See Further Details. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[in,out] | dA | REAL array, dimension (LDDA,N) Copy of A on GPU. |
[in] | ldda | INTEGER The leading dimension of the array dA. LDDA >= max(1,M). |
[out] | d | REAL array, dimension (NB) The diagonal elements of the first NB rows and columns of the reduced matrix. D(i) = A(i,i). |
[out] | e | REAL array, dimension (NB) The off-diagonal elements of the first NB rows and columns of the reduced matrix. |
[out] | tauq | REAL array dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix Q. See Further Details. |
[out] | taup | REAL array, dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix P. See Further Details. |
[out] | X | REAL array, dimension (LDX,NB) The m-by-nb matrix X required to update the unreduced part of A. |
[in] | ldx | INTEGER The leading dimension of the array X. LDX >= M. |
[out] | dX | REAL array, dimension (LDDX,NB) Copy of X on GPU. |
[in] | lddx | INTEGER The leading dimension of the array dX. LDDX >= M. |
[out] | Y | REAL array, dimension (LDY,NB) The n-by-nb matrix Y required to update the unreduced part of A. |
[in] | ldy | INTEGER The leading dimension of the array Y. LDY >= N. |
[out] | dY | REAL array, dimension (LDDY,NB) Copy of Y on GPU. |
[in] | lddy | INTEGER The leading dimension of the array dY. LDDY >= N. |
work | REAL array, dimension (LWORK) Workspace. | |
[in] | lwork | INTEGER The dimension of the array WORK. LWORK >= max( M, N ). |
[in] | queue | magma_queue_t Queue to execute in. |
The matrices Q and P are represented as products of elementary reflectors:
Q = H(1) H(2) . . . H(nb) and P = G(1) G(2) . . . G(nb)
Each H(i) and G(i) has the form:
H(i) = I - tauq * v * v' and G(i) = I - taup * u * u'
where tauq and taup are real scalars, and v and u are real vectors.
If m >= n, v(1:i-1) = 0, v(i) = 1, and v(i:m) is stored on exit in A(i:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
If m < n, v(1:i) = 0, v(i+1) = 1, and v(i+1:m) is stored on exit in A(i+2:m,i); u(1:i-1) = 0, u(i) = 1, and u(i:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
The elements of the vectors v and u together form the m-by-nb matrix V and the nb-by-n matrix U' which are needed, with X and Y, to apply the transformation to the unreduced part of the matrix, using a block update of the form: A := A - V*Y' - X*U'.
The contents of A on exit are illustrated by the following examples with nb = 2:
m = 6 and n = 5 (m > n): m = 5 and n = 6 (m < n): ( 1 1 u1 u1 u1 ) ( 1 u1 u1 u1 u1 u1 ) ( v1 1 1 u2 u2 ) ( 1 1 u2 u2 u2 u2 ) ( v1 v2 a a a ) ( v1 1 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a )
where a denotes an element of the original matrix which is unchanged, vi denotes an element of the vector defining H(i), and ui an element of the vector defining G(i).
magma_int_t magma_zlabrd_gpu | ( | magma_int_t | m, |
magma_int_t | n, | ||
magma_int_t | nb, | ||
magmaDoubleComplex * | A, | ||
magma_int_t | lda, | ||
magmaDoubleComplex_ptr | dA, | ||
magma_int_t | ldda, | ||
double * | d, | ||
double * | e, | ||
magmaDoubleComplex * | tauq, | ||
magmaDoubleComplex * | taup, | ||
magmaDoubleComplex * | X, | ||
magma_int_t | ldx, | ||
magmaDoubleComplex_ptr | dX, | ||
magma_int_t | lddx, | ||
magmaDoubleComplex * | Y, | ||
magma_int_t | ldy, | ||
magmaDoubleComplex_ptr | dY, | ||
magma_int_t | lddy, | ||
magmaDoubleComplex * | work, | ||
magma_int_t | lwork, | ||
magma_queue_t | queue ) |
ZLABRD reduces the first NB rows and columns of a complex general m by n matrix A to upper or lower bidiagonal form by an orthogonal transformation Q' * A * P, and returns the matrices X and Y which are needed to apply the transformation to the unreduced part of A.
If m >= n, A is reduced to upper bidiagonal form; if m < n, to lower bidiagonal form.
This is an auxiliary routine called by ZGEBRD.
[in] | m | INTEGER The number of rows in the matrix A. |
[in] | n | INTEGER The number of columns in the matrix A. |
[in] | nb | INTEGER The number of leading rows and columns of A to be reduced. |
[in,out] | A | COMPLEX_16 array, dimension (LDA,N) On entry, the m by n general matrix to be reduced. On exit, the first NB rows and columns of the matrix are overwritten; the rest of the array is unchanged. If m >= n, elements on and below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors; and elements above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. If m < n, elements below the diagonal in the first NB columns, with the array TAUQ, represent the orthogonal matrix Q as a product of elementary reflectors, and elements on and above the diagonal in the first NB rows, with the array TAUP, represent the orthogonal matrix P as a product of elementary reflectors. See Further Details. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,M). |
[in,out] | dA | COMPLEX_16 array, dimension (LDDA,N) Copy of A on GPU. |
[in] | ldda | INTEGER The leading dimension of the array dA. LDDA >= max(1,M). |
[out] | d | COMPLEX_16 array, dimension (NB) The diagonal elements of the first NB rows and columns of the reduced matrix. D(i) = A(i,i). |
[out] | e | COMPLEX_16 array, dimension (NB) The off-diagonal elements of the first NB rows and columns of the reduced matrix. |
[out] | tauq | COMPLEX_16 array dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix Q. See Further Details. |
[out] | taup | COMPLEX_16 array, dimension (NB) The scalar factors of the elementary reflectors which represent the orthogonal matrix P. See Further Details. |
[out] | X | COMPLEX_16 array, dimension (LDX,NB) The m-by-nb matrix X required to update the unreduced part of A. |
[in] | ldx | INTEGER The leading dimension of the array X. LDX >= M. |
[out] | dX | COMPLEX_16 array, dimension (LDDX,NB) Copy of X on GPU. |
[in] | lddx | INTEGER The leading dimension of the array dX. LDDX >= M. |
[out] | Y | COMPLEX_16 array, dimension (LDY,NB) The n-by-nb matrix Y required to update the unreduced part of A. |
[in] | ldy | INTEGER The leading dimension of the array Y. LDY >= N. |
[out] | dY | COMPLEX_16 array, dimension (LDDY,NB) Copy of Y on GPU. |
[in] | lddy | INTEGER The leading dimension of the array dY. LDDY >= N. |
work | COMPLEX_16 array, dimension (LWORK) Workspace. | |
[in] | lwork | INTEGER The dimension of the array WORK. LWORK >= max( M, N ). |
[in] | queue | magma_queue_t Queue to execute in. |
The matrices Q and P are represented as products of elementary reflectors:
Q = H(1) H(2) . . . H(nb) and P = G(1) G(2) . . . G(nb)
Each H(i) and G(i) has the form:
H(i) = I - tauq * v * v' and G(i) = I - taup * u * u'
where tauq and taup are complex scalars, and v and u are complex vectors.
If m >= n, v(1:i-1) = 0, v(i) = 1, and v(i:m) is stored on exit in A(i:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
If m < n, v(1:i) = 0, v(i+1) = 1, and v(i+1:m) is stored on exit in A(i+2:m,i); u(1:i-1) = 0, u(i) = 1, and u(i:n) is stored on exit in A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
The elements of the vectors v and u together form the m-by-nb matrix V and the nb-by-n matrix U' which are needed, with X and Y, to apply the transformation to the unreduced part of the matrix, using a block update of the form: A := A - V*Y' - X*U'.
The contents of A on exit are illustrated by the following examples with nb = 2:
m = 6 and n = 5 (m > n): m = 5 and n = 6 (m < n): ( 1 1 u1 u1 u1 ) ( 1 u1 u1 u1 u1 u1 ) ( v1 1 1 u2 u2 ) ( 1 1 u2 u2 u2 u2 ) ( v1 v2 a a a ) ( v1 1 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a ) ( v1 v2 a a a a ) ( v1 v2 a a a )
where a denotes an element of the original matrix which is unchanged, vi denotes an element of the vector defining H(i), and ui an element of the vector defining G(i).