gebrd#
Reduces a general matrix to bidiagonal form.
Description
gebrd
supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
The routine reduces a general \(m \times n\) matrix \(A\) to a bidiagonal matrix \(B\) by an orthogonal (unitary) transformation.
If \(m \ge n\), the reduction is given by \(A=QBP^H=\begin{pmatrix}B_1\\0\end{pmatrix}P^H=Q_1B_1P_H\)
where \(B_{1}\) is an \(n \times n\) upper diagonal matrix, \(Q\) and \(P\) are orthogonal or, for a complex \(A\), unitary matrices; \(Q_{1}\) consists of the first \(n\) columns of \(Q\).
If \(m < n\), the reduction is given by
\(A = QBP^H = Q\begin{pmatrix}B_1\\0\end{pmatrix}P^H = Q_1B_1P_1^H\),
where \(B_{1}\) is an \(m \times m\) lower diagonal matrix, \(Q\) and \(P\) are orthogonal or, for a complex \(A\), unitary matrices; \(P_{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 orgbr.
If the matrix \(A\) is complex,
to compute \(Q\) and \(P\) explicitly, call ungbr
gebrd (Buffer Version)#
Syntax
namespace oneapi::mkl::lapack {
void gebrd(cl::sycl::queue &queue, std::int64_t m, std::int64_t n, cl::sycl::buffer<T,1> &a, std::int64_t lda, cl::sycl::buffer<realT,1> &d, cl::sycl::buffer<realT,1> &e, cl::sycl::buffer<T,1> &tauq, cl::sycl::buffer<T,1> &taup, cl::sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}
Input Parameters
- queue
The queue where the routine should be executed.
- m
The number of rows in the matrix \(A\) (\(0 \le m\)).
- n
The number of columns in the matrix \(A\) (\(0 \le n\)).
- a
The buffer \(a\), size (
lda,*
). The buffera
contains the matrix \(A\). The second dimension ofa
must be at least \(\max(1, m)\).- lda
The leading dimension of \(a\).
- 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 gebrd_scratchpad_size function.
Output Parameters
- a
If \(m \ge 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, size at least \(\max(1, \min(m,n))\). Contains the diagonal elements of \(B\).
- e
Buffer, size at least \(\max(1, \min(m,n) - 1)\). Contains the off-diagonal elements of \(B\).
- tauq
Buffer, 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, size at least \(\max(1, \min(m, n))\). The scalar factors of the elementary reflectors which represent the orthogonal or unitary matrix \(P\).
- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
Throws
This routine shall throw the following exceptions if the associated condition is detected. An implementation may throw additional implementation-specific exception(s) in case of error conditions not covered here.
oneapi::mkl::unsupported_device
oneapi::mkl::lapack::invalid_argument
oneapi::mkl::lapack::computation_error
Exception is thrown in case of problems during calculations. The
info
code of the problem can be obtained by info() method of exception object:If
info=-i
, thei
-th parameter had an illegal value.If
info
equals to value passed as scratchpad size, and detail() returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return by detail() method of exception object.
gebrd (USM Version)#
Syntax
namespace oneapi::mkl::lapack {
cl::sycl::event gebrd(cl::sycl::queue &queue, std::int64_t m, std::int64_t n, T *a, std::int64_t lda, RealT *d, RealT *e, T *tauq, T *taup, T *scratchpad, std::int64_t scratchpad_size, const std::vector<cl::sycl::event> &events = {})
}
Input Parameters
- queue
The queue where the routine should be executed.
- m
The number of rows in the matrix \(A\) (\(0 \le m\)).
- n
The number of columns in the matrix \(A\) (\(0 \le n\)).
- a
Pointer to matrix \(A\). The second dimension of
a
must be at least \(\max(1, m)\).- lda
The leading dimension of
a
.- 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 gebrd_scratchpad_size function.
- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
If \(m \ge n\), the diagonal and first super-diagonal of a are overwritten by the upper bidiagonal matrix \(B\). The elements below the diagonal, with the array tauq, represent the orthogonal matrix \(Q\) as a product of elementary reflectors, and the elements above the first superdiagonal, with the array
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 array tauq, represent the orthogonal matrix \(Q\) as a product of elementary reflectors, and the elements above the diagonal, with the array
taup
, represent the orthogonal matrix \(P\) as a product of elementary reflectors.- d
Pointer to memory of size at least \(\max(1, \min(m,n))\). Contains the diagonal elements of \(B\).
- e
Pointer to memory of size at least \(\max(1, \min(m,n) - 1)\). Contains the off-diagonal elements of \(B\).
- tauq
Pointer to memory 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
Pointer to memory of size at least \(\max(1, \min(m, n))\). The scalar factors of the elementary reflectors which represent the orthogonal or unitary matrix \(P\).
- scratchpad
Pointer to scratchpad memory to be used by routine for storing intermediate results.
Throws
This routine shall throw the following exceptions if the associated condition is detected. An implementation may throw additional implementation-specific exception(s) in case of error conditions not covered here.
oneapi::mkl::unsupported_device
oneapi::mkl::lapack::invalid_argument
oneapi::mkl::lapack::computation_error
Exception is thrown in case of problems during calculations. The
info
code of the problem can be obtained by info() method of exception object:If
info=-i
, thei
-th parameter had an illegal value.If
info
equals to value passed as scratchpad size, and detail() returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return by detail() method of exception object.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: LAPACK Singular Value and Eigenvalue Problem Routines