ungqr#

Generates the complex unitary matrix \(Q\) of the QR factorization formed by geqrf.

Description

ungqr supports the following precisions.

T

std::complex<float>

std::complex<double>

The routine generates the whole or part of \(m \times m\) unitary matrix \(Q\) of the QR factorization formed by the routines geqrf.

Usually \(Q\) is determined from the QR factorization of an \(m \times p\) matrix \(A\) with \(m \ge p\). To compute the whole matrix \(Q\), use:

oneapi::mkl::lapack::ungqr(queue, m, m, p, a, lda, tau, scratchpad, scratchpad_size)

To compute the leading \(p\) columns of \(Q\) (which form an orthonormal basis in the space spanned by the columns of \(A\)):

oneapi::mkl::lapack::ungqr(queue, m, p, p, a, lda, tau, scratchpad, scratchpad_size)

To compute the matrix \(Q^{k}\) of the QR factorization of the leading \(k\) columns of the matrix \(A\):

oneapi::mkl::lapack::ungqr(queue, m, m, k, a, lda, tau, scratchpad, scratchpad_size)

To compute the leading \(k\) columns of \(Q^{k}\) (which form an orthonormal basis in the space spanned by the leading \(k\) columns of the matrix \(A\)):

oneapi::mkl::lapack::ungqr(queue, m, k, k, a, lda, tau, scratchpad, scratchpad_size)

ungqr (Buffer Version)#

Syntax

namespace oneapi::mkl::lapack {
  void ungqr(sycl::queue &queue, std::int64_t m, std::int64_t n, std::int64_t k, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<T,1> &tau, sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}

Input Parameters

queue

The queue where the routine should be executed.

m

The number of rows in the matrix \(A\) (\(0 \le m\)).

n

The number of columns in the matrix \(A\) (\(0 \le n\)).

k

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

a

The buffer a as returned by geqrf.

lda

The leading dimension of a (\(\text{lda} \le m\)).

tau

The buffer tau as returned by geqrf.

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 ungqr_scratchpad_size function.

Output Parameters

a

Overwritten by \(n\) leading columns of the \(m \times m\) orthogonal matrix \(Q\).

scratchpad

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

ungqr (USM Version)#

Syntax

namespace oneapi::mkl::lapack {
  sycl::event ungqr(sycl::queue &queue, std::int64_t m, std::int64_t n, std::int64_t k, T *a, std::int64_t lda, T *tau, T *scratchpad, std::int64_t scratchpad_size, const std::vector<sycl::event> &events = {})
}

Input Parameters

queue

The queue where the routine should be executed.

m

The number of rows in the matrix \(A\) (\(0 \le m\)).

n

The number of columns in the matrix \(A\) (\(0 \le n\)).

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 geqrf.

lda

The leading dimension of a (\(\text{lda} \le m\)).

tau

The pointer to tau as returned by geqrf.

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 ungqr_scratchpad_size function.

events

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

Output Parameters

a

Overwritten by \(n\) leading columns of the \(m \times m\) orthogonal matrix \(Q\).

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 Linear Equation Routines