herk#
Performs a Hermitian rank-k update.
Description
The herk
routines compute a rank-k update of a Hermitian matrix
C
by a general matrix A
. The operation is defined as:
where:
op(X
) is one of op(X
) = X
or op(X
) = X
H,
alpha
and beta
are real scalars,
C
is a Hermitian matrix and A
is a general matrix.
Here op(A
) is n
x k
, and C
is n
x n
.
herk
supports the following precisions:
T
Treal
std::complex<float>
float
std::complex<double>
double
herk (Buffer Version)#
Syntax
namespace oneapi::math::blas::column_major {
void herk(sycl::queue &queue,
oneapi::math::uplo upper_lower,
oneapi::math::transpose trans,
std::int64_t n,
std::int64_t k,
Treal alpha,
sycl::buffer<T,1> &a,
std::int64_t lda,
Treal beta,
sycl::buffer<T,1> &c,
std::int64_t ldc)
}
namespace oneapi::math::blas::row_major {
void herk(sycl::queue &queue,
oneapi::math::uplo upper_lower,
oneapi::math::transpose trans,
std::int64_t n,
std::int64_t k,
Treal alpha,
sycl::buffer<T,1> &a,
std::int64_t lda,
Treal 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 oneMath defined datatypes for more details.- trans
Specifies op(
A
), the transposition operation applied toA
. See oneMath defined datatypes for more details. Supported operations aretranspose::nontrans
andtranspose::conjtrans
.- n
The number of rows and columns in
C
.The value ofn
must be at least zero.- k
Number of columns in op(
A
).The value of
k
must be at least zero.- alpha
Real scaling factor for the rank-k update.
- a
Buffer holding input matrix
A
.trans
=transpose::nontrans
trans
=transpose::trans
ortranspose::conjtrans
Column major
A
is ann
-by-k
matrix so the arraya
must have size at leastlda
*k
.A
is ank
-by-n
matrix so the arraya
must have size at leastlda
*n
Row major
A
is ann
-by-k
matrix so the arraya
must have size at leastlda
*n
.A
is ank
-by-n
matrix so the arraya
must have size at leastlda
*k
.See Matrix Storage for more details.
- lda
The leading dimension of
A
. It must be positive.trans
=transpose::nontrans
trans
=transpose::trans
ortranspose::conjtrans
Column major
lda
must be at leastn
.lda
must be at leastk
.Row major
lda
must be at leastk
.lda
must be at leastn
.- 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
The output buffer, overwritten by
alpha
*op(A
)*op(A
)T +beta
*C
. The imaginary parts of the diagonal elements are set to zero.
Throws
This routine shall throw the following exceptions if the associated condition is detected. An implementation may throw additional implementation-specific exception(s) in case of error conditions not covered here.
oneapi::math::invalid_argument
oneapi::math::unsupported_device
herk (USM Version)#
Syntax
namespace oneapi::math::blas::column_major {
sycl::event herk(sycl::queue &queue,
oneapi::math::uplo upper_lower,
oneapi::math::transpose trans,
std::int64_t n,
std::int64_t k,
value_or_pointer<Treal> alpha,
const T *a,
std::int64_t lda,
value_or_pointer<Treal> beta,
T *c,
std::int64_t ldc,
const std::vector<sycl::event> &dependencies = {})
}
namespace oneapi::math::blas::row_major {
sycl::event herk(sycl::queue &queue,
oneapi::math::uplo upper_lower,
oneapi::math::transpose trans,
std::int64_t n,
std::int64_t k,
value_or_pointer<Treal> alpha,
const T *a,
std::int64_t lda,
value_or_pointer<Treal> 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 oneMath defined datatypes for more details.- trans
Specifies op(
A
), the transposition operation applied toA
. See oneMath defined datatypes for more details. Supported operations aretranspose::nontrans
andtranspose::conjtrans
.- n
The number of rows and columns in
C
.The value ofn
must be at least zero.- k
Number of columns in op(
A
).The value of
k
must be at least zero.- alpha
Real scaling factor for the rank-k update. See Scalar Arguments in BLAS for more details.
- a
Pointer to input matrix
A
.trans
=transpose::nontrans
trans
=transpose::trans
ortranspose::conjtrans
Column major
A
is ann
-by-k
matrix so the arraya
must have size at leastlda
*k
.A
is ank
-by-n
matrix so the arraya
must have size at leastlda
*n
Row major
A
is ann
-by-k
matrix so the arraya
must have size at leastlda
*n
.A
is ank
-by-n
matrix so the arraya
must have size at leastlda
*k
.See Matrix Storage for more details.
- lda
The leading dimension of
A
. It must be positive.trans
=transpose::nontrans
trans
=transpose::trans
ortranspose::conjtrans
Column major
lda
must be at leastn
.lda
must be at leastk
.Row major
lda
must be at leastk
.lda
must be at leastn
.- beta
Real scaling factor for matrix
C
. See Scalar Arguments in BLAS for more details.- 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
alpha
*op(A
)*op(A
)T +beta
*C
. The imaginary parts of the diagonal elements are set to zero.
Return Values
Output event to wait on to ensure computation is complete.
Throws
This routine shall throw the following exceptions if the associated condition is detected. An implementation may throw additional implementation-specific exception(s) in case of error conditions not covered here.
oneapi::math::invalid_argument
oneapi::math::unsupported_device
oneapi::math::device_bad_alloc
Parent topic: BLAS Level 3 Routines