.. _getrs_batch-buffer-strided-version: 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. .. contents:: :local: :depth: 1 Description *********** The routine solves for ``X``\ :sub:`i` the following systems of linear equations: - ``A``\ :sub:`i` \* ``X``\ :sub:`i` = ``B``\ :sub:`i` If ``trans = mkl::transpose::notrans`` - ``A``\ :sub:`i`\ :sup:`T` \* ``X``\ :sub:`i` = ``B``\ :sub:`i` If ``trans = mkl::transpose::trans`` - ``A``\ :sub:`i`\ :sup:`H` \* ``X``\ :sub:`i` = ``B``\ :sub:`i` If ``trans = mkl::transpose::conjtrans`` Before calling this routine you must call :ref:`getrf_batch-buffer-strided-version` to compute the LU factorization of ``A``\ :sub:`1`. API *** Syntax ------ .. code-block:: cpp 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 &a, std::int64_t lda, std::int64_t stride_a, cl::sycl::buffer &ipiv, std::int64_t stride_ipiv, cl::sycl::buffer &b, std::int64_t ldb, std::int64_t stride_b, std::int64_t batch_size, cl::sycl::buffer &scratchpad, std::int64_t scratchpad_size) } Function supports the following precisions and devices. .. list-table:: :header-rows: 1 * - T - Devices supported * - ``float`` - Host, CPU, and GPU * - ``double`` - Host, CPU, and GPU * - ``std::complex`` - Host, CPU, and GPU * - ``std::complex`` - 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``\ :sub:`i`\ \*\ ``X``\ :sub:`i` = ``B``\ :sub:`i` is solved for ``X``\ :sub:`i`. If trans = mkl::transpose::trans, then ``A``\ :sub:`i`\ :sup:`T`\ \*\ ``X``\ :sub:`i` = ``B``\ :sub:`i` is solved for ``X``\ :sub:`i`. If trans = mkl::transpose::conjtrans, then ``A``\ :sub:`i`\ :sup:`H`\ \*\ ``X``\ :sub:`i` = ``B``\ :sub:`i` is solved for ``X``\ :sub:`i`. n The order of the matrices ``A``\ :sub:`i` and the number of rows in matrices ``B``\ :sub:`i` (``0 ≤ n``). nrhs The number of right hand sides ``(0≤nrhs)``. a Array containing the factorizations of the matrices ``A``\ :sub:`i`, as returned by :ref:`getrf_batch-buffer-strided-version`. lda The leading dimension of ``A``\ :sub:`i`. stride_a The stride between the beginnings of matrices ``B``\ :sub:`i` inside the batch array ``b``. ipiv The ipiv array, as returned by :ref:`getrf_batch-buffer-strided-version`. stride_ipiv The stride between the beginnings of arrays ``ipiv``\ :sub:`i` inside the array ``ipiv``. b The array containing the matrices ``B``\ :sub:`i` whose columns are the right-hand sides for the systems of equations. ldb The leading dimensions of ``B``\ :sub:`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 :ref:`getrs_batch_scratchpad_size-strided-version` function. Output Parameters ----------------- b Overwritten by the solution matrices ``X``\ :sub:`i`. Exceptions ---------- .. tabularcolumns:: |\Y{0.3}|\Y{0.7}| .. list-table:: :header-rows: 1 * - 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``\ :sub:`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``\ :sub:`i` matrices using the infos() method of the exception object.