.. _ungqr-usm-version:

ungqr (USM Version)
===================

Generates the complex unitary matrix Q of the QR factorization formed by
the ``geqrf`` (USM Version) function. . This routine belongs to the
``oneapi::mkl::lapack`` namespace.


.. contents::
    :local:
    :depth: 1

Description
***********

The routine generates the whole or part of ``m``-by-``m`` unitary
matrix ``Q`` of the ``QR`` factorization formed by the :ref:`geqrf-usm-version`
function.


Usually ``Q`` is determined from the ``QR`` factorization of an ``m``
by ``p`` matrix ``A`` with ``mβ‰₯p``. To compute the whole matrix
``Q``, use:

.. code-block:: cpp

   ungqr(queue, m, m, p, a, lda, tau, ...)

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

.. code-block:: cpp

   ungqr(queue, m, p, p, a, lda, tau, ...)

To compute the matrix ``Q``\ :sup:`k` of the ``QR`` factorization of
the leading ``k`` columns of the matrix ``A``:

.. code-block:: cpp

   ungqr(queue, m, m, k, a, lda, tau, ...)

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

.. code-block:: cpp

   ungqr(queue, m, k, k, a, lda, tau, ...)


API
***


Syntax
------

.. code-block:: cpp

   namespace oneapi::mkl::lapack {
     cl::sycl::event ungqr(cl::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<cl::sycl::event> &events = {})
   }

``ungqr`` (USM version) supports the following precisions and devices:

.. list-table::
   :header-rows: 1

   * -  T
     -  Devices supported
   * -  ``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 matrix ``A`` (``0≀m``).


n
   The number of columns in the matrix ``A`` (``0≀n``).


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


a
   Pointer to the result of :ref:`geqrf-usm-version`.


lda
   The leading dimension of a (``lda≀m``).


tau
   Pointer to the result of :ref:`geqrf-usm-version`.


scratchpad
   Pointer to 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
   :ref:`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-by-m unitary matrix
   ``Q``.


Exceptions
----------

.. tabularcolumns:: |\Y{0.3}|\Y{0.7}|

.. list-table::
   :header-rows: 1

   * - Exception
     - Description

   * -     ``mkl::lapack::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 = -i``, the ``i``-th parameter had an illegal   value.

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


Return Values
-------------

Output event to wait on to ensure computation is complete.