MAGMA  1.6.3
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
double-complex precision

Functions

magma_int_t magma_znan_inf (magma_uplo_t uplo, magma_int_t m, magma_int_t n, const magmaDoubleComplex *A, magma_int_t lda, magma_int_t *cnt_nan, magma_int_t *cnt_inf)
 magma_znan_inf checks a matrix that is located on the CPU host for NAN (not-a-number) and INF (infinity) values. More...
 
magma_int_t magma_znan_inf_gpu (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magma_int_t *cnt_nan, magma_int_t *cnt_inf)
 magma_znan_inf checks a matrix that is located on the CPU host for NAN (not-a-number) and INF (infinity) values. More...
 
void magma_zprint (magma_int_t m, magma_int_t n, const magmaDoubleComplex *A, magma_int_t lda)
 magma_zprint prints a matrix that is located on the CPU host. More...
 
void magma_zprint_gpu (magma_int_t m, magma_int_t n, const magmaDoubleComplex *dA, magma_int_t ldda)
 magma_zprint_gpu prints a matrix that is located on the GPU device. More...
 
void magmablas_clat2z_q (magma_uplo_t uplo, magma_int_t n, magmaFloatComplex_const_ptr SA, magma_int_t ldsa, magmaDoubleComplex_ptr A, magma_int_t lda, magma_int_t *info, magma_queue_t queue)
 CLAT2Z_STREAM converts a single-complex matrix, SA, to a double-complex matrix, A. More...
 
void magmablas_clat2z (magma_uplo_t uplo, magma_int_t n, magmaFloatComplex_const_ptr SA, magma_int_t ldsa, magmaDoubleComplex_ptr A, magma_int_t lda, magma_int_t *info)
 
void magmablas_zclaswp_q (magma_int_t n, magmaDoubleComplex_ptr A, magma_int_t lda, magmaFloatComplex_ptr SA, magma_int_t m, const magma_int_t *ipiv, magma_int_t incx, magma_queue_t queue)
 Row i of A is cast to single precision in row ipiv[i] of SA (incx > 0), or row i of SA is cast to double precision in row ipiv[i] of A (incx < 0), for 0 <= i < M. More...
 
void magmablas_zclaswp (magma_int_t n, magmaDoubleComplex_ptr A, magma_int_t lda, magmaFloatComplex_ptr SA, magma_int_t m, const magma_int_t *ipiv, magma_int_t incx)
 
void magmablas_zgeadd_q (magma_int_t m, magma_int_t n, magmaDoubleComplex alpha, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb, magma_queue_t queue)
 ZGEADD adds two matrices, dB = alpha*dA + dB. More...
 
void magmablas_zgeadd (magma_int_t m, magma_int_t n, magmaDoubleComplex alpha, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb)
 
void magmablas_zgeadd_batched_q (magma_int_t m, magma_int_t n, magmaDoubleComplex alpha, magmaDoubleComplex_const_ptr const dAarray[], magma_int_t ldda, magmaDoubleComplex_ptr dBarray[], magma_int_t lddb, magma_int_t batchCount, magma_queue_t queue)
 ZGEADD adds two sets of matrices, dAarray[i] = alpha*dAarray[i] + dBarray[i], for i = 0, ..., batchCount-1. More...
 
void magmablas_zgeadd_batched (magma_int_t m, magma_int_t n, magmaDoubleComplex alpha, magmaDoubleComplex_const_ptr const dAarray[], magma_int_t ldda, magmaDoubleComplex_ptr dBarray[], magma_int_t lddb, magma_int_t batchCount)
 
void magmablas_zlacpy_q (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb, magma_queue_t queue)
 ZLACPY_Q copies all or part of a two-dimensional matrix dA to another matrix dB. More...
 
void magmablas_zlacpy (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb)
 
void magmablas_zlacpy_batched (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr const dAarray[], magma_int_t ldda, magmaDoubleComplex_ptr dBarray[], magma_int_t lddb, magma_int_t batchCount, magma_queue_t queue)
 ZLACPY_BATCHED copies all or part of each two-dimensional matrix dAarray[i] to matrix dBarray[i], for 0 <= i < batchcount. More...
 
void magmablas_zlacpy_sym_in_q (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magma_int_t *rows, magma_int_t *perm, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb, magma_queue_t queue)
 ZLACPY_Q copies all or part of a two-dimensional matrix dA to another matrix dB. More...
 
void magmablas_zlacpy_sym_in (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magma_int_t *rows, magma_int_t *perm, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb)
 
void magmablas_zlacpy_sym_out_q (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magma_int_t *rows, magma_int_t *perm, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb, magma_queue_t queue)
 ZLACPY_Q copies all or part of a two-dimensional matrix dA to another matrix dB. More...
 
void magmablas_zlacpy_sym_out (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magma_int_t *rows, magma_int_t *perm, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dB, magma_int_t lddb)
 
void magmablas_zlag2c_q (magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr A, magma_int_t lda, magmaFloatComplex_ptr SA, magma_int_t ldsa, magma_int_t *info, magma_queue_t queue)
 ZLAG2C_Q converts a double-complex matrix, A, to a single-complex matrix, SA. More...
 
void magmablas_zlag2c (magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr A, magma_int_t lda, magmaFloatComplex_ptr SA, magma_int_t ldsa, magma_int_t *info)
 
double magmablas_zlange (magma_norm_t norm, magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDouble_ptr dwork)
 ZLANGE returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real matrix A. More...
 
double magmablas_zlanhe (magma_norm_t norm, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDouble_ptr dwork)
 ZLANHE returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix A. More...
 
void magmablas_zlascl_q (magma_type_t type, magma_int_t kl, magma_int_t ku, double cfrom, double cto, magma_int_t m, magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_queue_t queue, magma_int_t *info)
 ZLASCL multiplies the M by N complex matrix A by the real scalar CTO/CFROM. More...
 
void magmablas_zlascl (magma_type_t type, magma_int_t kl, magma_int_t ku, double cfrom, double cto, magma_int_t m, magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *info)
 
void magmablas_zlascl2_q (magma_type_t type, magma_int_t m, magma_int_t n, magmaDouble_const_ptr dD, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_queue_t queue, magma_int_t *info)
 ZLASCL2 scales the M by N complex matrix A by the real diagonal matrix dD. More...
 
void magmablas_zlascl2 (magma_type_t type, magma_int_t m, magma_int_t n, magmaDouble_const_ptr dD, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *info)
 
void magmablas_zlascl_2x2_q (magma_type_t type, magma_int_t m, const magmaDoubleComplex *dW, magma_int_t lddw, magmaDoubleComplex *dA, magma_int_t ldda, magma_int_t *info, magma_queue_t queue)
 ZLASCL_2x2 scales the M by M complex matrix A by the 2-by-2 pivot. More...
 
void magmablas_zlascl_2x2 (magma_type_t type, magma_int_t m, magmaDoubleComplex *dW, magma_int_t lddw, magmaDoubleComplex *dA, magma_int_t ldda, magma_int_t *info)
 
void magmablas_zlascl_diag_q (magma_type_t type, magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dD, magma_int_t lddd, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *info, magma_queue_t queue)
 ZLASCL_DIAG scales the M by N complex matrix A by the real diagonal matrix dD. More...
 
