MAGMA  2.7.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
syr2k: Symmetric rank 2k update

\( C = \alpha A B^T + \alpha B A^T + \beta C \) where \( C \) is symmetric More...

Functions

void magmablas_csyr2k_batched (magma_uplo_t uplo, magma_trans_t trans, magma_int_t n, magma_int_t k, magmaFloatComplex alpha, magmaFloatComplex const *const *dA_array, magma_int_t ldda, magmaFloatComplex const *const *dB_array, magma_int_t lddb, magmaFloatComplex beta, magmaFloatComplex **dC_array, magma_int_t lddc, magma_int_t batchCount, magma_queue_t queue)
 CSYR2K performs one of the symmetric rank 2k operations. More...
 
void magmablas_csyr2k_vbatched (magma_uplo_t uplo, magma_trans_t trans, magma_int_t *n, magma_int_t *k, magmaFloatComplex alpha, magmaFloatComplex const *const *dA_array, magma_int_t *ldda, magmaFloatComplex const *const *dB_array, magma_int_t *lddb, magmaFloatComplex beta, magmaFloatComplex **dC_array, magma_int_t *lddc, magma_int_t batchCount, magma_queue_t queue)
 CSYR2K performs one of the symmetric rank 2k operations. More...
 
void magmablas_zsyr2k_batched (magma_uplo_t uplo, magma_trans_t trans, magma_int_t n, magma_int_t k, magmaDoubleComplex alpha, magmaDoubleComplex const *const *dA_array, magma_int_t ldda, magmaDoubleComplex const *const *dB_array, magma_int_t lddb, magmaDoubleComplex beta, magmaDoubleComplex **dC_array, magma_int_t lddc, magma_int_t batchCount, magma_queue_t queue)
 ZSYR2K performs one of the symmetric rank 2k operations. More...
 
void magmablas_zsyr2k_vbatched (magma_uplo_t uplo, magma_trans_t trans, magma_int_t *n, magma_int_t *k, magmaDoubleComplex alpha, magmaDoubleComplex const *const *dA_array, magma_int_t *ldda, magmaDoubleComplex const *const *dB_array, magma_int_t *lddb, magmaDoubleComplex beta, magmaDoubleComplex **dC_array, magma_int_t *lddc, magma_int_t batchCount, magma_queue_t queue)
 ZSYR2K performs one of the symmetric rank 2k operations. More...
 

Detailed Description

\( C = \alpha A B^T + \alpha B A^T + \beta C \) where \( C \) is symmetric

Function Documentation

void magmablas_csyr2k_batched ( magma_uplo_t  uplo,
magma_trans_t  trans,
magma_int_t  n,
magma_int_t  k,
magmaFloatComplex  alpha,
magmaFloatComplex const *const *  dA_array,
magma_int_t  ldda,
magmaFloatComplex const *const *  dB_array,
magma_int_t  lddb,
magmaFloatComplex  beta,
magmaFloatComplex **  dC_array,
magma_int_t  lddc,
magma_int_t  batchCount,
magma_queue_t  queue 
)

CSYR2K performs one of the symmetric rank 2k operations.

C := alpha*A*B**H + conjg( alpha )*B*A**H + beta*C,

or

C := alpha*A**H*B + conjg( alpha )*B**H*A + beta*C,

where alpha and beta are scalars with beta real, C is an n by n symmetric matrix and A and B are n by k matrices in the first case and k by n matrices in the second case.

Parameters

Parameters
[in]uplomagma_uplo_t. On entry, UPLO specifies whether the upper or lower triangular part of the array C is to be referenced as follows:
  • = MagmaUpper: Only the upper triangular part of C is to be referenced.
  • = MagmaLower: Only the lower triangular part of C is to be referenced.
[in]transmagma_trans_t. On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: C := alpha*A*B**H + conj( alpha )*B*A**H + beta*C.
  • = Magma_ConjTrans: C := alpha*A**H*B + conj( alpha )*B**H*A + beta*C.
