unmrq¶
Multiplies a complex matrix by the unitary matrix \(Q\) of the RQ factorization formed by gerqf.
Description
unmrq
supports the following precisions.
T
std::complex<float>
std::complex<double>
The routine multiplies a rectangular complex \(m \times n\) matrix \(C\) by \(Q\) or \(Q^H\), where \(Q\) is the complex unitary matrix defined as a product of \(k\) elementary reflectors \(H(i)\) of order \(n\): \(Q = H(1)^HH(2)^H ... H(k)^H\) as returned by the RQ factorization routine gerqf.
Depending on the parameters side
and trans
, the routine can form one of
the matrix products \(QC\), \(Q^HC\), \(CQ\), or \(CQ^H\)
(overwriting the result over \(C\)).
unmrq (Buffer Version)¶
Syntax
namespace oneapi::mkl::lapack {
void unmrq(cl::sycl::queue &queue, oneapi::mkl::side side, oneapi::mkl::transpose trans, std::int64_t m, std::int64_t n, std::int64_t k, cl::sycl::buffer<T,1> &a, std::int64_t lda, cl::sycl::buffer<T,1> &tau, cl::sycl::buffer<T,1> &c, std::int64_t ldc, cl::sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}
Input Parameters
- queue
The queue where the routine should be executed.
- side
If
side = oneapi::mkl::side::left
, \(Q\) or \(Q^{H}\) is applied to \(C\) from the left.If
side = oneapi::mkl::side::right
, \(Q\) or \(Q^{H}\) is applied to \(C\) from the right.- trans
If
trans = oneapi::mkl::transpose::nontrans
, the routine multiplies \(C\) by \(Q\).If
trans = oneapi::mkl::transpose::conjtrans
, the routine multiplies \(C\) by \(Q^{H}\).- m
The number of rows in the matrix \(C\) (\(0 \le m\)).
- n
The number of columns in the matrix \(C\) (\(0 \le n\)).
- k
The number of elementary reflectors whose product defines the matrix \(Q\)
If
side = oneapi::mkl::side::left
, \(0 \le k \le m\)If
side = oneapi::mkl::side::right
, \(0 \le k \le n\)- a
The buffer
a
as returned by gerqf. The second dimension ofa
must be at least \(\max(1,k)\).- lda
The leading dimension of
a
.- tau
The buffer
tau
as returned by gerqf.- c
The buffer
c
contains the matrix \(C\). The second dimension ofc
must be at least \(\max(1,n)\).- ldc
The leading dimension of
c
.- 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 unmrq_scratchpad_size function.
Output Parameters
- c
Overwritten by the product \(QC\), \(Q^{H}C\), \(CQ\), or \(CQ^H\) (as specified by
side
andtrans
).- 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 \(\text{info}=-i\), the \(i\)-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.
unmrq (USM Version)¶
Syntax
namespace oneapi::mkl::lapack {
cl::sycl::event unmrq(cl::sycl::queue &queue, oneapi::mkl::side side, oneapi::mkl::transpose trans, std::int64_t m, std::int64_t n, std::int64_t k, T *a, std::int64_t lda, T *tau, T *c, std::int64_t ldc, 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.
- side
If
side = oneapi::mkl::side::left
, \(Q\) or \(Q^{H}\) is applied to \(C\) from the left.If
side = oneapi::mkl::side::right
, \(Q\) or \(Q^{H}\) is applied to \(C\) from the right.- trans
If
trans = oneapi::mkl::transpose::nontrans
, the routine multiplies \(C\) by \(Q\).If
trans = oneapi::mkl::transpose::conjtrans
, the routine multiplies \(C\) by \(Q^{H}\).- m
The number of rows in the matrix \(C\) (\(0 \le m\)).
- n
The number of columns in the matrix \(C\) (\(0 \le n\)).
- k
The number of elementary reflectors whose product defines the matrix \(Q\)
If
side = oneapi::mkl::side::left
, \(0 \le k \le m\)If
side = oneapi::mkl::side::right
, \(0 \le k \le n\)- a
The pointer to
a
as returned by gerqf. The second dimension ofa
must be at least \(\max(1,k)\).- lda
The leading dimension of
a
.- tau
The pointer to
tau
as returned by gerqf.- c
The pointer
c
points to the matrix \(C\). The second dimension ofc
must be at least \(\max(1,n)\).- ldc
The leading dimension of
c
.- 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 unmrq_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- c
Overwritten by the product \(QC\), \(Q^{H}C\), \(CQ\), or \(CQ^{H}\) (as specified by
side
andtrans
).- 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 \(\text{info}=-i\), the \(i\)-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 Linear Equation Routines