void magmablas_zlascl_diag (magma_type_t type, magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dD, magma_int_t lddd, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t *info)
 
void magmablas_zlaset_q (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magmaDoubleComplex offdiag, magmaDoubleComplex diag, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_queue_t queue)
 ZLASET_Q initializes a 2-D array A to DIAG on the diagonal and OFFDIAG on the off-diagonals. More...
 
void magmablas_zlaset (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magmaDoubleComplex offdiag, magmaDoubleComplex diag, magmaDoubleComplex_ptr dA, magma_int_t ldda)
 
void magmablas_zlaset_band_q (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magma_int_t k, magmaDoubleComplex offdiag, magmaDoubleComplex diag, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_queue_t queue)
 ZLASET_BAND_STREAM initializes the main diagonal of dA to DIAG, and the K-1 sub- or super-diagonals to OFFDIAG. More...
 
void magmablas_zlaset_band (magma_uplo_t uplo, magma_int_t m, magma_int_t n, magma_int_t k, magmaDoubleComplex offdiag, magmaDoubleComplex diag, magmaDoubleComplex_ptr dA, magma_int_t ldda)
 
void magmablas_zlaswp_q (magma_int_t n, magmaDoubleComplex_ptr dAT, magma_int_t ldda, magma_int_t k1, magma_int_t k2, const magma_int_t *ipiv, magma_int_t inci, magma_queue_t queue)
 

Purpose:

ZLASWP performs a series of row interchanges on the matrix A. More...
 
void magmablas_zlaswp (magma_int_t n, magmaDoubleComplex_ptr dAT, magma_int_t ldda, magma_int_t k1, magma_int_t k2, const magma_int_t *ipiv, magma_int_t inci)
 
void magmablas_zlaswpx_q (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldx, magma_int_t ldy, magma_int_t k1, magma_int_t k2, const magma_int_t *ipiv, magma_int_t inci, magma_queue_t queue)
 

Purpose:

ZLASWPX performs a series of row interchanges on the matrix A. More...
 
void magmablas_zlaswpx (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldx, magma_int_t ldy, magma_int_t k1, magma_int_t k2, const magma_int_t *ipiv, magma_int_t inci)
 
void magmablas_zlaswp2_q (magma_int_t n, magmaDoubleComplex_ptr dAT, magma_int_t ldda, magma_int_t k1, magma_int_t k2, magmaInt_const_ptr d_ipiv, magma_int_t inci, magma_queue_t queue)
 

Purpose:

ZLASWP2 performs a series of row interchanges on the matrix A. More...
 
void magmablas_zlaswp2 (magma_int_t n, magmaDoubleComplex_ptr dAT, magma_int_t ldda, magma_int_t k1, magma_int_t k2, magmaInt_const_ptr d_ipiv, magma_int_t inci)
 
void magmablas_zlaswp_sym_q (magma_int_t n, magmaDoubleComplex *dA, magma_int_t lda, magma_int_t k1, magma_int_t k2, const magma_int_t *ipiv, magma_int_t inci, magma_queue_t queue)
 

Purpose:

ZLASWPX performs a series of row interchanges on the matrix A. More...
 
void magmablas_zlaswp_sym (magma_int_t n, magmaDoubleComplex *dA, magma_int_t lda, magma_int_t k1, magma_int_t k2, const magma_int_t *ipiv, magma_int_t inci)
 
void magmablas_zlat2c_q (magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex_const_ptr A, magma_int_t lda, magmaFloatComplex_ptr SA, magma_int_t ldsa, magma_int_t *info, magma_queue_t queue)
 ZLAT2C converts a double-complex matrix, A, to a single-complex matrix, SA. More...
 
void magmablas_zlat2c (magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex_const_ptr A, magma_int_t lda, magmaFloatComplex_ptr SA, magma_int_t ldsa, magma_int_t *info)
 
void magmablas_zswapdblk_q (magma_int_t n, magma_int_t nb, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t inca, magmaDoubleComplex_ptr dB, magma_int_t lddb, magma_int_t incb, magma_queue_t queue)
 zswapdblk swaps diagonal blocks of size nb x nb between matrices dA and dB on the GPU. More...
 
void magmablas_zswapdblk (magma_int_t n, magma_int_t nb, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t inca, magmaDoubleComplex_ptr dB, magma_int_t lddb, magma_int_t incb)
 
void magmablas_zswapdblk_batched_q (magma_int_t n, magma_int_t nb, magmaDoubleComplex **dA_array, magma_int_t ldda, magma_int_t inca, magmaDoubleComplex **dB_array, magma_int_t lddb, magma_int_t incb, magma_int_t batchCount, magma_queue_t queue)
 zswapdblk swaps diagonal blocks of size nb x nb between matrices dA and dB on the GPU. More...
 
void magmablas_zswapdblk_batched (magma_int_t n, magma_int_t nb, magmaDoubleComplex **dA_array, magma_int_t ldda, magma_int_t inca, magmaDoubleComplex **dB_array, magma_int_t lddb, magma_int_t incb, magma_int_t batchCount)
 
void magmablas_zsymmetrize_q (magma_uplo_t uplo, magma_int_t m, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_queue_t queue)
 ZSYMMETRIZE copies lower triangle to upper triangle, or vice-versa, to make dA a general representation of a symmetric matrix. More...
 
void magmablas_zsymmetrize (magma_uplo_t uplo, magma_int_t m, magmaDoubleComplex_ptr dA, magma_int_t ldda)
 
void magmablas_zsymmetrize_tiles_q (magma_uplo_t uplo, magma_int_t m, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t ntile, magma_int_t mstride, magma_int_t nstride, magma_queue_t queue)
 ZSYMMETRIZE_TILES copies lower triangle to upper triangle, or vice-versa, to make some blocks of dA into general representations of a symmetric block. More...
 
void magmablas_zsymmetrize_tiles (magma_uplo_t uplo, magma_int_t m, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_int_t ntile, magma_int_t mstride, magma_int_t nstride)
 
void magmablas_ztranspose_q (magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dAT, magma_int_t lddat, magma_queue_t queue)
 ztranspose_q copies and transposes a matrix dA to matrix dAT. More...
 
void magmablas_ztranspose (magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dAT, magma_int_t lddat)
 
void magmablas_ztranspose_batched_q (magma_int_t m, magma_int_t n, magmaDoubleComplex **dA_array, magma_int_t ldda, magmaDoubleComplex **dAT_array, magma_int_t lddat, magma_int_t batchCount, magma_queue_t queue)
 ztranspose_batched_q copies and transposes a matrix dA_array[i] to matrix dAT_array[i]. More...
 
void magmablas_ztranspose_batched (magma_int_t m, magma_int_t n, magmaDoubleComplex **dA_array, magma_int_t ldda, magmaDoubleComplex **dAT_array, magma_int_t lddat, magma_int_t batchCount)
 
void magmablas_ztranspose_conj_q (magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dAT, magma_int_t lddat, magma_queue_t queue)
 ztranspose_conj_q copies and conjugate-transposes a matrix dA to matrix dAT. More...
 
void magmablas_ztranspose_conj (magma_int_t m, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_ptr dAT, magma_int_t lddat)
 
