MAGMA  2.7.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
lahru: Partial factorization; used by gehrd

Functions

magma_int_t magma_clahru (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, magmaFloatComplex *A, magma_int_t lda, magmaFloatComplex_ptr dA, magma_int_t ldda, magmaFloatComplex_ptr dY, magma_int_t lddy, magmaFloatComplex_ptr dV, magma_int_t lddv, magmaFloatComplex_ptr dT, magmaFloatComplex_ptr dwork, magma_queue_t queue)
 CLAHRU is an auxiliary MAGMA routine that is used in CGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 
magma_int_t magma_clahru_m (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, magmaFloatComplex *A, magma_int_t lda, struct cgehrd_data *data)
 CLAHRU is an auxiliary MAGMA routine that is used in CGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 
magma_int_t magma_dlahru (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, double *A, magma_int_t lda, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dY, magma_int_t lddy, magmaDouble_ptr dV, magma_int_t lddv, magmaDouble_ptr dT, magmaDouble_ptr dwork, magma_queue_t queue)
 DLAHRU is an auxiliary MAGMA routine that is used in DGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 
magma_int_t magma_dlahru_m (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, double *A, magma_int_t lda, struct dgehrd_data *data)
 DLAHRU is an auxiliary MAGMA routine that is used in DGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 
magma_int_t magma_slahru (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, float *A, magma_int_t lda, magmaFloat_ptr dA, magma_int_t ldda, magmaFloat_ptr dY, magma_int_t lddy, magmaFloat_ptr dV, magma_int_t lddv, magmaFloat_ptr dT, magmaFloat_ptr dwork, magma_queue_t queue)
 SLAHRU is an auxiliary MAGMA routine that is used in SGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 
magma_int_t magma_slahru_m (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, float *A, magma_int_t lda, struct sgehrd_data *data)
 SLAHRU is an auxiliary MAGMA routine that is used in SGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 
magma_int_t magma_zlahru (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, magmaDoubleComplex *A, magma_int_t lda, magmaDoubleComplex_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dY, magma_int_t lddy, magmaDoubleComplex_ptr dV, magma_int_t lddv, magmaDoubleComplex_ptr dT, magmaDoubleComplex_ptr dwork, magma_queue_t queue)
 ZLAHRU is an auxiliary MAGMA routine that is used in ZGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 
magma_int_t magma_zlahru_m (magma_int_t n, magma_int_t ihi, magma_int_t k, magma_int_t nb, magmaDoubleComplex *A, magma_int_t lda, struct zgehrd_data *data)
 ZLAHRU is an auxiliary MAGMA routine that is used in ZGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels. More...
 

Detailed Description

Function Documentation

magma_int_t magma_clahru ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
magmaFloatComplex *  A,
magma_int_t  lda,
magmaFloatComplex_ptr  dA,
magma_int_t  ldda,
magmaFloatComplex_ptr  dY,
magma_int_t  lddy,
magmaFloatComplex_ptr  dV,
magma_int_t  lddv,
magmaFloatComplex_ptr  dT,
magmaFloatComplex_ptr  dwork,
magma_queue_t  queue 
)

CLAHRU is an auxiliary MAGMA routine that is used in CGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in cgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]ACOMPLEX array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dACOMPLEX array on the GPU, dimension (LDDA,N-K). On entry, the N-by-(N-K) general matrix to be updated. On exit, the 1st K rows (matrix Am) of A are updated by applying an orthogonal transformation from the right Am = Am (I-V T V'), and sub-matrix Ag is updated by Ag = (I - V T V') Ag (I - V T V(NB+1:)' ) where Q = I - V T V' represent the orthogonal matrix (as a product of elementary reflectors V) used to reduce the current panel of A to upper Hessenberg form. After Am is updated Am(:,1:NB) is sent to the CPU. See Further Details below.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,N).
[in,out]dY(workspace) COMPLEX array on the GPU, dimension (LDDY, NB). On entry the (N-K)-by-NB Y = A V. It is used internally as workspace, so its value is changed on exit.
[in]lddyINTEGER The leading dimension of the array dY. LDDY >= max(1,N).
[in,out]dV(workspace) COMPLEX array on the GPU, dimension (LDDV, NB). On entry the (N-K)-by-NB matrix V of elementary reflectors used to reduce the current panel of A to upper Hessenberg form. The rest K-by-NB part is used as workspace. V is unchanged on exit.
[in]lddvINTEGER The leading dimension of the array dV. LDDV >= max(1,N).
[in]dTCOMPLEX array on the GPU, dimension (NB, NB). On entry the NB-by-NB upper trinagular matrix defining the orthogonal Hessenberg reduction transformation matrix for the current panel. The lower triangular part are 0s.
dwork(workspace) COMPLEX array on the GPU, dimension N*NB.
[in]queuemagma_queue_t Queue to execute in.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.