[in]nINTEGER. On entry, N specifies the order of the matrix C. N must be at least zero.
[in]kINTEGER. On entry with TRANS = MagmaNoTrans, k specifies the number of columns of the matrices A and B, and on entry with TRANS = Magma_ConjTrans, k specifies the number of rows of the matrices A and B. k must be at least zero.
[in]alphaCOMPLEX. On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of DIMENSION ( ldda, ka ), where ka is k when TRANS = MagmaNoTrans, and is n otherwise. Before entry with TRANS = MagmaNoTrans, the leading n by k part of the array A must contain the matrix A, otherwise the leading k by n part of the array A must contain the matrix A.
[in]lddaINTEGER. On entry, ldda specifies the first dimension of A as declared in the calling (sub) program. When TRANS = MagmaNoTrans then ldda must be at least max( 1, n ), otherwise ldda must be at least max( 1, k ).
[in]dB_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of DIMENSION ( ldb, kb ), where kb is k when TRANS = MagmaNoTrans, and is n otherwise. Before entry with TRANS = MagmaNoTrans, the leading n by k part of the array B must contain the matrix B, otherwise the leading k by n part of the array B must contain the matrix B.
[in]lddbINTEGER On entry, lddb specifies the first dimension of B as declared in the calling (sub) program. When TRANS = MagmaNoTrans then lddb must be at least max( 1, n ), otherwise lddb must be at least max( 1, k ). Unchanged on exit.
[in]betaCOMPLEX. On entry, BETA specifies the scalar beta.
[in,out]dC_arrayArray of pointers, dimension (batchCount). Each is COMPLEX array of DIMENSION ( lddc, n ). Before entry with UPLO = MagmaUpper, the leading n by n upper triangular part of the array C must contain the upper triangular part of the symmetric matrix and the strictly lower triangular part of C is not referenced. On exit, the upper triangular part of the array C is overwritten by the upper triangular part of the updated matrix. Before entry with UPLO = MagmaLower, the leading n by n lower triangular part of the array C must contain the lower triangular part of the symmetric matrix and the strictly upper triangular part of C is not referenced. On exit, the lower triangular part of the array C is overwritten by the lower triangular part of the updated matrix. Note that the imaginary parts of the diagonal elements need not be set, they are assumed to be zero, and on exit they are set to zero.
[in]lddcINTEGER. On entry, lddc specifies the first dimension of each array C as declared in the calling (sub) program. lddc must be at least max( 1, n ).
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_csyr2k_vbatched ( magma_uplo_t  uplo,
magma_trans_t  trans,
magma_int_t *  n,
magma_int_t *  k,
magmaFloatComplex  alpha,
magmaFloatComplex const *const *  dA_array,
magma_int_t *  ldda,
magmaFloatComplex const *const *  dB_array,
magma_int_t *  lddb,
magmaFloatComplex  beta,
magmaFloatComplex **  dC_array,
magma_int_t *  lddc,
magma_int_t  batchCount,
magma_queue_t  queue 
)

CSYR2K performs one of the symmetric rank 2k operations.

C := alpha*A*B**T + conjg( alpha )*B*A**T + beta*C,

or

C := alpha*A**T*B + conjg( alpha )*B**T*A + beta*C,

where alpha and beta are scalars with beta real, C is an n by n symmetric matrix and A and B are n by k matrices in the first case and k by n matrices in the second case.

Parameters

Parameters
[in]uplomagma_uplo_t. On entry, UPLO specifies whether the upper or lower triangular part of the array C is to be referenced as follows:
  • = MagmaUpper: Only the upper triangular part of C is to be referenced.
  • = MagmaLower: Only the lower triangular part of C is to be referenced.
[in]transmagma_trans_t. On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: C := alpha*A*B**T + conj( alpha )*B*A**T + beta*C.
  • = Magma_ConjTrans: C := alpha*A**T*B + conj( alpha )*B**T*A + beta*C.