void magmablas_ztranspose_conj_batched_q (magma_int_t m, magma_int_t n, magmaDoubleComplex **dA_array, magma_int_t ldda, magmaDoubleComplex **dAT_array, magma_int_t lddat, magma_int_t batchCount, magma_queue_t queue)
 ztranspose_conj_batched_q copies and conjugate-transposes a matrix dA_array[i] to matrix dAT_array[i]. More...
 
void magmablas_ztranspose_conj_batched (magma_int_t m, magma_int_t n, magmaDoubleComplex **dA_array, magma_int_t ldda, magmaDoubleComplex **dAT_array, magma_int_t lddat, magma_int_t batchCount)
 
void magmablas_ztranspose_conj_inplace_q (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_queue_t queue)
 ztranspose_conj_inplace_q conjugate-transposes a square N-by-N matrix in-place. More...
 
void magmablas_ztranspose_conj_inplace (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda)
 
void magmablas_ztranspose_inplace_q (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, magma_queue_t queue)
 ztranspose_inplace_q transposes a square N-by-N matrix in-place. More...
 
void magmablas_ztranspose_inplace (magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda)
 

Detailed Description

Function Documentation

magma_int_t magma_znan_inf ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
const magmaDoubleComplex *  A,
magma_int_t  lda,
magma_int_t *  cnt_nan,
magma_int_t *  cnt_inf 
)

magma_znan_inf checks a matrix that is located on the CPU host for NAN (not-a-number) and INF (infinity) values.

NAN is created by 0/0 and similar. INF is created by x/0 and similar, where x != 0.

Parameters
[in]uplomagma_uplo_t Specifies what part of the matrix A to check.
  • = MagmaUpper: Upper triangular part of A
  • = MagmaLower: Lower triangular part of A
  • = MagmaFull: All of A
[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]ACOMPLEX_16 array, dimension (LDA,N), on the CPU host. The M-by-N matrix to be printed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
[out]cnt_nanINTEGER* If non-NULL, on exit contains the number of NAN values in A.
[out]cnt_infINTEGER* If non-NULL, on exit contains the number of INF values in A.
Returns
  • >= 0: Returns number of NAN + number of INF values.
  • < 0: If it returns -i, the i-th argument had an illegal value, or another error occured, such as memory allocation failed.
magma_int_t magma_znan_inf_gpu ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magma_int_t *  cnt_nan,
magma_int_t *  cnt_inf 
)

magma_znan_inf checks a matrix that is located on the CPU host for NAN (not-a-number) and INF (infinity) values.

NAN is created by 0/0 and similar. INF is created by x/0 and similar, where x != 0.

Parameters
[in]uplomagma_uplo_t Specifies what part of the matrix A to check.
  • = MagmaUpper: Upper triangular part of A
  • = MagmaLower: Lower triangular part of A
  • = MagmaFull: All of A
[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]dACOMPLEX_16 array, dimension (LDDA,N), on the GPU device. The M-by-N matrix to be printed.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]cnt_nanINTEGER* If non-NULL, on exit contains the number of NAN values in A.
[out]cnt_infINTEGER* If non-NULL, on exit contains the number of INF values in A.
Returns
  • >= 0: Returns number of NAN + number of INF values.
  • < 0: If it returns -i, the i-th argument had an illegal value, or another error occured, such as memory allocation failed.
void magma_zprint ( magma_int_t  m,
magma_int_t  n,
const magmaDoubleComplex *  A,
magma_int_t  lda 
)

magma_zprint prints a matrix that is located on the CPU host.

The output is intended to be Matlab compatible, to be useful in debugging.

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]ACOMPLEX_16 array, dimension (LDA,N), on the CPU host. The M-by-N matrix to be printed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,M).
void magma_zprint_gpu ( magma_int_t  m,
magma_int_t  n,
const magmaDoubleComplex *  dA,
magma_int_t  ldda 
)

magma_zprint_gpu prints a matrix that is located on the GPU device.

Internally, it allocates CPU memory and copies the matrix to the CPU. The output is intended to be Matlab compatible, to be useful in debugging.

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]dACOMPLEX_16 array, dimension (LDDA,N), on the GPU device. The M-by-N matrix to be printed.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
void magmablas_clat2z ( magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex_const_ptr  SA,
magma_int_t  ldsa,
magmaDoubleComplex_ptr  A,
magma_int_t  lda,
magma_int_t *  info 
)
void magmablas_clat2z_q ( magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex_const_ptr  SA,
magma_int_t  ldsa,
magmaDoubleComplex_ptr  A,
magma_int_t  lda,
magma_int_t *  info,
magma_queue_t  queue 
)

CLAT2Z_STREAM converts a single-complex matrix, SA, to a double-complex matrix, A.

Note that while it is possible to overflow while converting from double to single, it is not possible to overflow when converting from single to double.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix A to be converted.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
[in]nINTEGER The number of columns of the matrix A. n >= 0.
[in]ACOMPLEX_16 array, dimension (LDA,n) On entry, the n-by-n coefficient matrix A.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,n).
[out]SACOMPLEX array, dimension (LDSA,n) On exit, if INFO=0, the n-by-n coefficient matrix SA; if INFO > 0, the content of SA is unspecified.
[in]ldsaINTEGER The leading dimension of the array SA. LDSA >= max(1,n).
[out]infoINTEGER
  • = 0: successful exit.
  • < 0: if INFO = -i, the i-th argument had an illegal value
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zclaswp ( magma_int_t  n,
magmaDoubleComplex_ptr  A,
magma_int_t  lda,
magmaFloatComplex_ptr  SA,
magma_int_t  m,
const magma_int_t *  ipiv,
magma_int_t  incx 
)
void magmablas_zclaswp_q ( magma_int_t  n,
magmaDoubleComplex_ptr  A,
magma_int_t  lda,
magmaFloatComplex_ptr  SA,
magma_int_t  m,
const magma_int_t *  ipiv,
magma_int_t  incx,
magma_queue_t  queue 
)

Row i of A is cast to single precision in row ipiv[i] of SA (incx > 0), or row i of SA is cast to double precision in row ipiv[i] of A (incx < 0), for 0 <= i < M.

Parameters
[in]nINTEGER. On entry, N specifies the number of columns of the matrix A.
[in,out]ADOUBLE PRECISION array on the GPU, dimension (LDA,N) On entry, the M-by-N matrix to which the row interchanges will be applied. TODO update docs
[in]ldaINTEGER. LDA specifies the leading dimension of A.
[in,out]SAREAL array on the GPU, dimension (LDA,N) On exit, the single precision, permuted matrix. TODO update docs
[in]mThe number of rows to be interchanged.
[in]ipivINTEGER array on the GPU, dimension (M) The vector of pivot indices. Row i of A is cast to single precision in row ipiv[i] of SA, for 0 <= i < m.
[in]incxINTEGER If INCX is negative, the pivots are applied in reverse order, otherwise in straight-forward order.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zgeadd ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  alpha,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb 
)
void magmablas_zgeadd_batched ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  alpha,
magmaDoubleComplex_const_ptr const  dAarray[],
magma_int_t  ldda,
magmaDoubleComplex_ptr  dBarray[],
magma_int_t  lddb,
magma_int_t  batchCount 
)
void magmablas_zgeadd_batched_q ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  alpha,
magmaDoubleComplex_const_ptr const  dAarray[],
magma_int_t  ldda,
magmaDoubleComplex_ptr  dBarray[],
magma_int_t  lddb,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ZGEADD adds two sets of matrices, dAarray[i] = alpha*dAarray[i] + dBarray[i], for i = 0, ..., batchCount-1.

