orgqr_batch

Generates the orthogonal/unitary matrix \(Q_i\) of the QR factorizations for a group of general matrices.

Description

orgqr_batch supports the following precisions.

T

float

double

orgqr_batch (Buffer Version)

Description

The buffer version of orgqr_batch supports only the strided API.

Strided API

The routine generates the wholes or parts of \(m \times n\) orthogonal matrices \(Q_i\) of the batch of QR factorizations formed by the Strided API of the geqrf_batch (Buffer Version) function.
Usually \(Q_i\) is determined from the QR factorization of an \(m \times p\) matrix \(A_i\) with \(m \ge p\).
To compute the whole matrices \(Q_i\), use:
orgqr_batch(queue, m, m, p, a, ...)
To compute the leading \(p\) columns of \(Q_i\) (which form an orthonormal basis in the space spanned by the columns of \(A_i\)):
orgqr_batch(queue, m, p, p, a, ...)
To compute the matrices \(Q_i^k\) of the QR factorizations of leading \(k\) columns of the matrices \(A_i\):
orgqr_batch(queue, m, m, k, a, ...)
To compute the leading \(k\) columns of \(Q_i^k\) (which form an orthonormal basis in the space spanned by leading \(k\) columns of the matrices \(A_i\)):
orgqr_batch(queue, m, k, k, a, ...)

Syntax

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

Input Parameters

queue

Device queue where calculations will be performed.

m

Number of rows in the matrices \(A_i\) (\(0 \le m\)).

n

Number of columns in the matrices \(A_i\) (\(0 \le n\)).

k

Number of elementary reflectors whose product defines the matrices \(Q_i\) (\(0 \le k \le n\)).

a

Array resulting after call to the Strided API of the geqrf_batch (Buffer Version) function.

lda

Leading dimension of \(A_i\) (\(\text{lda} \le m\)).

stride_a

The stride between the beginnings of matrices \(A_i\) inside the batch array a.

tau

Array resulting from call to the Strided API of the geqrf_batch (Buffer Version) function.

stride_tau

Stride between the beginnings of arrays \(\tau_i\) 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 then the value returned by the Strided API of the orgqr_batch_scratchpad_size function.

Output Parameters

a

Batch of \(n\) leading columns of the \(m \times m\) orthogonal matrices \(Q_i\).

orgqr_batch (USM Version)

Description

The USM version of orgqr_batch supports the group API and strided API.

Group API

The routine generates the wholes or parts of \(m \times n\) orthogonal matrices \(Q_i\) of the batch of QR factorizations formed by the Group API of the onemkl_lapack_geqrf_batch_usm function.
Usually \(Q_i\) is determined from the QR factorization of an \(m \times p\) matrix \(A_i\) with \(m \ge p\).
To compute the whole matrices \(Q_i\), use:
orgqr_batch(queue, m, m, p, a, ...)
To compute the leading \(p\) columns of \(Q_i\) (which form an orthonormal basis in the space spanned by the columns of \(A_i\)):
orgqr_batch(queue, m, p, p, a, ...)
To compute the matrices \(Q_i^k\) of the QR factorizations of leading \(k\) columns of the matrices \(A_i\):
orgqr_batch(queue, m, m, k, a, ...)
To compute the leading \(k\) columns of \(Q_i^k\) (which form an orthonormal basis in the space spanned by leading \(k\) columns of the matrices \(A_i\)):
orgqr_batch(queue, m, k, k, a, ...)

Syntax

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

Input Parameters

queue

Device queue where calculations will be performed.

m

Array of group_count \(m_g\) parameters as previously supplied to group version of geqrf_batch function.

n

Array of group_count \(n_g\) parameters as previously supplied to group version of geqrf_batch function.

k

Array of group_count \(k_g\) parameters as previously supplied to the Group API of the onemkl_lapack_geqrf_batch_usm function. The number of elementary reflectors whose product defines the matrices \(Q_i\) (\(0 \le k_g \le n_g\)).

a

Array resulting after call to the Group API of the onemkl_lapack_geqrf_batch_usm function.

lda

Array of leading dimensions of \(A_i\) as previously supplied to the Group API of the onemkl_lapack_geqrf_batch_usm function.

tau

Array resulting after call to the Group API of the onemkl_lapack_geqrf_batch_usm function.

group_count

Number of groups of parameters. Must be at least 0.

group_sizes

Array of group_count integers. Array element with index \(g\) specifies the number of problems to solve for each of the groups of parameters \(g\). So the total number of problems to solve, batch_size, is a sum of all parameter group sizes.

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 then the value returned by Group API of the orgqr_batch_scratchpad_size function.

events

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

Output Parameters

a

\(n_g\) leading columns of the \(m_g \times m_g\) orthogonal matrices \(Q_i\), where \(g\) is an index of group of parameters corresponding to \(Q_i\).

Return Values

Output event to wait on to ensure computation is complete.

Strided API

The routine generates the wholes or parts of \(m \times n\) orthogonal matrices \(Q_i\) of the batch of QR factorizations formed by the Strided API of the onemkl_lapack_geqrf_batch_usm function.
Usually \(Q_i\) is determined from the QR factorization of an \(m \times p\) matrix \(A_i\) with \(m \ge p\).
To compute the whole matrices \(Q_i\), use:
orgqr_batch(queue, m, m, p, a, ...)
To compute the leading \(p\) columns of \(Q_i\) (which form an orthonormal basis in the space spanned by the columns of \(A_i\)):
orgqr_batch(queue, m, p, p, a, ...)
To compute the matrices \(Q_i^k\) of the QR factorizations of leading \(k\) columns of the matrices \(A_i\):
orgqr_batch(queue, m, m, k, a, ...)
To compute the leading \(k\) columns of \(Q_i^k\) (which form an orthonormal basis in the space spanned by leading \(k\) columns of the matrices \(A_i\)):
orgqr_batch(queue, m, k, k, a, ...)

Syntax

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

Input Parameters

queue

Device queue where calculations will be performed.

m

Number of rows in the matrices \(A_i\) (\(0 \le m\)).

n

Number of columns in the matrices \(A_i\) (\(0 \le n\)).

k

Number of elementary reflectors whose product defines the matrices \(Q_i\) (\(0 \le k \le n\)).

a

Array resulting after call to the Strided API of the onemkl_lapack_geqrf_batch_usm function.

lda

Leading dimension of \(A_i\) (\(\text{lda} \le m\)).

stride_a

The stride between the beginnings of matrices \(A_i\) inside the batch array a.

tau

Array resulting from call to the Strided API of the onemkl_lapack_geqrf_batch_usm function.

stride_tau

Stride between the beginnings of arrays \(\tau_i\) 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 then the value returned by the Strided API of the orgqr_batch_scratchpad_size function.

events

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

Output Parameters

a

Batch of \(n\) leading columns of the \(m \times m\) orthogonal matrices \(Q_i\).

Return Values

Output event to wait on to ensure computation is complete.

Parent topic: LAPACK-like Extensions Routines