A newer version of this document is available. Customers should click here to go to the newest version.
orgqr_batch (Buffer Strided Version)
Generates the real orthogonal matrices Qi of the batch QR factorizations formed by geqrf_batch (Buffer Strided Version). This routine belongs to the oneapi::mkl::lapack namespace.
Description
The routine generates the wholes or parts of the orthogonal matrices Qi of the batch of QR factorizations formed by the routine geqrf_batch (Buffer Strided Version). Usually, Qi is determined from the QR factorization of an m-by-p matrix Ai with m ≥ p.
To compute the whole matrices Qi, use:
orgqr_batch(queue, m, m, p, a, ...)
To compute the leading p columns of Qi (which form an orthonormal basis in the space spanned by the columns of Ai):
orgqr_batch(queue, m, p, p, a, ...)
To compute the matrices Qik of the QR factorizations of leading k columns of the matrices Ai:
orgqr_batch(queue, m, m, k, a, ...)
To compute the leading k columns of Qik (which form an orthonormal basis in the space spanned by leading k columns of the matrices Ai):
orgqr_batch(queue, m, k, k, a, ...)
API
Syntax
namespace oneapi::mkl::lapack {
  void orgqr_batch(sycl::queue &queue,
  std::int64_t m,
  std::int64_t n,
  std::int64_t k,
  sycl::buffer<T> &a,
  std::int64_t lda,
  std::int64_t stride_a,
  sycl::buffer<T> &tau,
  std::int64_t stride_tau,
  std::int64_t batch_size,
  sycl::buffer<T> &scratchpad,
  std::int64_t scratchpad_size)
} 
   This function supports the following precisions and devices:
T  |  
        Devices supported  |  
       
|---|---|
float  |  
        CPU and GPU  |  
       
double  |  
        CPU and GPU  |  
       
Input Parameters
- queue
 -  
     
Device queue where calculations will be performed.
 - m
 -  
     
The number of rows in the matrices Ai (m ≥ 0).
 - n
 -  
     
The number of columns in the matrices Ai (n ≥ 0).
 - k
 -  
     
the number of elementary reflectors whose product defines the matrices Qi (0 ≤ k ≤ n) .
 - a
 -  
     
Array resulting after a call to geqrf_batch (Buffer Strided Version).
 - lda
 -  
     
The leading dimension of Ai (lda ≥ max(1,m)).
 - stride_a
 -  
     
The stride between the beginnings of matrices Ai inside the batch array a (stride_a ≥ max(1, lda * n)).
 - tau
 -  
     
Array resulting after a call to geqrf_batch (Buffer Strided Version).
 - stride_tau
 -  
     
The stride between the beginnings of arrays taui inside the array tau (stride_tau ≥ max(1, min(m,n))).
 - batch_size
 -  
     
The number of problems in a batch (batch_size ≥ 0).
 - 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 than the value returned by orgqr_batch_scratchpad_size (Strided Version).
 
Output Parameters
- a
 -  
     
Array data is overwritten by a batch of n leading columns of the m-by-m orthogonal matrices Qi.
 
Exceptions
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.  |