Parameters
[in]mINTEGER The number of rows of each matrix dAarray[i]. M >= 0.
[in]nINTEGER The number of columns of each matrix dAarray[i]. N >= 0.
[in]alphaCOMPLEX_16 The scalar alpha.
[in]dAarrayarray on GPU, dimension(batchCount), of pointers to arrays, with each array a COMPLEX_16 array, dimension (LDDA,N) The m by n matrices dAarray[i].
[in]lddaINTEGER The leading dimension of each array dAarray[i]. LDDA >= max(1,M).
[in,out]dBarrayarray on GPU, dimension(batchCount), of pointers to arrays, with each array a COMPLEX_16 array, dimension (LDDB,N) The m by n matrices dBarray[i].
[in]lddbINTEGER The leading dimension of each array dBarray[i]. LDDB >= max(1,M).
[in]batchCountINTEGER The number of matrices to add; length of dAarray and dBarray. batchCount >= 0.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zgeadd_q ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  alpha,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magma_queue_t  queue 
)

ZGEADD adds two matrices, dB = alpha*dA + dB.

Parameters
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]alphaCOMPLEX_16 The scalar alpha.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The m by n matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[in,out]dBCOMPLEX_16 array, dimension (LDDB,N) The m by n matrix dB.
[in]lddbINTEGER The leading dimension of the array dB. LDDB >= max(1,M).
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlacpy ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb 
)
void magmablas_zlacpy_batched ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr const  dAarray[],
magma_int_t  ldda,
magmaDoubleComplex_ptr  dBarray[],
magma_int_t  lddb,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ZLACPY_BATCHED copies all or part of each two-dimensional matrix dAarray[i] to matrix dBarray[i], for 0 <= i < batchcount.

Parameters
[in]uplomagma_uplo_t Specifies the part of each matrix dA to be copied to dB.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part Otherwise: All of each matrix dA
[in]mINTEGER The number of rows of each matrix dA. M >= 0.
[in]nINTEGER The number of columns of each matrix dA. N >= 0.
[in]dAarrayCOMPLEX_16* array, dimension (batchCount) Array of pointers to the matrices dA, where each dA is of dimension (LDDA,N). The M-by-N matrix dA. If UPLO = MagmaUpper, only the upper triangle or trapezoid is accessed; if UPLO = MagmaLower, only the lower triangle or trapezoid is accessed.
[in]lddaINTEGER The leading dimension of each array dA. LDDA >= max(1,M).
[out]dBarrayCOMPLEX_16* array, dimension (batchCount) Array of pointers to the matrices dB, where each dB is of dimension (LDDB,N). The M-by-N matrix dB. On exit, dB = dA in the locations specified by UPLO.
[in]lddbINTEGER The leading dimension of each array dB. LDDB >= max(1,M).
[in]batchCountNumber of matrices in dAarray and dBarray.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlacpy_q ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magma_queue_t  queue 
)

ZLACPY_Q copies all or part of a two-dimensional matrix dA to another matrix dB.

This is the same as ZLACPY, but adds queue argument.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix dA to be copied to dB.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
  • = MagmaFull: All of the matrix dA
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The M-by-N matrix dA. If UPLO = MagmaUpper, only the upper triangle or trapezoid is accessed; if UPLO = MagmaLower, only the lower triangle or trapezoid is accessed.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[out]dBCOMPLEX_16 array, dimension (LDDB,N) The M-by-N matrix dB. On exit, dB = dA in the locations specified by UPLO.
[in]lddbINTEGER The leading dimension of the array dB. LDDB >= max(1,M).
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlacpy_sym_in ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magma_int_t *  rows,
magma_int_t *  perm,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb 
)
void magmablas_zlacpy_sym_in_q ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magma_int_t *  rows,
magma_int_t *  perm,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magma_queue_t  queue 
)

ZLACPY_Q copies all or part of a two-dimensional matrix dA to another matrix dB.

This is the same as ZLACPY, but adds queue argument.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix dA to be copied to dB.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
  • = MagmaFull: All of the matrix dA
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The M-by-N matrix dA. If UPLO = MagmaUpper, only the upper triangle or trapezoid is accessed; if UPLO = MagmaLower, only the lower triangle or trapezoid is accessed.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[out]dBCOMPLEX_16 array, dimension (LDDB,N) The M-by-N matrix dB. On exit, dB = dA in the locations specified by UPLO.
[in]lddbINTEGER The leading dimension of the array dB. LDDB >= max(1,M).
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlacpy_sym_out ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magma_int_t *  rows,
magma_int_t *  perm,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb 
)
void magmablas_zlacpy_sym_out_q ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magma_int_t *  rows,
magma_int_t *  perm,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magma_queue_t  queue 
)

ZLACPY_Q copies all or part of a two-dimensional matrix dA to another matrix dB.

This is the same as ZLACPY, but adds queue argument.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix dA to be copied to dB.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
  • = MagmaFull: All of the matrix dA
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The M-by-N matrix dA. If UPLO = MagmaUpper, only the upper triangle or trapezoid is accessed; if UPLO = MagmaLower, only the lower triangle or trapezoid is accessed.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[out]dBCOMPLEX_16 array, dimension (LDDB,N) The M-by-N matrix dB. On exit, dB = dA in the locations specified by UPLO.
[in]lddbINTEGER The leading dimension of the array dB. LDDB >= max(1,M).
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlag2c ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  A,
magma_int_t  lda,
magmaFloatComplex_ptr  SA,
magma_int_t  ldsa,
magma_int_t *  info 
)
void magmablas_zlag2c_q ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  A,
magma_int_t  lda,
magmaFloatComplex_ptr  SA,
magma_int_t  ldsa,
magma_int_t *  info,
magma_queue_t  queue 
)

ZLAG2C_Q converts a double-complex matrix, A, to a single-complex matrix, SA.

RMAX is the overflow for the single-complex arithmetic. ZLAG2C checks that all the entries of A are between -RMAX and RMAX. If not, the conversion is aborted and a flag is raised.

This is the same as ZLAG2C, but adds queue argument.

Parameters
[in]mINTEGER The number of lines of the matrix A. m >= 0.
[in]nINTEGER The number of columns of the matrix A. n >= 0.
[in]ACOMPLEX_16 array, dimension (LDA,n) On entry, the m-by-n coefficient matrix A.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,m).
[out]SACOMPLEX array, dimension (LDSA,n) On exit, if INFO=0, the m-by-n coefficient matrix SA; if INFO > 0, the content of SA is unspecified.
[in]ldsaINTEGER The leading dimension of the array SA. LDSA >= max(1,m).
[out]infoINTEGER
  • = 0: successful exit.
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • = 1: an entry of the matrix A is greater than the COMPLEX overflow threshold, in this case, the content of SA on exit is unspecified.
[in]queuemagma_queue_t Queue to execute in.
double magmablas_zlange ( magma_norm_t  norm,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dwork 
)

ZLANGE returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real matrix A.

Description

ZLANGE returns the value

