gerqf#

Computes the RQ factorization of a general \(m \times n\) matrix.

Description

gerqf supports the following precisions.

T

float

double

std::complex<float>

std::complex<double>

The routine forms the RQ factorization of a general \(m \times n\) matrix \(A\). No pivoting is performed. The routine does not form the matrix \(Q\) explicitly. Instead, \(Q\) is represented as a product of \(\min(m, n)\) elementary reflectors. Routines are provided to work with \(Q\) in this representation

gerqf (Buffer Version)#

Syntax

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

Input Parameters

queue

Device queue where calculations will be performed.

m

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

n

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

a

Buffer holding input matrix \(A\). The second dimension of a must be at least \(\max(1, n)\).

lda

The leading dimension of a, at least \(\max(1, m)\).

scratchpad

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

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 the gerqf_scratchpad_size function.

Output Parameters

a

Output buffer, overwritten by the factorization data as follows:

If \(m \le n\), the upper triangle of the subarray a(1:m, n-m+1:n) contains the \(m \times m\) upper triangular matrix \(R\); if \(m \ge n\), the elements on and above the \((m-n)\)-th subdiagonal contain the \(m \times n\) upper trapezoidal matrix \(R\)

In both cases, the remaining elements, with the array tau, represent the orthogonal/unitary matrix \(Q\) as a product of \(\min(m,n)\) elementary reflectors.

tau

Array, size at least \(\min(m,n)\).

Contains scalars that define elementary reflectors for the matrix \(Q\) in its decomposition in a product of elementary reflectors.

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::mkl::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

oneapi::mkl::unsupported_device

oneapi::mkl::lapack::invalid_argument

oneapi::mkl::lapack::computation_error

Exception is thrown in case of problems during calculations. The info code of the problem can be obtained by info() method of exception object:

If info = -i, the \(i\)-th parameter had an illegal value.

If info equals to value passed as scratchpad size, and detail() returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return by detail() method of exception object.

gerqf (USM Version)#

Syntax

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

Input Parameters

queue

Device queue where calculations will be performed.

m

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

n

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

a

Buffer holding input matrix \(A\). The second dimension of a must be at least \(\max(1, n)\).

lda

The leading dimension of a, at least \(\max(1, m)\).

scratchpad

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

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 the gerqf_scratchpad_size function.

events

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

Output Parameters

a

Output buffer, overwritten by the factorization data as follows:

If \(m \le n\), the upper triangle of the subarray a(1:m, n-m+1:n) contains the \(m \times m\) upper triangular matrix \(R\); if \(m \ge n\), the elements on and above the \((m-n)\)-th subdiagonal contain the \(m \times n\) upper trapezoidal matrix \(R\)

In both cases, the remaining elements, with the array tau, represent the orthogonal/unitary matrix \(Q\) as a product of \(\min(m,n)\) elementary reflectors.

tau

Array, size at least \(\min(m,n)\).

Contains scalars that define elementary reflectors for the matrix \(Q\) in its decomposition in a product of elementary reflectors.

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::mkl::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

oneapi::mkl::unsupported_device

oneapi::mkl::lapack::invalid_argument

oneapi::mkl::lapack::computation_error

Exception is thrown in case of problems during calculations. The info code of the problem can be obtained by info() method of exception object:

If info = -i, the \(i\)-th parameter had an illegal value.

If info equals to value passed as scratchpad size, and detail() returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return by detail() method of exception object.

Return Values

Output event to wait on to ensure computation is complete.

Parent topic: LAPACK Linear Equation Routines