geqrf_batch (Buffer Strided Version)

Computes the batch of QR factorizations of a general m-by-n matrices. This routine belongs to the oneapi::mkl::lapack namespace.

Description

The routine forms the QiRi factorizations of a general m-by-n matrices Ai. No pivoting is performed.

The routine does not form the matrix Qi explicitly. Instead, Qi is represented as a product of min(m, n) elementary reflectors. Routines are provided to work with Qi in this representation.

API

Syntax

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

Function supports the following precisions and devices.

T

Devices supported

float

Host, CPU, and GPU

double

Host, CPU, and GPU

std::complex<float>

Host, CPU, and GPU

std::complex<double>

Host, CPU, and GPU

Input Parameters

queue

Device queue where calculations will be performed.

m

The number of rows in the matrices Ai (0 m).

n

The number of columns in the matrices Ai (0 n).

a

Array holding input matrices Ai.

lda

The leading dimension of Ai .

stride_a

The stride between the beginnings of matrices Ai inside the batch array a.

stride_tau

The stride between the beginnings of arrays taui inside the array tau.

batch_size

Specifies the number of problems in a batch.

scratchpad

Scratchpad memory to be used by 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 geqrf_batch_scratchpad_size (Strided Version).

Output Parameters

a

Overwritten by the factorization data as follows:

The elements on and above the diagonal of the array contain the min(m,n)-by-n upper trapezoidal matrices Ri (Ri is upper triangular if m≥n); the elements below the diagonal, with the array taui, present the orthogonal matrix Qi as a product of min(m,n) elementary reflectors.

tau

Array to store batch of taui, each of size min(m,n), containing scalars that define elementary reflectors for the matrices Qi in its decomposition in a product of elementary reflectors.

Exceptions

Exception

Description

mkl::lapack::batch_exception

This exception is thrown when problems occur during calculations. You can obtain the info code of the problem using the info() method of the exception object:

If info = -n, the n-th parameter had an illegal value.

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