ZLANGE = ( max(abs(A(i,j))), NORM = 'M' or 'm' ( ( norm1(A), NORM = '1', 'O' or 'o' ( ( normI(A), NORM = 'I' or 'i' ( ( normF(A), NORM = 'F', 'f', 'E' or 'e' ** not yet supported

where norm1 denotes the one norm of a matrix (maximum column sum), normI denotes the infinity norm of a matrix (maximum row sum) and normF denotes the Frobenius norm of a matrix (square root of sum of squares). Note that max(abs(A(i,j))) is not a consistent matrix norm.

Parameters
[in]normCHARACTER*1 Specifies the value to be returned in ZLANGE as described above.
[in]mINTEGER The number of rows of the matrix A. M >= 0. When M = 0, ZLANGE is set to zero.
[in]nINTEGER The number of columns of the matrix A. N >= 0. When N = 0, ZLANGE is set to zero.
[in]dADOUBLE PRECISION array on the GPU, dimension (LDDA,N) The m by n matrix A.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(M,1).
dwork(workspace) DOUBLE PRECISION array on the GPU, dimension (LWORK).
double magmablas_zlanhe ( magma_norm_t  norm,
magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dwork 
)

ZLANHE returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix A.

ZLANHE = ( max(abs(A(i,j))), NORM = 'M' or 'm' ( ( norm1(A), NORM = '1', 'O' or 'o' ** supported only for (PRECISION_s || PRECISION_d || PRECISION_c || CUDA_ARCH >= 200) ( ( normI(A), NORM = 'I' or 'i' ** supported only for (PRECISION_s || PRECISION_d || PRECISION_c || CUDA_ARCH >= 200) ( ( normF(A), NORM = 'F', 'f', 'E' or 'e' ** not yet supported

where norm1 denotes the one norm of a matrix (maximum column sum), normI denotes the infinity norm of a matrix (maximum row sum) and normF denotes the Frobenius norm of a matrix (square root of sum of squares). Note that max(abs(A(i,j))) is not a consistent matrix norm.

On error, returns ZLANHE < 0: if ZLANHE = -i, the i-th argument had an illegal value.

Arguments:

Parameters
[in]normCHARACTER*1 Specifies the value to be returned in ZLANHE as described above.
[in]uplomagma_uplo_t Specifies whether the upper or lower triangular part of the Hermitian matrix A is to be referenced.
  • = MagmaUpper: Upper triangular part of A is referenced
  • = MagmaLower: Lower triangular part of A is referenced
[in]nINTEGER The order of the matrix A. N >= 0. When N = 0, ZLANHE is set to zero.
[in]dACOMPLEX*16 array on the GPU, dimension (LDDA,N) The Hermitian matrix A. If UPLO = MagmaUpper, the leading n by n upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading n by n lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. Note that the imaginary parts of the diagonal elements need not be set and are assumed to be zero.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(N,1).
dwork(workspace) DOUBLE PRECISION array on the GPU, dimension (MAX(1,LWORK)), where LWORK >= N. NOTE: this is different than LAPACK, where WORK is required only for norm1 and normI. Here max-norm also requires work.
void magmablas_zlascl ( magma_type_t  type,
magma_int_t  kl,
magma_int_t  ku,
double  cfrom,
double  cto,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t *  info 
)
void magmablas_zlascl2 ( magma_type_t  type,
magma_int_t  m,
magma_int_t  n,
magmaDouble_const_ptr  dD,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t *  info 
)
void magmablas_zlascl2_q ( magma_type_t  type,
magma_int_t  m,
magma_int_t  n,
magmaDouble_const_ptr  dD,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_queue_t  queue,
magma_int_t *  info 
)

ZLASCL2 scales the M by N complex matrix A by the real diagonal matrix dD.

TYPE specifies that A may be full, upper triangular, lower triangular.

Parameters
[in]typemagma_type_t TYPE indices the storage type of the input matrix A. = MagmaFull: full matrix. = MagmaLower: lower triangular matrix. = MagmaUpper: upper triangular matrix. Other formats that LAPACK supports, MAGMA does not currently support.
[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]dDDOUBLE PRECISION vector, dimension (M) The diagonal matrix containing the scalar factors. Stored as a vector.
[in,out]dACOMPLEX*16 array, dimension (LDDA,N) The matrix to be scaled by dD. See TYPE for the storage type.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlascl_2x2 ( magma_type_t  type,
magma_int_t  m,
magmaDoubleComplex *  dW,
magma_int_t  lddw,
magmaDoubleComplex *  dA,
magma_int_t  ldda,
magma_int_t *  info 
)
void magmablas_zlascl_2x2_q ( magma_type_t  type,
magma_int_t  m,
const magmaDoubleComplex *  dW,
magma_int_t  lddw,
magmaDoubleComplex *  dA,
magma_int_t  ldda,
magma_int_t *  info,
magma_queue_t  queue 
)

ZLASCL_2x2 scales the M by M complex matrix A by the 2-by-2 pivot.

TYPE specifies that A may be upper or lower triangular.

Parameters
[in]typemagma_type_t TYPE indices the storage type of the input matrix A. = MagmaLower: lower triangular matrix. = MagmaUpper: upper triangular matrix. Other formats that LAPACK supports, MAGMA does not currently support.
[in]mINTEGER The number of rows of the matrix A. M >= 0.
[in]dWDOUBLE PRECISION vector, dimension (2*lddw) The matrix containing the 2-by-2 pivot.
[in]lddwINTEGER The leading dimension of the array W. LDDA >= max(1,M).
[in,out]dACOMPLEX*16 array, dimension (LDDA,N) The matrix to be scaled by dW. See TYPE for the storage type.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value.
void magmablas_zlascl_diag ( magma_type_t  type,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dD,
magma_int_t  lddd,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t *  info 
)
void magmablas_zlascl_diag_q ( magma_type_t  type,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dD,
magma_int_t  lddd,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t *  info,
magma_queue_t  queue 
)

ZLASCL_DIAG scales the M by N complex matrix A by the real diagonal matrix dD.

TYPE specifies that A may be full, upper triangular, lower triangular.

Parameters
[in]typemagma_type_t TYPE indices the storage type of the input matrix A. = MagmaFull: full matrix. = MagmaLower: lower triangular matrix. = MagmaUpper: upper triangular matrix. Other formats that LAPACK supports, MAGMA does not currently support.
[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]dDDOUBLE PRECISION vector, dimension (LDDD,M) The matrix storing the scaling factor on its diagonal.
[in]ldddINTEGER The leading dimension of the array D.
[in,out]dACOMPLEX*16 array, dimension (LDDA,N) The matrix to be scaled by dD. See TYPE for the storage type.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value.
void magmablas_zlascl_q ( magma_type_t  type,
magma_int_t  kl,
magma_int_t  ku,
double  cfrom,
double  cto,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_queue_t  queue,
magma_int_t *  info 
)

ZLASCL multiplies the M by N complex matrix A by the real scalar CTO/CFROM.

This is done without over/underflow as long as the final result CTO*A(I,J)/CFROM does not over/underflow. TYPE specifies that A may be full, upper triangular, lower triangular.

Parameters
[in]typemagma_type_t TYPE indices the storage type of the input matrix A. = MagmaFull: full matrix. = MagmaLower: lower triangular matrix. = MagmaUpper: upper triangular matrix. Other formats that LAPACK supports, MAGMA does not currently support.
[in]klINTEGER Unused, for LAPACK compatability.
[in]kuKU is INTEGER Unused, for LAPACK compatability.
[in]cfromDOUBLE PRECISION
[in]ctoDOUBLE PRECISION
The matrix A is multiplied by CTO/CFROM. A(I,J) is computed without over/underflow if the final result CTO*A(I,J)/CFROM can be represented without over/underflow. CFROM must be nonzero. CFROM and CTO must not be NAN.
[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]dACOMPLEX*16 array, dimension (LDDA,N) The matrix to be multiplied by CTO/CFROM. See TYPE for the storage type.
[in]lddaINTEGER The leading dimension of the array A. LDDA >= max(1,M).
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlaset ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  offdiag,
magmaDoubleComplex  diag,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda 
)
void magmablas_zlaset_band ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magma_int_t  k,
magmaDoubleComplex  offdiag,
magmaDoubleComplex  diag,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda 
)
void magmablas_zlaset_band_q ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magma_int_t  k,
magmaDoubleComplex  offdiag,
magmaDoubleComplex  diag,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_queue_t  queue 
)

ZLASET_BAND_STREAM initializes the main diagonal of dA to DIAG, and the K-1 sub- or super-diagonals to OFFDIAG.

This is the same as ZLASET_BAND, but adds queue argument.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix dA to be set.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]kINTEGER The number of diagonals to set, including the main diagonal. K >= 0. Currently, K <= 1024 due to CUDA restrictions (max. number of threads per block).
[in]offdiagCOMPLEX_16 Off-diagonal elements in the band are set to OFFDIAG.
[in]diagCOMPLEX_16 All the main diagonal elements are set to DIAG.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The M-by-N matrix dA. If UPLO = MagmaUpper, only the upper triangle or trapezoid is accessed; if UPLO = MagmaLower, only the lower triangle or trapezoid is accessed. On exit, A(i,j) = ALPHA, 1 <= i <= m, 1 <= j <= n where i != j, abs(i-j) < k; and A(i,i) = BETA, 1 <= i <= min(m,n)
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[in]queuemagma_queue_t Stream to execute ZLASET in.
void magmablas_zlaset_q ( magma_uplo_t  uplo,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  offdiag,
magmaDoubleComplex  diag,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_queue_t  queue 
)

ZLASET_Q initializes a 2-D array A to DIAG on the diagonal and OFFDIAG on the off-diagonals.

This is the same as ZLASET, but adds queue argument.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix dA to be set.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
  • = MagmaFull: All of the matrix dA
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]offdiagCOMPLEX_16 The scalar OFFDIAG. (In LAPACK this is called ALPHA.)
[in]diagCOMPLEX_16 The scalar DIAG. (In LAPACK this is called BETA.)
[in]dACOMPLEX_16 array, dimension (LDDA,N) The M-by-N matrix dA. If UPLO = MagmaUpper, only the upper triangle or trapezoid is accessed; if UPLO = MagmaLower, only the lower triangle or trapezoid is accessed. On exit, A(i,j) = OFFDIAG, 1 <= i <= m, 1 <= j <= n, i != j; and A(i,i) = DIAG, 1 <= i <= min(m,n)
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlaswp ( magma_int_t  n,
magmaDoubleComplex_ptr  dAT,
magma_int_t  ldda,
magma_int_t  k1,
magma_int_t  k2,
const magma_int_t *  ipiv,
magma_int_t  inci 
)
void magmablas_zlaswp2 ( magma_int_t  n,
magmaDoubleComplex_ptr  dAT,
magma_int_t  ldda,
magma_int_t  k1,
magma_int_t  k2,
magmaInt_const_ptr  d_ipiv,
magma_int_t  inci 
)
void magmablas_zlaswp2_q ( magma_int_t  n,
magmaDoubleComplex_ptr  dAT,
magma_int_t  ldda,
magma_int_t  k1,
magma_int_t  k2,
magmaInt_const_ptr  d_ipiv,
magma_int_t  inci,
magma_queue_t  queue 
)

Purpose:

ZLASWP2 performs a series of row interchanges on the matrix A.

One row interchange is initiated for each of rows K1 through K2 of A.

Unlike LAPACK, here A is stored row-wise (hence dAT). ** Otherwise, this is identical to LAPACK's interface.

Here, d_ipiv is passed in GPU memory.

Arguments:

Parameters
[in]nINTEGER The number of columns of the matrix A.
[in,out]dATCOMPLEX*16 array on GPU, stored row-wise, dimension (LDDA,*) On entry, the matrix of column dimension N to which the row interchanges will be applied. On exit, the permuted matrix.
[in]lddaINTEGER The leading dimension of the array A. (I.e., stride between elements in a column.)
[in]k1INTEGER The first element of IPIV for which a row interchange will be done. (One based index.)
[in]k2INTEGER The last element of IPIV for which a row interchange will be done. (One based index.)
[in]d_ipivINTEGER array, on GPU, dimension (K2*abs(INCI)) The vector of pivot indices. Only the elements in positions K1 through K2 of IPIV are accessed. IPIV(K) = L implies rows K and L are to be interchanged.
[in]inciINTEGER The increment between successive values of IPIV. Currently, IPIV > 0. TODO: If IPIV is negative, the pivots are applied in reverse order.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlaswp_q ( magma_int_t  n,
magmaDoubleComplex_ptr  dAT,
magma_int_t  ldda,
magma_int_t  k1,
magma_int_t  k2,
const magma_int_t *  ipiv,
magma_int_t  inci,
magma_queue_t  queue 
)

Purpose:

ZLASWP performs a series of row interchanges on the matrix A.

One row interchange is initiated for each of rows K1 through K2 of A.

Unlike LAPACK, here A is stored row-wise (hence dAT). ** Otherwise, this is identical to LAPACK's interface.

Arguments:

Parameters
[in]nINTEGER The number of columns of the matrix A.
[in,out]dATCOMPLEX*16 array on GPU, stored row-wise, dimension (LDDA,N) On entry, the matrix of column dimension N to which the row interchanges will be applied. On exit, the permuted matrix.
[in]lddaINTEGER The leading dimension of the array A. ldda >= n.
[in]k1INTEGER The first element of IPIV for which a row interchange will be done. (Fortran one-based index: 1 <= k1 .)
[in]k2INTEGER The last element of IPIV for which a row interchange will be done. (Fortran one-based index: 1 <= k2 .)
[in]ipivINTEGER array, on CPU, dimension (K2*abs(INCI)) The vector of pivot indices. Only the elements in positions K1 through K2 of IPIV are accessed. IPIV(K) = L implies rows K and L are to be interchanged.
[in]inciINTEGER The increment between successive values of IPIV. Currently, INCI > 0. TODO: If INCI is negative, the pivots are applied in reverse order.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlaswp_sym ( magma_int_t  n,
magmaDoubleComplex *  dA,
magma_int_t  lda,
magma_int_t  k1,
magma_int_t  k2,
const magma_int_t *  ipiv,
magma_int_t  inci 
)
void magmablas_zlaswp_sym_q ( magma_int_t  n,
magmaDoubleComplex *  dA,
magma_int_t  lda,
magma_int_t  k1,
magma_int_t  k2,
const magma_int_t *  ipiv,
magma_int_t  inci,
magma_queue_t  queue 
)

Purpose:

ZLASWPX performs a series of row interchanges on the matrix A.

One row interchange is initiated for each of rows K1 through K2 of A.

Unlike LAPACK, here A is stored either row-wise or column-wise, depending on ldx and ldy. ** Otherwise, this is identical to LAPACK's interface.

Arguments:

Parameters
[in]nINTEGER The number of columns of the matrix A.
[in,out]dACOMPLEX*16 array on GPU, dimension (*,*) On entry, the matrix of column dimension N to which the row interchanges will be applied. On exit, the permuted matrix.
[in]ldxINTEGER Stride between elements in same column.
[in]ldyINTEGER Stride between elements in same row. For A stored row-wise, set ldx=lda and ldy=1. For A stored column-wise, set ldx=1 and ldy=lda.
[in]k1INTEGER The first element of IPIV for which a row interchange will be done. (One based index.)
[in]k2INTEGER The last element of IPIV for which a row interchange will be done. (One based index.)
[in]ipivINTEGER array, on CPU, dimension (K2*abs(INCI)) The vector of pivot indices. Only the elements in positions K1 through K2 of IPIV are accessed. IPIV(K) = L implies rows K and L are to be interchanged.
[in]inciINTEGER The increment between successive values of IPIV. Currently, IPIV > 0. TODO: If IPIV is negative, the pivots are applied in reverse order.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlaswpx ( magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldx,
magma_int_t  ldy,
magma_int_t  k1,
magma_int_t  k2,
const magma_int_t *  ipiv,
magma_int_t  inci 
)
void magmablas_zlaswpx_q ( magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldx,
magma_int_t  ldy,
magma_int_t  k1,
magma_int_t  k2,
const magma_int_t *  ipiv,
magma_int_t  inci,
magma_queue_t  queue 
)

Purpose:

ZLASWPX performs a series of row interchanges on the matrix A.

One row interchange is initiated for each of rows K1 through K2 of A.

Unlike LAPACK, here A is stored either row-wise or column-wise, depending on ldx and ldy. ** Otherwise, this is identical to LAPACK's interface.

Arguments:

Parameters
[in]nINTEGER The number of columns of the matrix A.
[in,out]dACOMPLEX*16 array on GPU, dimension (*,*) On entry, the matrix of column dimension N to which the row interchanges will be applied. On exit, the permuted matrix.
[in]ldxINTEGER Stride between elements in same column.
[in]ldyINTEGER Stride between elements in same row. For A stored row-wise, set ldx=ldda and ldy=1. For A stored column-wise, set ldx=1 and ldy=ldda.
[in]k1INTEGER The first element of IPIV for which a row interchange will be done. (One based index.)
[in]k2INTEGER The last element of IPIV for which a row interchange will be done. (One based index.)
[in]ipivINTEGER array, on CPU, dimension (K2*abs(INCI)) The vector of pivot indices. Only the elements in positions K1 through K2 of IPIV are accessed. IPIV(K) = L implies rows K and L are to be interchanged.
[in]inciINTEGER The increment between successive values of IPIV. Currently, IPIV > 0. TODO: If IPIV is negative, the pivots are applied in reverse order.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zlat2c ( magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex_const_ptr  A,
magma_int_t  lda,
magmaFloatComplex_ptr  SA,
magma_int_t  ldsa,
magma_int_t *  info 
)
void magmablas_zlat2c_q ( magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex_const_ptr  A,
magma_int_t  lda,
magmaFloatComplex_ptr  SA,
magma_int_t  ldsa,
magma_int_t *  info,
magma_queue_t  queue 
)

ZLAT2C converts a double-complex matrix, A, to a single-complex matrix, SA.

RMAX is the overflow for the single-complex arithmetic. ZLAT2C checks that all the entries of A are between -RMAX and RMAX. If not, the conversion is aborted and a flag is raised.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix A to be converted.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
[in]nINTEGER The number of columns of the matrix A. n >= 0.
[in]ACOMPLEX_16 array, dimension (LDA,n) On entry, the n-by-n coefficient matrix A.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,n).
[out]SACOMPLEX array, dimension (LDSA,n) On exit, if INFO=0, the n-by-n coefficient matrix SA; if INFO > 0, the content of SA is unspecified.
[in]ldsaINTEGER The leading dimension of the array SA. LDSA >= max(1,n).
[out]infoINTEGER
  • = 0: successful exit.
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • = 1: an entry of the matrix A is greater than the COMPLEX overflow threshold, in this case, the content of SA on exit is unspecified.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zswapdblk ( magma_int_t  n,
magma_int_t  nb,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t  inca,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magma_int_t  incb 
)
void magmablas_zswapdblk_batched ( magma_int_t  n,
magma_int_t  nb,
magmaDoubleComplex **  dA_array,
magma_int_t  ldda,
magma_int_t  inca,
magmaDoubleComplex **  dB_array,
magma_int_t  lddb,
magma_int_t  incb,
magma_int_t  batchCount 
)
void magmablas_zswapdblk_batched_q ( magma_int_t  n,
magma_int_t  nb,
magmaDoubleComplex **  dA_array,
magma_int_t  ldda,
magma_int_t  inca,
magmaDoubleComplex **  dB_array,
magma_int_t  lddb,
magma_int_t  incb,
magma_int_t  batchCount,
magma_queue_t  queue 
)

zswapdblk swaps diagonal blocks of size nb x nb between matrices dA and dB on the GPU.

It swaps nblocks = ceil(n/nb) blocks. For i = 1 .. nblocks, submatrices dA( i*nb*inca, i*nb ) and dB( i*nb*incb, i*nb ) are swapped.

Parameters
[in]nINTEGER The number of columns of the matrices dA and dB. N >= 0.
[in]nbINTEGER The size of diagonal blocks. NB > 0 and NB <= maximum threads per CUDA block (512 or 1024).
[in,out]dACOMPLEX_16 array, dimension (LDDA,N) The matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= (nblocks - 1)*nb*inca + nb.
[in]incaINTEGER The row increment between diagonal blocks of dA. inca >= 0. For example, inca = 1 means blocks are stored on the diagonal at dA(i*nb, i*nb), inca = 0 means blocks are stored side-by-side at dA(0, i*nb).
[in,out]dBCOMPLEX_16 array, dimension (LDDB,N) The matrix dB.
[in]lddbINTEGER The leading dimension of the array db. LDDB >= (nblocks - 1)*nb*incb + nb.
[in]incbINTEGER The row increment between diagonal blocks of dB. incb >= 0. See inca.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zswapdblk_q ( magma_int_t  n,
magma_int_t  nb,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t  inca,
magmaDoubleComplex_ptr  dB,
magma_int_t  lddb,
magma_int_t  incb,
magma_queue_t  queue 
)

zswapdblk swaps diagonal blocks of size nb x nb between matrices dA and dB on the GPU.

It swaps nblocks = n/nb blocks. For i = 1 .. nblocks, submatrices dA( i*nb*inca, i*nb ) and dB( i*nb*incb, i*nb ) are swapped.

Parameters
[in]nINTEGER The number of columns of the matrices dA and dB. N >= 0.
[in]nbINTEGER The size of diagonal blocks. NB > 0 and NB <= maximum threads per CUDA block (512 or 1024).
[in,out]dACOMPLEX_16 array, dimension (LDDA,N) The matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= (nblocks - 1)*nb*inca + nb.
[in]incaINTEGER The row increment between diagonal blocks of dA. inca >= 0. For example, inca = 1 means blocks are stored on the diagonal at dA(i*nb, i*nb), inca = 0 means blocks are stored side-by-side at dA(0, i*nb).
[in,out]dBCOMPLEX_16 array, dimension (LDDB,N) The matrix dB.
[in]lddbINTEGER The leading dimension of the array db. LDDB >= (nblocks - 1)*nb*incb + nb.
[in]incbINTEGER The row increment between diagonal blocks of dB. incb >= 0. See inca.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zsymmetrize ( magma_uplo_t  uplo,
magma_int_t  m,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda 
)
void magmablas_zsymmetrize_q ( magma_uplo_t  uplo,
magma_int_t  m,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_queue_t  queue 
)

ZSYMMETRIZE copies lower triangle to upper triangle, or vice-versa, to make dA a general representation of a symmetric matrix.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix dA that is valid on input.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in,out]dACOMPLEX_16 array, dimension (LDDA,N) The m by m matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1,M).
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zsymmetrize_tiles ( magma_uplo_t  uplo,
magma_int_t  m,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t  ntile,
magma_int_t  mstride,
magma_int_t  nstride 
)
void magmablas_zsymmetrize_tiles_q ( magma_uplo_t  uplo,
magma_int_t  m,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_int_t  ntile,
magma_int_t  mstride,
magma_int_t  nstride,
magma_queue_t  queue 
)

ZSYMMETRIZE_TILES copies lower triangle to upper triangle, or vice-versa, to make some blocks of dA into general representations of a symmetric block.

This processes NTILE blocks, typically the diagonal blocks. Each block is offset by mstride rows and nstride columns from the previous block.

Parameters
[in]uplomagma_uplo_t Specifies the part of the matrix dA that is valid on input.
  • = MagmaUpper: Upper triangular part
  • = MagmaLower: Lower triangular part
[in]mINTEGER The number of rows & columns of each square block of dA. M >= 0.
[in,out]dACOMPLEX_16 array, dimension (LDDA,N) The matrix dA. N = m + nstride*(ntile-1).
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= max(1, m + mstride*(ntile-1)).
[in]ntileINTEGER Number of blocks to symmetrize. ntile >= 0.
[in]mstrideINTEGER Row offset from start of one block to start of next block. mstride >= 0. Either (mstride >= m) or (nstride >= m), to prevent m-by-m tiles from overlapping.
[in]nstrideINTEGER Column offset from start of one block to start of next block. nstride >= 0.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_ztranspose ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dAT,
magma_int_t  lddat 
)
void magmablas_ztranspose_batched ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex **  dA_array,
magma_int_t  ldda,
magmaDoubleComplex **  dAT_array,
magma_int_t  lddat,
magma_int_t  batchCount 
)
void magmablas_ztranspose_batched_q ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex **  dA_array,
magma_int_t  ldda,
magmaDoubleComplex **  dAT_array,
magma_int_t  lddat,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ztranspose_batched_q copies and transposes a matrix dA_array[i] to matrix dAT_array[i].

Same as ztranspose_batched, but adds queue argument.

Parameters
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]dA_arrayCOMPLEX_16* array, dimension (batchCount) array of pointers to the matrices dA, where each dA is of dimension (LDDA,N) The M-by-N matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= M.
[in]dAT_arrayCOMPLEX_16* array, dimension (batchCount) array of pointers to the matrices dAT, where each dAT is of dimension (LDDAT,M) The N-by-M matrix dAT.
[in]lddatINTEGER The leading dimension of the array dAT. LDDAT >= N.
[in]queuemagma_queue_t Queue to execute in.
[in]batchCountNumber of matrices in dA_array and dAT_array
void magmablas_ztranspose_conj ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dAT,
magma_int_t  lddat 
)
void magmablas_ztranspose_conj_batched ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex **  dA_array,
magma_int_t  ldda,
magmaDoubleComplex **  dAT_array,
magma_int_t  lddat,
magma_int_t  batchCount 
)
void magmablas_ztranspose_conj_batched_q ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex **  dA_array,
magma_int_t  ldda,
magmaDoubleComplex **  dAT_array,
magma_int_t  lddat,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ztranspose_conj_batched_q copies and conjugate-transposes a matrix dA_array[i] to matrix dAT_array[i].

Same as ztranspose_conj_batched, but adds queue argument.

Parameters
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]dA_arrayCOMPLEX_16* array, dimension (batchCount) array of pointers to the matrices dA, where each dA is of dimension (LDDA,N) The M-by-N matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= M.
[in]dAT_arrayCOMPLEX_16* array, dimension (batchCount) array of pointers to the matrices dAT, where each dAT is of dimension (LDDAT,M) The N-by-M matrix dAT.
[in]lddatINTEGER The leading dimension of the array dAT. LDDAT >= N.
[in]queuemagma_queue_t Queue to execute in.
[in]batchCountNumber of matrices in dA_array and dAT_array
void magmablas_ztranspose_conj_inplace ( magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda 
)
void magmablas_ztranspose_conj_inplace_q ( magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_queue_t  queue 
)

