.. _gebrd: gebrd ===== Reduces a general matrix to bidiagonal form. This routine belongs to the ``oneapi::mkl::lapack`` namespace. .. contents:: :local: :depth: 1 Description *********** The routine reduces a general ``m``-by-``n`` matrix ``A`` to a bidiagonal matrix ``B`` by an orthogonal (unitary) transformation. If ``m≥n``, the reduction is given by .. math:: A = QBP^H = \left( \begin{array}{ccc} B_1 \\ 0 \end{array} \right) P^H = Q_1 B_1 P_H where ``B``\ :sub:`1` is an ``n``-by-``n`` upper diagonal matrix, ``Q`` and ``P`` are orthogonal or, for a complex ``A``, unitary matrices; ``Q``\ :sub:`1` consists of the first ``n`` columns of ``Q``. If ``m < n``, the reduction is given by ``A = Q*B*PH = Q*(B10)*PH = Q1*B1*P1H``, where ``B``\ :sub:`1` is an ``m``-by-``m`` lower diagonal matrix, ``Q`` and ``P`` are orthogonal or, for a complex ``A``, unitary matrices; ``P``\ :sub:`1` consists of the first ``m`` columns of ``P``. The routine does not form the matrices ``Q`` and ``P`` explicitly, but represents them as products of elementary reflectors. Routines are provided to work with the matrices ``Q`` and ``P`` in this representation: If the matrix ``A`` is real, - to compute ``Q`` and ``P`` explicitly, call :ref:`orgbr`. If the matrix ``A`` is complex, - to compute ``Q`` and ``P`` explicitly, call :ref:`ungbr`. API *** Syntax ------ .. code-block:: cpp namespace oneapi::mkl::lapack { void gebrd(cl::sycl::queue &queue, std::int64_t m, std::int64_t n, cl::sycl::buffer &a, std::int64_t lda, cl::sycl::buffer &d, cl::sycl::buffer &e, cl::sycl::buffer &tauq, cl::sycl::buffer &taup, cl::sycl::buffer &scratchpad, std::int64_t scratchpad_size) } ``gebrd`` supports the following precision and devices. .. list-table:: :header-rows: 1 * - T - Devices Supported * - ``float`` - Host, CPU, GPU * - ``double`` - Host, CPU, GPU * - ``std::complex`` - Host, CPU, GPU * - ``std::complex`` - Host, CPU, 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``). a The buffer holding matrix ``A``. The second dimension of ``a`` must be at least ``max(1, m)``. lda The leading dimension of ``a``. scratchpad Buffer holding scratchpad memory to be used by the routine for storing intermediate results. scratchpad_size Size of scratchpad memory as a number offloating point elements of typeT. Size should not be less then the valuereturned by the :ref:`gebrd_scratchpad_size`\ function. Output Parameters ----------------- a If ``m≥n``, the diagonal and first super-diagonal of a are overwritten by the upper bidiagonal matrix ``B``. The elements below the diagonal, with the buffer tauq, represent the orthogonal matrix ``Q`` as a product of elementary reflectors, and the elements above the first superdiagonal, with the buffer taup, represent the orthogonal matrix ``P`` as a product of elementary reflectors. If ``m