orgtr

Generates the real orthogonal matrix \(Q\) determined by sytrd.

Description

orgtr supports the following precisions.

T

float

double

The routine explicitly generates the \(n \times n\) orthogonal matrix \(Q\) formed by sytrd when reducing a real symmetric matrix \(A\) to tridiagonal form: \(A = QTQ^T\). Use this routine after a call to sytrd.

orgtr (Buffer Version)

Syntax

namespace oneapi::mkl::lapack {
  void orgtr(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<T,1> &tau, sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}

Input Parameters

queue

The queue where the routine should be executed.

upper_lower

Must be uplo::upper or uplo::lower. Uses the same upper_lower as supplied to sytrd.

n

The order of the matrix \(Q\) \((0 \le n)\).

a

The buffer a as returned by sytrd. The second dimension of a must be at least \(\max(1,n)\).

lda

The leading dimension of a \((n \le \text{lda})\).

tau

The buffer tau as returned by sytrd. The dimension of tau must be at least \(\max(1, n-1)\).

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 orgtr_scratchpad_size function.

Output Parameters

a

Overwritten by the orthogonal matrix \(Q\).

scratchpad

Buffer holding scratchpad memory to be used by routine for storing intermediate results.

orgtr (USM Version)

Syntax

namespace oneapi::mkl::lapack {
  sycl::event orgtr(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, T *a, std::int64_t lda, T *tau, T *scratchpad, std::int64_t scratchpad_size, const std::vector<sycl::event> &events = {})
}

Input Parameters

queue

The queue where the routine should be executed.

upper_lower

Must be uplo::upper or uplo::lower. Uses the same upper_lower as supplied to sytrd.

n

The order of the matrix \(Q\) \((0 \le n)\).

a

The pointer to a as returned by sytrd. The second dimension of a must be at least \(\max(1,n)\).

lda

The leading dimension of a \((n \le \text{lda})\).

tau

The pointer to tau as returned by sytrd. The dimension of tau must be at least \(\max(1, n-1)\).

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 orgtr_scratchpad_size function.

events

List of events to wait for before starting computation. Defaults to empty list.

Output Parameters

a

Overwritten by the orthogonal matrix \(Q\).

scratchpad

Pointer to scratchpad memory to be used by routine for storing intermediate results.

Return Values

Output event to wait on to ensure computation is complete.

Parent topic: LAPACK Singular Value and Eigenvalue Problem Routines