Developer Reference

Contents

getrs_batch (Buffer Strided Version)

Solves a batch of linear equations with an LU-factored square coefficient matrices, with multiple right-hand sides. This routine belongs to the
oneapi::mkl::lapack
namespace.

Description

The routine solves for
X
i
the following systems of linear equations:
  • A
    i
    *
    X
    i
    =
    B
    i
    If
    trans = mkl::transpose::notrans
  • A
    i
    T
    *
    X
    i
    =
    B
    i
    If
    trans = mkl::transpose::trans
  • A
    i
    H
    *
    X
    i
    =
    B
    i
    If
    trans = mkl::transpose::conjtrans
Before calling this routine you must call getrf_batch (Buffer Strided Version) to compute the LU factorization of
A
1
.

API

Syntax
namespace oneapi::mkl::lapack { void getrs_batch(cl::sycl::queue &queue, mkl::transpose trans, std::int64_t n, std::int64_t nrhs, cl::sycl::buffer<T> &a, std::int64_t lda, std::int64_t stride_a, cl::sycl::buffer<std::int64_t> &ipiv, std::int64_t stride_ipiv, cl::sycl::buffer<T> &b, std::int64_t ldb, std::int64_t stride_b, std::int64_t batch_size, cl::sycl::buffer<T> &scratchpad, std::int64_t scratchpad_size) }
Function supports the following precisions and devices.
T
Devices supported
float
Host, CPU, and GPU
double
Host, CPU, and GPU
std::complex<float>
Host, CPU, and GPU
std::complex<double>
Host, CPU, and GPU
Input Parameters
queue
Device queue where calculations will be performed.
trans
Indicates the form of the equations:
If trans = mkl::transpose::nontrans, then
A
i
*
X
i
=
B
i
is solved for
X
i
.
If trans = mkl::transpose::trans, then
A
i
T
*
X
i
=
B
i
is solved for
X
i
.
If trans = mkl::transpose::conjtrans, then
A
i
H
*
X
i
=
B
i
is solved for
X
i
.
n
The order of the matrices
A
i
and the number of rows in matrices
B
i
(
0 ≤ n
).
nrhs
The number of right hand sides
(0≤nrhs)
.
a
Array containing the factorizations of the matrices
A
i
, as returned by getrf_batch (Buffer Strided Version).
lda
The leading dimension of
A
i
.
stride_a
The stride between the beginnings of matrices
B
i
inside the batch array
b
.
ipiv
The ipiv array, as returned by getrf_batch (Buffer Strided Version).
stride_ipiv
The stride between the beginnings of arrays
ipiv
i
inside the array
ipiv
.
b
The array containing the matrices
B
i
whose columns are the right-hand sides for the systems of equations.
ldb
The leading dimensions of
B
i
.
batch_size
Specifies the number of problems in a batch.
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 then the value returned by stride version of getrs_batch_scratchpad_size (Strided Version) function.
Output Parameters
b
Overwritten by the solution matrices
X
i
.
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.
If
info
is zero, then the diagonal element of some of
U
i
is zero, and the solve could not be completed. The indexes of such matrices in the batch can be obtained with the ids() method of the exception object. You can obtain the indexes of the first zero diagonal elements in these
U
i
matrices using the infos() method of the exception object.

Product and Performance Information

1

Performance varies by use, configuration and other factors. Learn more at www.Intel.com/PerformanceIndex.