.. _ungqr_batch-group-version:

ungqr_batch (Group Version)
===========================


Generates the complex unitary matrices Q\ :sub:`i` of the batch of QR
factorizations formed by ``geqrf_batch``. This routine belongs to the
``oneapi::mkl::lapack`` namespace.


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

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

The routine generates the wholes or parts of ``m``-by-``m`` unitary
matrix ``Q``\ :sub:`i` of the ``QR`` factorization formed by the
:ref:`geqrf_batch-group-version`
function.


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

.. code-block:: cpp

   ungqr_batch(queue, m, m, p, a, ...)


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


.. code-block:: cpp

   ungqr_batch(queue, m, p, p, a, ...)

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

.. code-block:: cpp

   ungqr_batch(queue, m, m, k, a, ...)

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

.. code-block:: cpp

   ungqr_batch(queue, m, k, k, a, ...)


API
***


Syntax
------

.. code-block:: cpp

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

Function 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
   Array of ``group_count`` parameters m\ :sub:`g` as previously
   supplied to :ref:`geqrf_batch-group-version`.


n
   Array of ``group_count`` parameters n\ :sub:`g` as previously
   supplied to :ref:`geqrf_batch-group-version`.


k
   Array of ``group_count`` parameters ``k``\ :sub:`g` as previously
   supplied to :ref:`geqrf_batch-group-version`. The
   number of elementary reflectors whose product defines the matrices
   ``Q``\ :sub:`i` (0 ≀ ``k``\ :sub:`g` ≀ ``n``\ :sub:`g`) .


a
   Array resulting after a call to :ref:`geqrf_batch-group-version`.


lda
   Array of leading dimension of ``A``\ :sub:`i` as previously
   supplied to :ref:`geqrf_batch-group-version`.


tau
   Array resulting after a call to :ref:`geqrf_batch-group-version`.


group_count
   Specifies the 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
   :ref:`ungqr_batch_scratchpad_size-group-version`.


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


Output Parameters
-----------------

a
   Matrices pointed to by array ``a`` are overwritten by
   ``n``\ :sub:`g` leading columns of the
   ``m``\ :sub:`g`-by-``m``\ :sub:`g` unitary matrices
   ``Q``\ :sub:`i`, where ``g`` is an index of group of parameters
   corresponding to ``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.