her2k#
Performs a Hermitian rank-2k update.
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:
where A is n x k and B is k x n.
If trans = transpose::conjtrans, the operation is defined as:
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 supports the following precisions:
T
T_real
std::complex<float>
float
std::complex<double>
double
her2k (Buffer Version)#
Syntax
namespace oneapi::mkl::blas::column_major {
void 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)
}
namespace oneapi::mkl::blas::row_major {
void 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::nontransandtranspose::conjtrans.- n
The number of rows and columns in
C. The value ofnmust be at least zero.- k
The inner dimension of matrix multiplications. The value of
kmust be at least equal to zero.- alpha
Complex scaling factor for the rank-2k update.
- a
Buffer holding input matrix
A.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*k.Ais ank-by-nmatrix so the arrayamust have size at leastlda*nRow major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*n.Ais ank-by-nmatrix so the arrayamust have size at leastlda*k.See Matrix Storage for more details.
- lda
The leading dimension of
A. It must be positive.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
ldamust be at leastn.ldamust be at leastk.Row major
ldamust be at leastk.ldamust be at leastn.- b
Buffer holding input matrix
B.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*n.Bis ann-by-kmatrix so the arraybmust have size at leastldb*kRow major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*k.Bis ann-by-kmatrix so the arraybmust have size at leastldb*n.See Matrix Storage for more details.
- ldb
The leading dimension of
B. It must be positive.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
ldbmust be at leastk.ldbmust be at leastn.Row major
ldbmust be at leastn.ldbmust be at leastk.- beta
Real scaling factor for matrix
C.- c
Buffer holding input/output matrix
C. Must have size at leastldc*n. See Matrix Storage for more details.- ldc
Leading dimension of
C. Must be positive and at leastn.
Output Parameters
- c
Output buffer, overwritten by the updated
Cmatrix.
her2k (USM Version)#
Syntax
namespace oneapi::mkl::blas::column_major {
sycl::event 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 std::vector<sycl::event> &dependencies = {})
}
namespace oneapi::mkl::blas::row_major {
sycl::event 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 std::vector<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::nontransandtranspose::conjtrans.- n
The number of rows and columns in
C. The value ofnmust be at least zero.- k
The inner dimension of matrix multiplications. The value of
kmust be at least equal to zero.- alpha
Complex scaling factor for the rank-2k update.
- a
Pointer to input matrix
A.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*k.Ais ank-by-nmatrix so the arrayamust have size at leastlda*nRow major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*n.Ais ank-by-nmatrix so the arrayamust have size at leastlda*k.See Matrix Storage for more details.
- lda
The leading dimension of
A. It must be positive.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
ldamust be at leastn.ldamust be at leastk.Row major
ldamust be at leastk.ldamust be at leastn.- b
Pointer to input matrix
B.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*n.Bis ann-by-kmatrix so the arraybmust have size at leastldb*kRow major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*k.Bis ann-by-kmatrix so the arraybmust have size at leastldb*n.See Matrix Storage for more details.
- ldb
The leading dimension of
B. It must be positive.trans=transpose::nontranstrans=transpose::transortranspose::conjtransColumn major
ldbmust be at leastk.ldbmust be at leastn.Row major
ldbmust be at leastn.ldbmust be at leastk.- beta
Real scaling factor for matrix
C.- c
Pointer to input/output matrix
C. Must have size at leastldc*n. See Matrix 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
Cmatrix.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: BLAS Level 3 Routines