.. _orgqr_batch-usm-strided-version: orgqr_batch (USM Strided Version) ================================= Generates the orthogonal matrices ``Q``\ :sub:`i` of the batch QR factorizations formed by ``geqrf_batch`` (USM Strided Version). This routine belongs to the ``oneapi::mkl::lapack`` namespace. .. contents:: :local: :depth: 1 Description *********** The routine generates the wholes or parts of the orthogonal matrices ``Q``\ :sub:`1` of the batch of QR factorizations formed by the routine :ref:`geqrf_batch-usm-strided-version`. Usually, ``Q``\ :sub:`1` is determined from the QR factorization of an ``m``-by-``p`` matrix ``A``\ :sub:`1` with ``m`` ≥\ ``p``. To compute the whole matrices ``Q``\ :sub:`1`, use: .. code-block:: cpp orgqr_batch(queue, m, m, p, a, ...) To compute the leading ``p`` columns of ``Q``\ :sub:`1` (which form an orthonormal basis in the space spanned by the columns of ``A``\ :sub:`1`): .. code-block:: cpp orgqr_batch(queue, m, p, p, a, ...) To compute the matrices ``Q``\ :sub:`1`\ :sup:`k` of the QR factorizations of leading ``k`` columns of the matrices ``A``\ :sub:`1`: .. code-block:: cpp orgqr_batch(queue, m, m, k, a, ...) To compute the leading ``k`` columns of ``Q``\ :sub:`1`\ :sup:`k` (which form an orthonormal basis in the space spanned by leading ``k`` columns of the matrices ``A``\ :sub:`1`): .. code-block:: cpp orgqr_batch(queue, m, k, k, a, ...) API *** Syntax ------ .. code-block:: cpp namespace oneapi::mkl::lapack { cl::sycl::event orgqr_batch(cl::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 &events = {}) } Function supports the following precisions and devices. .. list-table:: :header-rows: 1 * - T - Devices supported * - ``float`` - Host, CPU, and GPU * - ``double`` - Host, CPU, and GPU Input Parameters ---------------- queue Device queue where calculations will be performed. m The number of rows in the matrices ``A``\ :sub:`i` (``0 ≤ m``). n The number of columns in the matrices ``A``\ :sub:`i` (``0 ≤ n``). k the number of elementary reflectors whose product defines the matrices ``Q``\ :sub:`i` (0 ≤ ``k`` ≤ ``n``) . a Array resulting after a call to :ref:`geqrf_batch-usm-strided-version`. lda The leading dimension of ``A``\ :sub:`i` (``lda`` ≤ ``m``). stride_a The stride between the beginnings of matrices ``A``\ :sub:`i` inside the batch array ``a``. tau Array resulting after a call to :ref:`geqrf_batch-usm-strided-version`. stride_tau The stride between the beginnings of arrays ``tau``\ :sub:`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 than the value returned by :ref:`orgqr_batch_scratchpad_size-strided-version`. events List of events to wait for before starting computation. Defaults to empty list. Output Parameters ----------------- a Array data is overwritten by a batch of ``n`` leading columns of the ``m``-by-``m`` orthogonal matrices ``Q``\ :sub:`i`. Exceptions ---------- .. tabularcolumns:: |\Y{0.3}|\Y{0.7}| .. list-table:: :header-rows: 1 * - 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 then value returned by the detail() method of the exception object. Return Values ------------- Output event to wait on to ensure computation is complete.