unmtr#

Multiplies a complex matrix by the complex unitary matrix Q determined by hetrd.

Description

unmtr supports the following precisions.

T

std::complex<float>

std::complex<double>

The routine multiplies a complex matrix \(C\) by \(Q\) or \(Q^{H}\), where \(Q\) is the unitary matrix \(Q\) formed by hetrd when reducing a complex Hermitian matrix \(A\) to tridiagonal form: \(A = QTQ^H\). Use this routine after a call to hetrd.

Depending on the parameters side and trans, the routine can form one of the matrix products \(QC\), \(Q^{H}C\), \(CQ\), or \(CQ^{H}\) (overwriting the result on \(C\)).

unmtr (Buffer Version)#

Syntax

namespace oneapi::mkl::lapack {
  void unmtr(sycl::queue &queue, oneapi::mkl::side side, oneapi::mkl::uplo upper_lower, oneapi::mkl::transpose trans, std::int64_t m, std::int64_t n, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<T,1> &tau, sycl::buffer<T,1> &c, std::int64_t ldc, sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}

Input Parameters

In the descriptions below, r denotes the order of \(Q\):

\(r\)=\(m\)

if side = side::left

\(r\)=\(n\)

if side = side::right

queue

The queue where the routine should be executed.

side

Must be either side::left or side::right.

If side=side::left, \(Q\) or \(Q^{H}\) is applied to \(C\) from the left.

If side=side::right, \(Q\) or \(Q^{H}\) is applied to \(C\) from the right.

upper_lower

Must be either uplo::upper or uplo::lower. Uses the same upper_lower as supplied to hetrd.

trans

Must be either transpose::nontrans or transpose::conjtrans.

If trans=transpose::nontrans, the routine multiplies \(C\) by \(Q\).

If trans=transpose::conjtrans, the routine multiplies \(C\) by \(Q^{H}\).

m

The number of rows in the matrix \(C\) (\(m \ge 0\)).

n

The number of columns the matrix \(C\) (\(n \ge 0\)).

k

The number of elementary reflectors whose product defines the matrix \(Q\) (\(0 \le k \le n\)).

a

The buffer a as returned by hetrd.

lda

The leading dimension of a \((\max(1,r) \le \text{lda})\).

tau

The buffer tau as returned by hetrd. The dimension of tau must be at least \(\max(1,r-1)\).

c

The buffer c contains the matrix \(C\). The second dimension of c must be at least \(\max(1,n)\).

ldc

The leading dimension of c \((\max(1,n) \le \text{ldc})\).

scratchpad_size

Size of scratchpad memory as a number of floating point elements of type T. Size should not be less than the value returned by unmtr_scratchpad_size function.

Output Parameters

c

Overwritten by the product \(QC\), \(Q^{H}C\), \(CQ\), or \(CQ^{H}\) (as specified by side and trans).

scratchpad

Buffer holding scratchpad memory to be used by routine for storing intermediate results.

unmtr (USM Version)#

Syntax

namespace oneapi::mkl::lapack {
  sycl::event unmtr(sycl::queue &queue, oneapi::mkl::side side, oneapi::mkl::uplo upper_lower, oneapi::mkl::transpose trans, std::int64_t m, std::int64_t n, T *a, std::int64_t lda, T *tau, T *c, std::int64_t ldc, T *scratchpad, std::int64_t scratchpad_size, const std::vector<sycl::event> &events = {})
}

Input Parameters

In the descriptions below, r denotes the order of \(Q\):

\(r\)=\(m\)

if side = side::left

\(r\)=\(n\)

if side = side::right

queue

The queue where the routine should be executed.

side

Must be either side::left or side::right.

If side=side::left, \(Q\) or \(Q^{H}\) is applied to \(C\) from the left.

If side=side::right, \(Q\) or \(Q^{H}\) is applied to \(C\) from the right.

upper_lower

Must be either uplo::upper or uplo::lower. Uses the same upper_lower as supplied to hetrd.

trans

Must be either transpose::nontrans or transpose::conjtrans.

If trans=transpose::nontrans, the routine multiplies \(C\) by \(Q\).

If trans=transpose::conjtrans, the routine multiplies \(C\) by \(Q^{H}\).

m

The number of rows in the matrix \(C\) (\(m \ge 0\)).

n

The number of columns the matrix \(C\) (\(n \ge 0\)).

k

The number of elementary reflectors whose product defines the matrix \(Q\) (\(0 \le k \le n\)).

a

The pointer to a as returned by hetrd.

lda

The leading dimension of a \((\max(1,r) \le \text{lda})\).

tau

The pointer to tau as returned by hetrd. The dimension of tau must be at least \(\max(1,r-1)\).

c

The array c contains the matrix \(C\). The second dimension of c must be at least \(\max(1,n)\).

ldc

The leading dimension of c \((\max(1,n) \le \text{ldc})\).

scratchpad_size

Size of scratchpad memory as a number of floating point elements of type T. Size should not be less than the value returned by unmtr_scratchpad_size function.

events

List of events to wait for before starting computation. Defaults to empty list.

Output Parameters

c

Overwritten by the product \(QC\), \(Q^{H}C\), \(CQ\), or \(CQ^{H}\) (as specified by side and trans).

scratchpad

Pointer to scratchpad memory to be used by routine for storing intermediate results.

Return Values

Output event to wait on to ensure computation is complete.

Parent topic: LAPACK Singular Value and Eigenvalue Problem Routines