magma_int_t magma_clahru_m ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
magmaFloatComplex *  A,
magma_int_t  lda,
struct cgehrd_data data 
)

CLAHRU is an auxiliary MAGMA routine that is used in CGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in cgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]ACOMPLEX array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dataStructure with pointers to dA, dT, dV, dW, dY which are distributed across multiple GPUs.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.

magma_int_t magma_dlahru ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
double *  A,
magma_int_t  lda,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dY,
magma_int_t  lddy,
magmaDouble_ptr  dV,
magma_int_t  lddv,
magmaDouble_ptr  dT,
magmaDouble_ptr  dwork,
magma_queue_t  queue 
)

DLAHRU is an auxiliary MAGMA routine that is used in DGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in dgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]ADOUBLE PRECISION array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dADOUBLE PRECISION array on the GPU, dimension (LDDA,N-K). On entry, the N-by-(N-K) general matrix to be updated. On exit, the 1st K rows (matrix Am) of A are updated by applying an orthogonal transformation from the right Am = Am (I-V T V'), and sub-matrix Ag is updated by Ag = (I - V T V') Ag (I - V T V(NB+1:)' ) where Q = I - V T V' represent the orthogonal matrix (as a product of elementary reflectors V) used to reduce the current panel of A to upper Hessenberg form. After Am is updated Am(:,1:NB) is sent to the CPU. See Further Details below.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,N).
[in,out]dY(workspace) DOUBLE PRECISION array on the GPU, dimension (LDDY, NB). On entry the (N-K)-by-NB Y = A V. It is used internally as workspace, so its value is changed on exit.
[in]lddyINTEGER The leading dimension of the array dY. LDDY >= max(1,N).
[in,out]dV(workspace) DOUBLE PRECISION array on the GPU, dimension (LDDV, NB). On entry the (N-K)-by-NB matrix V of elementary reflectors used to reduce the current panel of A to upper Hessenberg form. The rest K-by-NB part is used as workspace. V is unchanged on exit.
[in]lddvINTEGER The leading dimension of the array dV. LDDV >= max(1,N).
[in]dTDOUBLE PRECISION array on the GPU, dimension (NB, NB). On entry the NB-by-NB upper trinagular matrix defining the orthogonal Hessenberg reduction transformation matrix for the current panel. The lower triangular part are 0s.
dwork(workspace) DOUBLE PRECISION array on the GPU, dimension N*NB.
[in]queuemagma_queue_t Queue to execute in.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.

magma_int_t magma_dlahru_m ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
double *  A,
magma_int_t  lda,
struct dgehrd_data data 
)

DLAHRU is an auxiliary MAGMA routine that is used in DGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in dgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]ADOUBLE PRECISION array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dataStructure with pointers to dA, dT, dV, dW, dY which are distributed across multiple GPUs.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.

magma_int_t magma_slahru ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
float *  A,
magma_int_t  lda,
magmaFloat_ptr  dA,
magma_int_t  ldda,
magmaFloat_ptr  dY,
magma_int_t  lddy,
magmaFloat_ptr  dV,
magma_int_t  lddv,
magmaFloat_ptr  dT,
magmaFloat_ptr  dwork,
magma_queue_t  queue 
)

SLAHRU is an auxiliary MAGMA routine that is used in SGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in sgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]AREAL array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dAREAL array on the GPU, dimension (LDDA,N-K). On entry, the N-by-(N-K) general matrix to be updated. On exit, the 1st K rows (matrix Am) of A are updated by applying an orthogonal transformation from the right Am = Am (I-V T V'), and sub-matrix Ag is updated by Ag = (I - V T V') Ag (I - V T V(NB+1:)' ) where Q = I - V T V' represent the orthogonal matrix (as a product of elementary reflectors V) used to reduce the current panel of A to upper Hessenberg form. After Am is updated Am(:,1:NB) is sent to the CPU. See Further Details below.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,N).
[in,out]dY(workspace) REAL array on the GPU, dimension (LDDY, NB). On entry the (N-K)-by-NB Y = A V. It is used internally as workspace, so its value is changed on exit.
[in]lddyINTEGER The leading dimension of the array dY. LDDY >= max(1,N).
[in,out]dV(workspace) REAL array on the GPU, dimension (LDDV, NB). On entry the (N-K)-by-NB matrix V of elementary reflectors used to reduce the current panel of A to upper Hessenberg form. The rest K-by-NB part is used as workspace. V is unchanged on exit.
[in]lddvINTEGER The leading dimension of the array dV. LDDV >= max(1,N).
[in]dTREAL array on the GPU, dimension (NB, NB). On entry the NB-by-NB upper trinagular matrix defining the orthogonal Hessenberg reduction transformation matrix for the current panel. The lower triangular part are 0s.
dwork(workspace) REAL array on the GPU, dimension N*NB.
[in]queuemagma_queue_t Queue to execute in.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.

