.. _trtrs: trtrs ===== Solves a system of linear equations with a triangular coefficient matrix, with multiple right-hand sides. This routine belongs to the ``oneapi::mkl::lapack`` namespace. .. contents:: :local: :depth: 1 Description *********** The routine solves for ``X`` the following systems of linear equations with a triangular matrix ``A``, with multiple right-hand sides stored in ``B``: .. list-table:: :header-rows: 0 * - A*X = B - - if ``transa`` =\ ``transpose::nontrans``, * - \ ``AT*X = B``\ - - if ``transa`` =\ ``transpose::trans``, * - A\ :sup:`H`\ ``*X`` = B - - if ``transa`` =\ ``transpose::conjtrans`` (for complex matrices only). API *** Syntax ------ .. code-block:: cpp namespace oneapi::mkl::lapack { void trtrs(cl::sycl::queue &queue, mkl::uplo uplo, mkl::transpose trans, mkl::diag diag, std::int64_t n, std::int64_t nrhs, cl::sycl::buffer &a, std::int64_t lda, cl::sycl::buffer &b, std::int64_t ldb, cl::sycl::buffer &scratchpad, std::int64_t scratchpad_size) } ``trtrs`` 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. uplo Indicates whether ``A`` is upper or lower triangular: If uplo = ``uplo::upper``, then ``A`` is upper triangular. If uplo = ``uplo::lower``, then ``A`` is lower triangular. trans If transa = ``transpose::nontrans``, then ``A``\ \*\ ``X`` = ``B`` is solved for ``X``. If transa = ``transpose::trans``, then ``A``\ :sup:`T`\ \*\ ``X`` = ``B`` is solved for ``X``. If transa = ``transpose::conjtrans``, then ``A``\ :sup:`H`\ \*\ ``X`` = ``B`` is solved for ``X``. diag If diag = ``diag::nonunit``, then ``A`` is not a unit triangular matrix. If diag = ``diag::unit``, then ``A`` is unit triangular: diagonal elements of ``A`` are assumed to be 1 and not referenced in the array a. n The order of ``A``; the number of rows in ``B``; n\ ``≥ 0``. nrhs The number of right-hand sides; nrhs\ ``≥ 0``. a Array containing the matrix ``A``. The second dimension of a must be at least ``max(1,n)``. lda The leading dimension of ``a``; lda\ ``≥ max(1, n)``. b Array containing the matrix ``B`` whose columns are the right-hand sides for the systems of equations. The second dimension of b at least ``max(1,nrhs)``. ldb The leading dimension of b; ldb\ ``≥ max(1, n)``. scratchpad Buffer holding scratchpad memory to be used by the 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 the :ref:`trtrs_scratchpad_size` function. Output Parameters ----------------- b Overwritten by the solution matrix ``X``. Exceptions ---------- .. tabularcolumns:: |\Y{0.3}|\Y{0.7}| .. list-table:: :header-rows: 1 * - Exception - Description * - ``mkl::lapack::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 = -i``, the ``i``-th parameter had an illegal value. If ``info`` is equal to the value passed as scratchpad size, and detail() returns non zero, then the passed scratchpad has an insufficient size, and the required size should not be less than the value returned by the detail() method of the exception object. Return Values ------------- info Buffer containing error information. If ``info`` = 0, the execution is successful. If ``info`` = -``i``, the ``i``-th parameter had an illegal value. Known Limitations ------------------ GPU support for this function does not include error reporting through the ``info`` parameter.