ztranspose_conj_inplace_q conjugate-transposes a square N-by-N matrix in-place.

Same as ztranspose_conj_inplace, but adds queue argument.

Parameters
[in]nINTEGER The number of rows & columns of the matrix dA. N >= 0.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The N-by-N matrix dA. On exit, dA(j,i) = dA_original(i,j), for 0 <= i,j < N.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= N.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_ztranspose_conj_q ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dAT,
magma_int_t  lddat,
magma_queue_t  queue 
)

ztranspose_conj_q copies and conjugate-transposes a matrix dA to matrix dAT.

Same as ztranspose_conj, but adds queue argument.

Parameters
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The M-by-N matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= M.
[in]dATCOMPLEX_16 array, dimension (LDDAT,M) The N-by-M matrix dAT.
[in]lddatINTEGER The leading dimension of the array dAT. LDDAT >= N.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_ztranspose_inplace ( magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda 
)
void magmablas_ztranspose_inplace_q ( magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
magma_queue_t  queue 
)

ztranspose_inplace_q transposes a square N-by-N matrix in-place.

Same as ztranspose_inplace, but adds queue argument.

Parameters
[in]nINTEGER The number of rows & columns of the matrix dA. N >= 0.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The N-by-N matrix dA. On exit, dA(j,i) = dA_original(i,j), for 0 <= i,j < N.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= N.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_ztranspose_q ( magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_ptr  dAT,
magma_int_t  lddat,
magma_queue_t  queue 
)

ztranspose_q copies and transposes a matrix dA to matrix dAT.

Same as ztranspose, but adds queue argument.

Parameters
[in]mINTEGER The number of rows of the matrix dA. M >= 0.
[in]nINTEGER The number of columns of the matrix dA. N >= 0.
[in]dACOMPLEX_16 array, dimension (LDDA,N) The M-by-N matrix dA.
[in]lddaINTEGER The leading dimension of the array dA. LDDA >= M.
[in]dATCOMPLEX_16 array, dimension (LDDAT,M) The N-by-M matrix dAT.
[in]lddatINTEGER The leading dimension of the array dAT. LDDAT >= N.
[in]queuemagma_queue_t Queue to execute in.