magma_int_t magma_slahru_m ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
float *  A,
magma_int_t  lda,
struct sgehrd_data data 
)

SLAHRU is an auxiliary MAGMA routine that is used in SGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in sgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]AREAL array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dataStructure with pointers to dA, dT, dV, dW, dY which are distributed across multiple GPUs.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.

magma_int_t magma_zlahru ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
magmaDoubleComplex *  A,
magma_int_t  lda,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dY,
magma_int_t  lddy,
magmaDoubleComplex_ptr  dV,
magma_int_t  lddv,
magmaDoubleComplex_ptr  dT,
magmaDoubleComplex_ptr  dwork,
magma_queue_t  queue 
)

ZLAHRU is an auxiliary MAGMA routine that is used in ZGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in zgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]ACOMPLEX_16 array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dACOMPLEX_16 array on the GPU, dimension (LDDA,N-K). On entry, the N-by-(N-K) general matrix to be updated. On exit, the 1st K rows (matrix Am) of A are updated by applying an orthogonal transformation from the right Am = Am (I-V T V'), and sub-matrix Ag is updated by Ag = (I - V T V') Ag (I - V T V(NB+1:)' ) where Q = I - V T V' represent the orthogonal matrix (as a product of elementary reflectors V) used to reduce the current panel of A to upper Hessenberg form. After Am is updated Am(:,1:NB) is sent to the CPU. See Further Details below.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,N).
[in,out]dY(workspace) COMPLEX_16 array on the GPU, dimension (LDDY, NB). On entry the (N-K)-by-NB Y = A V. It is used internally as workspace, so its value is changed on exit.
[in]lddyINTEGER The leading dimension of the array dY. LDDY >= max(1,N).
[in,out]dV(workspace) COMPLEX_16 array on the GPU, dimension (LDDV, NB). On entry the (N-K)-by-NB matrix V of elementary reflectors used to reduce the current panel of A to upper Hessenberg form. The rest K-by-NB part is used as workspace. V is unchanged on exit.
[in]lddvINTEGER The leading dimension of the array dV. LDDV >= max(1,N).
[in]dTCOMPLEX_16 array on the GPU, dimension (NB, NB). On entry the NB-by-NB upper trinagular matrix defining the orthogonal Hessenberg reduction transformation matrix for the current panel. The lower triangular part are 0s.
dwork(workspace) COMPLEX_16 array on the GPU, dimension N*NB.
[in]queuemagma_queue_t Queue to execute in.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.

magma_int_t magma_zlahru_m ( magma_int_t  n,
magma_int_t  ihi,
magma_int_t  k,
magma_int_t  nb,
magmaDoubleComplex *  A,
magma_int_t  lda,
struct zgehrd_data data 
)

ZLAHRU is an auxiliary MAGMA routine that is used in ZGEHRD to update the trailing sub-matrices after the reductions of the corresponding panels.

See further details below.

Parameters
[in]nINTEGER The order of the matrix A. N >= 0.
[in]ihiINTEGER Last row to update. Same as IHI in zgehrd.
[in]kINTEGER Number of rows of the matrix Am (see details below)
[in]nbINTEGER Block size
[out]ACOMPLEX_16 array, dimension (LDA,N-K) On entry, the N-by-(N-K) general matrix to be updated. The computation is done on the GPU. After Am is updated on the GPU only Am(1:NB) is transferred to the CPU - to update the corresponding Am matrix. See Further Details below.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in,out]dataStructure with pointers to dA, dT, dV, dW, dY which are distributed across multiple GPUs.

Further Details

This implementation follows the algorithm and notations described in:

S. Tomov and J. Dongarra, "Accelerating the reduction to upper Hessenberg form through hybrid GPU-based computing," University of Tennessee Computer Science Technical Report, UT-CS-09-642 (also LAPACK Working Note 219), May 24, 2009.

The difference is that here Am is computed on the GPU. M is renamed Am, G is renamed Ag.