Visible to Intel only — GUID: GUID-4EB9CEC6-04E1-4905-AF4C-009267A2AF88
Visible to Intel only — GUID: GUID-4EB9CEC6-04E1-4905-AF4C-009267A2AF88
gebrd
Reduces a general matrix to bidiagonal form. This routine belongs to the oneapi::mkl::lapack namespace.
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
where B1 is an n-by-n upper diagonal matrix, Q and P are orthogonal or, for a complex A, unitary matrices; Q1 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 B1 is an m-by-m lower diagonal matrix, Q and P are orthogonal or, for a complex A, unitary matrices; P1 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 orgbr.
If the matrix A is complex,
to compute Q and P explicitly, call ungbr.
API
Syntax
namespace oneapi::mkl::lapack { void gebrd(sycl::queue &queue, std::int64_t m, std::int64_t n, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<realT,1> &d, sycl::buffer<realT,1> &e, sycl::buffer<T,1> &tauq, sycl::buffer<T,1> &taup, sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size) }
gebrd supports the following precision and devices.
T |
Devices Supported |
---|---|
float |
CPU, GPU |
double |
CPU, GPU |
std::complex<float> |
CPU, GPU |
std::complex<double> |
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 gebrd_scratchpad_sizefunction.
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<n, the diagonal and first sub-diagonal of a are overwritten by the lower bidiagonal matrix B. The elements below the first subdiagonal, with the buffer tauq, represent the orthogonal matrix Q as a product of elementary reflectors, and the elements above the diagonal, with the buffer taup, represent the orthogonal matrix P as a product of elementary reflectors.
- d
-
Buffer holding array of size at least max(1, min(m,n)). Contains the diagonal elements of B.
- e
-
Buffer holding array of size at least max(1, min(m,n) - 1). Contains the off-diagonal elements of B.
- tauq
-
Buffer holding array of size at least max(1, min(m, n)). The scalar factors of the elementary reflectors which represent the orthogonal or unitary matrix Q.
- taup
-
Buffer holding array of size at least max(1, min(m, n)). The scalar factors of the elementary reflectors which represent the orthogonal or unitary matrix P.
Exceptions
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()i method of the exception object. |