[in]nArray of integers, size(batchCount + 1). On entry, each INTEGER N specifies the order of the corresponding matrix C. N must be at least zero. The last element of the array is used internally by the routine.
[in]kArray of integers, size(batchCount + 1). On entry with TRANS = MagmaNoTrans, each INTEGER K specifies the number of columns of the corresponding matrices A and B, and on entry with TRANS = Magma_ConjTrans, K specifies the number of rows of the corresponding matrices A and B. k must be at least zero.
[in]alphaCOMPLEX. On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of DIMENSION ( LDDA, Ka ), where Ka is K when TRANS = MagmaNoTrans, and is N otherwise. Before entry with TRANS = MagmaNoTrans, the leading N by K part of the array A must contain the matrix A, otherwise the leading K by N part of the array A must contain the matrix A.
[in]lddaArray of integers, size(batchCount + 1). On entry, each INTEGER LDDA specifies the first dimension of the corresponding matrix A as declared in the calling (sub) program. When TRANS = MagmaNoTrans then LDDA must be at least max( 1, N ), otherwise LDDA must be at least max( 1, K ). The last element of the array is used internally by the routine.
[in]dB_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of DIMENSION ( LDDB, Kb ), where Kb is K when TRANS = MagmaNoTrans, and is N otherwise. Before entry with TRANS = MagmaNoTrans, the leading N by K part of the array B must contain the matrix B, otherwise the leading K by N part of the array B must contain the matrix B.
[in]lddbArray of integers, size(batchCount + 1). On entry, each INTEGER LDDB specifies the first dimension of the corresponding matrix B as declared in the calling (sub) program. When TRANS = MagmaNoTrans then LDDB must be at least max( 1, N ), otherwise LDDB must be at least max( 1, K ). Unchanged on exit. The last element of the array is used internally by the routine.
[in]betaCOMPLEX. On entry, BETA specifies the scalar beta.
[in,out]dC_arrayArray of pointers, dimension (batchCount). Each is COMPLEX array of DIMENSION ( lddc, n ). Before entry with UPLO = MagmaUpper, the leading n by n upper triangular part of the array C must contain the upper triangular part of the symmetric matrix and the strictly lower triangular part of C is not referenced. On exit, the upper triangular part of the array C is overwritten by the upper triangular part of the updated matrix. Before entry with UPLO = MagmaLower, the leading n by n lower triangular part of the array C must contain the lower triangular part of the symmetric matrix and the strictly upper triangular part of C is not referenced. On exit, the lower triangular part of the array C is overwritten by the lower triangular part of the updated matrix. Note that the imaginary parts of the diagonal elements need not be set, they are assumed to be zero, and on exit they are set to zero.
[in]lddcArray of integers, size(batchCount + 1). On entry, each INTEGER LDDC specifies the first dimension of each array C as declared in the calling (sub) program. LDDC must be at least max( 1, N ). The last element of the array is used internally by the routine.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zsyr2k_batched ( magma_uplo_t  uplo,
magma_trans_t  trans,
magma_int_t  n,
magma_int_t  k,
magmaDoubleComplex  alpha,
magmaDoubleComplex const *const *  dA_array,
magma_int_t  ldda,
magmaDoubleComplex const *const *  dB_array,
magma_int_t  lddb,
magmaDoubleComplex  beta,
magmaDoubleComplex **  dC_array,
magma_int_t  lddc,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ZSYR2K performs one of the symmetric rank 2k operations.

C := alpha*A*B**H + conjg( alpha )*B*A**H + beta*C,

or

C := alpha*A**H*B + conjg( alpha )*B**H*A + beta*C,

where alpha and beta are scalars with beta real, C is an n by n symmetric matrix and A and B are n by k matrices in the first case and k by n matrices in the second case.

Parameters

Parameters
[in]uplomagma_uplo_t. On entry, UPLO specifies whether the upper or lower triangular part of the array C is to be referenced as follows:
  • = MagmaUpper: Only the upper triangular part of C is to be referenced.
  • = MagmaLower: Only the lower triangular part of C is to be referenced.
[in]transmagma_trans_t. On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: C := alpha*A*B**H + conj( alpha )*B*A**H + beta*C.
  • = Magma_ConjTrans: C := alpha*A**H*B + conj( alpha )*B**H*A + beta*C.
[in]nINTEGER. On entry, N specifies the order of the matrix C. N must be at least zero.
[in]kINTEGER. On entry with TRANS = MagmaNoTrans, k specifies the number of columns of the matrices A and B, and on entry with TRANS = Magma_ConjTrans, k specifies the number of rows of the matrices A and B. k must be at least zero.
[in]alphaCOMPLEX*16. On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX*16 array of DIMENSION ( ldda, ka ), where ka is k when TRANS = MagmaNoTrans, and is n otherwise. Before entry with TRANS = MagmaNoTrans, the leading n by k part of the array A must contain the matrix A, otherwise the leading k by n part of the array A must contain the matrix A.
[in]lddaINTEGER. On entry, ldda specifies the first dimension of A as declared in the calling (sub) program. When TRANS = MagmaNoTrans then ldda must be at least max( 1, n ), otherwise ldda must be at least max( 1, k ).
[in]dB_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX*16 array of DIMENSION ( ldb, kb ), where kb is k when TRANS = MagmaNoTrans, and is n otherwise. Before entry with TRANS = MagmaNoTrans, the leading n by k part of the array B must contain the matrix B, otherwise the leading k by n part of the array B must contain the matrix B.
[in]lddbINTEGER On entry, lddb specifies the first dimension of B as declared in the calling (sub) program. When TRANS = MagmaNoTrans then lddb must be at least max( 1, n ), otherwise lddb must be at least max( 1, k ). Unchanged on exit.
[in]betaCOMPLEX*16. On entry, BETA specifies the scalar beta.
[in,out]dC_arrayArray of pointers, dimension (batchCount). Each is COMPLEX*16 array of DIMENSION ( lddc, n ). Before entry with UPLO = MagmaUpper, the leading n by n upper triangular part of the array C must contain the upper triangular part of the symmetric matrix and the strictly lower triangular part of C is not referenced. On exit, the upper triangular part of the array C is overwritten by the upper triangular part of the updated matrix. Before entry with UPLO = MagmaLower, the leading n by n lower triangular part of the array C must contain the lower triangular part of the symmetric matrix and the strictly upper triangular part of C is not referenced. On exit, the lower triangular part of the array C is overwritten by the lower triangular part of the updated matrix. Note that the imaginary parts of the diagonal elements need not be set, they are assumed to be zero, and on exit they are set to zero.
[in]lddcINTEGER. On entry, lddc specifies the first dimension of each array C as declared in the calling (sub) program. lddc must be at least max( 1, n ).
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.
void magmablas_zsyr2k_vbatched ( magma_uplo_t  uplo,
magma_trans_t  trans,
magma_int_t *  n,
magma_int_t *  k,
magmaDoubleComplex  alpha,
magmaDoubleComplex const *const *  dA_array,
magma_int_t *  ldda,
magmaDoubleComplex const *const *  dB_array,
magma_int_t *  lddb,
magmaDoubleComplex  beta,
magmaDoubleComplex **  dC_array,
magma_int_t *  lddc,
magma_int_t  batchCount,
magma_queue_t  queue 
)

ZSYR2K performs one of the symmetric rank 2k operations.

C := alpha*A*B**T + conjg( alpha )*B*A**T + beta*C,

or

C := alpha*A**T*B + conjg( alpha )*B**T*A + beta*C,

where alpha and beta are scalars with beta real, C is an n by n symmetric matrix and A and B are n by k matrices in the first case and k by n matrices in the second case.

Parameters

Parameters
[in]uplomagma_uplo_t. On entry, UPLO specifies whether the upper or lower triangular part of the array C is to be referenced as follows:
  • = MagmaUpper: Only the upper triangular part of C is to be referenced.
  • = MagmaLower: Only the lower triangular part of C is to be referenced.
[in]transmagma_trans_t. On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: C := alpha*A*B**T + conj( alpha )*B*A**T + beta*C.
  • = Magma_ConjTrans: C := alpha*A**T*B + conj( alpha )*B**T*A + beta*C.
[in]nArray of integers, size(batchCount + 1). On entry, each INTEGER N specifies the order of the corresponding matrix C. N must be at least zero. The last element of the array is used internally by the routine.
[in]kArray of integers, size(batchCount + 1). On entry with TRANS = MagmaNoTrans, each INTEGER K specifies the number of columns of the corresponding matrices A and B, and on entry with TRANS = Magma_ConjTrans, K specifies the number of rows of the corresponding matrices A and B. k must be at least zero.
[in]alphaCOMPLEX*16. On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX*16 array of DIMENSION ( LDDA, Ka ), where Ka is K when TRANS = MagmaNoTrans, and is N otherwise. Before entry with TRANS = MagmaNoTrans, the leading N by K part of the array A must contain the matrix A, otherwise the leading K by N part of the array A must contain the matrix A.
[in]lddaArray of integers, size(batchCount + 1). On entry, each INTEGER LDDA specifies the first dimension of the corresponding matrix A as declared in the calling (sub) program. When TRANS = MagmaNoTrans then LDDA must be at least max( 1, N ), otherwise LDDA must be at least max( 1, K ). The last element of the array is used internally by the routine.
[in]dB_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX*16 array of DIMENSION ( LDDB, Kb ), where Kb is K when TRANS = MagmaNoTrans, and is N otherwise. Before entry with TRANS = MagmaNoTrans, the leading N by K part of the array B must contain the matrix B, otherwise the leading K by N part of the array B must contain the matrix B.
[in]lddbArray of integers, size(batchCount + 1). On entry, each INTEGER LDDB specifies the first dimension of the corresponding matrix B as declared in the calling (sub) program. When TRANS = MagmaNoTrans then LDDB must be at least max( 1, N ), otherwise LDDB must be at least max( 1, K ). Unchanged on exit. The last element of the array is used internally by the routine.
[in]betaCOMPLEX*16. On entry, BETA specifies the scalar beta.
[in,out]dC_arrayArray of pointers, dimension (batchCount). Each is COMPLEX*16 array of DIMENSION ( lddc, n ). Before entry with UPLO = MagmaUpper, the leading n by n upper triangular part of the array C must contain the upper triangular part of the symmetric matrix and the strictly lower triangular part of C is not referenced. On exit, the upper triangular part of the array C is overwritten by the upper triangular part of the updated matrix. Before entry with UPLO = MagmaLower, the leading n by n lower triangular part of the array C must contain the lower triangular part of the symmetric matrix and the strictly upper triangular part of C is not referenced. On exit, the lower triangular part of the array C is overwritten by the lower triangular part of the updated matrix. Note that the imaginary parts of the diagonal elements need not be set, they are assumed to be zero, and on exit they are set to zero.
[in]lddcArray of integers, size(batchCount + 1). On entry, each INTEGER LDDC specifies the first dimension of each array C as declared in the calling (sub) program. LDDC must be at least max( 1, N ). The last element of the array is used internally by the routine.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.