her2k¶
Performs a Hermitian rank-2k update.
her2k
supports the following precisions:
T
T_real
std::complex<float>
float
std::complex<double>
double
Description
The her2k
routines perform a rank-2k update of an n
x n
Hermitian matrix C
by general matrices A
and B
. If
trans
= transpose::nontrans
. The operation is defined as
C <- alpha*A*B H + conjg(alpha)*B*A H + beta*C
where A
is n
x k
and B
is k
x n
.
If trans
= transpose::conjtrans
, the operation is defined as:
C <- alpha*B*A H + conjg(alpha)*A*B H + beta*C
where A
is k
x n
and B
is n
x k
.
In both cases:
alpha
is a complex scalar and beta
is a real scalar.
C
is a Hermitian matrix and A, B
are general matrices.
The inner dimension of both matrix multiplications is k
.
her2k (Buffer Version)¶
Syntax
-
void
onemkl::blas
::
her2k
(sycl::queue &queue, onemkl::uplo upper_lower, onemkl::transpose trans, std::int64_t n, std::int64_t k, T alpha, sycl::buffer<T, 1> &a, std::int64_t lda, sycl::buffer<T, 1> &b, std::int64_t ldb, T_real beta, sycl::buffer<T, 1> &c, std::int64_t ldc)¶
Input Parameters
- queue
The queue where the routine should be executed.
- upper_lower
Specifies whether
A
’s data is stored in its upper or lower triangle. See oneMKL defined datatypes for more details.- trans
Specifies the operation to apply, as described above. Supported operations are
transpose::nontrans
andtranspose::conjtrans
.- n
The number of rows and columns in
C
. The value ofn
must be at least zero.- k
The inner dimension of matrix multiplications. The value of
k
must be at least equal to zero.- alpha
Complex scaling factor for the rank-2
k
update.- a
Buffer holding input matrix
A
. Iftrans
=transpose::nontrans
,A
is ann
-by-k
matrix so the arraya
must have size at leastlda
*k
. Otherwise,A
is ank
-by-n
matrix so the arraya
must have size at leastlda
*n
. See Matrix and Vector Storage for more details.- lda
Leading dimension of
A
. Must be at leastn
iftrans
=transpose::nontrans
, and at leastk
otherwise. Must be positive.- beta
Real scaling factor for matrix
C
.- b
Buffer holding input matrix
B
. Iftrans
=transpose::nontrans
,B
is ank
-by-n
matrix so the arrayb
must have size at leastldb
*n
. Otherwise,B
is ann
-by-k
matrix so the arrayb
must have size at leastldb
*k
. See Matrix and Vector Storage for more details.- ldb
Leading dimension of
B
. Must be at leastk
iftrans
=transpose::nontrans
, and at leastn
otherwise. Must be positive.- c
Buffer holding input/output matrix
C
. Must have size at leastldc
*n
. See Matrix and Vector Storage for more details.- ldc
Leading dimension of
C
. Must be positive and at leastn
.
Output Parameters
- c
Output buffer, overwritten by the updated
C
matrix.
her2k (USM Version)¶
Syntax
-
sycl::event
onemkl::blas
::
her2k
(sycl::queue &queue, onemkl::uplo upper_lower, onemkl::transpose trans, std::int64_t n, std::int64_t k, T alpha, const T *a, std::int64_t lda, const T *b, std::int64_t ldb, T_real beta, T *c, std::int64_t ldc, const sycl::vector_class<sycl::event> &dependencies = {})¶
Input Parameters
- queue
The queue where the routine should be executed.
- upper_lower
Specifies whether
A
’s data is stored in its upper or lower triangle. See oneMKL defined datatypes for more details.- trans
Specifies the operation to apply, as described above. Supported operations are
transpose::nontrans
andtranspose::conjtrans
.- n
The number of rows and columns in
C
. The value ofn
must be at least zero.- k
The inner dimension of matrix multiplications. The value of
k
must be at least equal to zero.- alpha
Complex scaling factor for the rank-2
k
update.- a
Pointer to input matrix
A
. Iftrans
=transpose::nontrans
,A
is ann
-by-k
matrix so the arraya
must have size at leastlda
*k
. Otherwise,A
is ank
-by-n
matrix so the arraya
must have size at leastlda
*n
. See Matrix and Vector Storage for more details.- lda
Leading dimension of
A
. Must be at leastn
iftrans
=transpose::nontrans
, and at leastk
otherwise. Must be positive.- beta
Real scaling factor for matrix
C
.- b
Pointer to input matrix
B
. Iftrans
=transpose::nontrans
,B
is ank
-by-n
matrix so the arrayb
must have size at leastldb
*n
. Otherwise,B
is ann
-by-k
matrix so the arrayb
must have size at leastldb
*k
. See Matrix and Vector Storage for more details.- ldb
Leading dimension of
B
. Must be at leastk
iftrans
=transpose::nontrans
, and at leastn
otherwise. Must be positive.- c
Pointer to input/output matrix
C
. Must have size at leastldc
*n
. See Matrix and Vector Storage for more details.- ldc
Leading dimension of
C
. Must be positive and at leastn
.- dependencies
List of events to wait for before starting computation, if any. If omitted, defaults to no dependencies.
Output Parameters
- c
Pointer to the output matrix, overwritten by the updated
C
matrix.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: BLAS Level 3 Routines