|
TRIQS/nda 1.3.0
Multi-dimensional array library for C++
|
Interface to parts of the LAPACK library.
Classes | |
| class | nda::lapack::gelss_worker< T > |
| Worker class for solving linear least squares problems. More... | |
| class | nda::lapack::gelss_worker_hermitian |
| Specialized worker class for solving linear least squares problems while enforcing a certain hermitian symmetry. More... | |
Functions | |
| template<MemoryMatrix A, MemoryArray B, MemoryVector S> requires (have_same_value_type_v<A, B> and mem::have_host_compatible_addr_space<A, B, S> and is_blas_lapack_v<get_value_t<A>>) | |
| int | nda::lapack::gelss (A &&a, B &&b, S &&s, double rcond, int &rank) |
| Interface to the LAPACK gelss routine. | |
| template<MemoryMatrix A, MemoryVector JPVT, MemoryVector TAU> requires (mem::have_host_compatible_addr_space<A, JPVT, TAU> and have_same_value_type_v<A, TAU> and is_blas_lapack_v<get_value_t<A>>) | |
| int | nda::lapack::geqp3 (A &&a, JPVT &&jpvt, TAU &&tau) |
| Interface to the LAPACK geqp3 routine. | |
| template<MemoryMatrix A, MemoryVector S, MemoryMatrix U, MemoryMatrix VT> requires (have_same_value_type_v<A, U, VT> and mem::have_compatible_addr_space<A, S, U, VT> and is_blas_lapack_v<get_value_t<A>> and std::same_as<double, get_value_t<S>>) | |
| int | nda::lapack::gesvd (A &&a, S &&s, U &&u, VT &&vt) |
| Interface to the LAPACK gesvd routine. | |
| template<MemoryMatrix A, MemoryVector IPIV> requires (mem::have_compatible_addr_space<A, IPIV> and is_blas_lapack_v<get_value_t<A>> and std::is_same_v<get_value_t<IPIV>, int>) | |
| int | nda::lapack::getrf (A &&a, IPIV &&ipiv) |
| Interface to the LAPACK getrf routine. | |
| template<MemoryMatrix A, MemoryVector IPIV> requires (mem::have_host_compatible_addr_space<A, IPIV> and is_blas_lapack_v<get_value_t<A>> and std::is_same_v<get_value_t<IPIV>, int>) | |
| int | nda::lapack::getri (A &&a, IPIV const &ipiv) |
| Interface to the LAPACK getri routine. | |
| template<Matrix A, MemoryArray B, MemoryVector IPIV> requires ((MemoryMatrix<A> or is_conj_array_expr<A>) and have_same_value_type_v<A, B> and mem::have_compatible_addr_space<A, B, IPIV> and is_blas_lapack_v<get_value_t<A>>) | |
| int | nda::lapack::getrs (A const &a, B &&b, IPIV const &ipiv) |
| Interface to the LAPACK getrs routine. | |
| template<MemoryVector DL, MemoryVector D, MemoryVector DU, MemoryArray B> requires (have_same_value_type_v<DL, D, DU, B> and mem::have_host_compatible_addr_space<DL, D, DU, B> and is_blas_lapack_v<get_value_t<DL>>) | |
| int | nda::lapack::gtsv (DL &&dl, D &&d, DU &&du, B &&b) |
| Interface to the LAPACK gtsv routine. | |
| template<MemoryMatrix A, MemoryVector W> requires (mem::have_host_compatible_addr_space<A> and std::same_as<std::complex<double>, get_value_t<A>> and std::same_as<double, get_value_t<W>>) | |
| int | nda::lapack::heev (A &&a, W &&w, char jobz='V') |
| Interface to the LAPACK heev routine. | |
| template<MemoryMatrix A, MemoryMatrix B, MemoryVector W> requires (mem::have_host_compatible_addr_space<A, B> and std::same_as<std::complex<double>, get_value_t<A>> and have_same_value_type_v<A, B> and std::same_as<double, get_value_t<W>>) | |
| int | nda::lapack::hegv (A &&a, B &&b, W &&w, char jobz='V', int itype=1) |
| Interface to the LAPACK hegv routine. | |
| template<MemoryMatrix A, MemoryVector TAU> requires (mem::have_host_compatible_addr_space<A> and std::is_same_v<double, get_value_t<A>> and have_same_value_type_v<A, TAU>) | |
| int | nda::lapack::orgqr (A &&a, TAU &&tau) |
| Interface to the LAPACK orgqr routine. | |
| template<MemoryMatrix A, MemoryVector W> requires (mem::have_host_compatible_addr_space<A> and std::same_as<double, get_value_t<A>> and have_same_value_type_v<A, W>) | |
| int | nda::lapack::syev (A &&a, W &&w, char jobz='V') |
| Interface to the LAPACK syev routine. | |
| template<MemoryMatrix A, MemoryMatrix B, MemoryVector W> requires (mem::have_host_compatible_addr_space<A, B> and std::same_as<double, get_value_t<A>> and have_same_value_type_v<A, B, W>) | |
| int | nda::lapack::sygv (A &&a, B &&b, W &&w, char jobz='V', int itype=1) |
| Interface to the LAPACK sygv routine. | |
| template<MemoryMatrix A, MemoryVector TAU> requires (mem::have_host_compatible_addr_space<A> and std::is_same_v<std::complex<double>, get_value_t<A>> and have_same_value_type_v<A, TAU>) | |
| int | nda::lapack::ungqr (A &&a, TAU &&tau) |
| Interface to the LAPACK ungqr routine. | |
| int nda::lapack::gelss | ( | A && | a, |
| B && | b, | ||
| S && | s, | ||
| double | rcond, | ||
| int & | rank ) |
#include <nda/lapack/gelss.hpp>
Interface to the LAPACK gelss routine.
Computes the minimum norm solution to a complex linear least squares problem:
\[ \min_{\mathbf{x}} | \mathbf{b} - \mathbf{A x} |_2 \]
using the singular value decomposition (SVD) of \( \mathbf{A} \). \( \mathbf{A} \) is an \( m \times n \) matrix which may be rank-deficient.
Several right hand side vectors \( \mathbf{b} \) and solution vectors \( \mathbf{x} \) can be handled in a single call; they are stored as the columns of the \( m \times n_{\mathrm{rhs}} \) right hand side matrix \(\mathbf{B} \) and the \( n \times n_{\mathrm{rhs}} \) solution matrix \( \mathbf{X} \).
The effective rank of \( \mathbf{A} \) is determined by treating as zero those singular values which are less than \( r_{\mathrm{cond}} \) times the largest singular value.
| A | nda::MemoryMatrix type. |
| B | nda::MemoryArray type. |
| S | nda::MemoryVector type. |
| a | Input/output matrix. On entry, the \( m \times n \) matrix \( \mathbf{A} \). On exit, the first \(\min(m,n) \) rows of \( \mathbf{A} \) are overwritten with its right singular vectors, stored rowwise. |
| b | Input/output array. On entry, the \( m \times n_{\mathrm{rhs}} \) right hand side matrix \( \mathbf{B} \). On exit, \( \mathbf{B} \) is overwritten by the \( n \times n_{\mathrm{rhs}} \) solution matrix \(\mathbf{X} \). If \( m \geq n \) and if the effective rank is equal \( n \), the residual sum-of-squares for the solution in the ith column is given by the sum of squares of the modulus of elements \( n + 1 \) to \( m \) in that column. |
| s | Output vector. The singular values of \( \mathbf{A} \) in decreasing order. The condition number of \(\mathbf{A} \) in the 2-norm is \( s_1 / s_{min(m,n)} \). |
| rcond | It is used to determine the effective rank of \( \mathbf{A} \). Singular values \( s_i \leq r_{\mathrm{cond}} s_1 \) are treated as zero. If \( r_{\mathrm{cond}} < 0 \), machine precision is used instead. |
| rank | Output variable. The effective rank of \( \mathbf{A} \), i.e. the number of singular values which are greater than \( r_{\mathrm{cond}} s_1 \). |
| int nda::lapack::geqp3 | ( | A && | a, |
| JPVT && | jpvt, | ||
| TAU && | tau ) |
#include <nda/lapack/geqp3.hpp>
Interface to the LAPACK geqp3 routine.
Computes a QR factorization with column pivoting of a matrix \( \mathbf{A} \):
\[ \mathbf{A P} = \mathbf{Q R} \]
using Level 3 BLAS.
The matrix \( \mathbf{Q} \) is represented as a product of elementary reflectors
\[ \mathbf{Q} = \mathbf{H}(1) \mathbf{H}(2) \ldots \mathbf{H}(k) \; , \]
where \( k = \min(m,n) \).
Each \( \mathbf{H}(i) \) has the form
\[ \mathbf{H}(i) = \mathbf{I} - \tau_i * \mathbf{v}_i \mathbf{v}_i^H \]
where \( \tau_i \) is a real/complex scalar, and \( \mathbf{v}_i \) is a real/complex vector with
| A | nda::MemoryMatrix type. |
| JPVT | nda::MemoryVector type. |
| TAU | nda::MemoryVector type. |
| a | Input/output matrix. On entry, the \( m \times n \) matrix \( \mathbf{A} \). On exit, the upper triangle of the array contains the \( \min(m,n) \times n \) upper trapezoidal matrix \( \mathbf{R} \); the elements below the diagonal, together with the array \( \mathbf{\tau} \), represent the unitary matrix \(\mathbf{Q} \) as a product of \( \min(m,n) \) elementary reflectors. |
| jpvt | Input/output vector. On entry, if the jth element is \( \neq 0 \), the jth column of \( \mathbf{A} \) is permuted to the front of \( \mathbf{A P} \) (a leading column); if the jth element is equal 0, the jth column of \( \mathbf{A} \) is a free column. On exit, if the jth element is equal \( k \), then the jth column of \( \mathbf{A P} \) was the the kth column of \( \mathbf{A} \). |
| tau | Output vector. The scalar factors \( \tau_i \) of the elementary reflectors \( \mathbf{H}(i) \). |
| int nda::lapack::gesvd | ( | A && | a, |
| S && | s, | ||
| U && | u, | ||
| VT && | vt ) |
#include <nda/lapack/gesvd.hpp>
Interface to the LAPACK gesvd routine.
Computes the singular value decomposition (SVD) of an \( m \times n \) matrix \( \mathbf{A} \). The SVD is written as
\[ \mathbf{A} = \mathbf{U} \mathbf{S} \mathbf{V}^H \; , \]
where \( \mathbf{S} \) is an \( m \times n \) matrix which is zero except for its \( \min(m,n) \) diagonal elements, \( \mathbf{U} \) is an \( m \times m \) unitary matrix, and \( \mathbf{V} \) is an \( n \times n \) unitary matrix. The diagonal elements of \( \mathbf{S} \) are the singular values of \( \mathbf{A} \); they are real and non-negative, and are returned in descending order. The first \( min(m,n) \) columns of \(\mathbf{U} \) and \( \mathbf{V} \) are the left and right singular vectors of \( \mathbf{A} \).
Note that the routine returns \( \mathbf{V}^H \), not \( \mathbf{V} \).
| A | nda::MemoryMatrix type. |
| S | nda::MemoryVector type. |
| U | nda::MemoryMatrix type. |
| VT | nda::MemoryMatrix type. |
| a | Input/output matrix. On entry, the \( m \times n \) matrix \( \mathbf{A} \). On exit, the contents of \( \mathbf{A} \) are destroyed. |
| s | Output vector. The singular values of \( \mathbf{A} \), sorted so that \( s_i \geq s_{i+1} \). |
| u | Output matrix. It contains the \( m \times m \) unitary matrix \( \mathbf{U} \). |
| vt | Output matrix. It contains the \( n \times n \) unitary matrix \( \mathbf{V}^H \). |
| int nda::lapack::getrf | ( | A && | a, |
| IPIV && | ipiv ) |
#include <nda/lapack/getrf.hpp>
Interface to the LAPACK getrf routine.
Computes an LU factorization of a general \( m \times n \) matrix \( \mathbf{A} \) using partial pivoting with row interchanges.
The factorization has the form
\[ \mathbf{A} = \mathbf{P L U} \]
where \( \mathbf{P} \) is a permutation matrix, \( \mathbf{L} \) is lower triangular with unit diagonal elements (lower trapezoidal if \( m > n \)), and \( \mathbf{U} \) is upper triangular (upper trapezoidal if \(m < n \)).
This is the right-looking Level 3 BLAS version of the algorithm.
| A | nda::MemoryMatrix type. |
| IPIV | nda::MemoryVector type. |
| a | Input/output matrix. On entry, the \( m \times n \) matrix to be factored. On exit, the factors \(\mathbf{L} \) and \( \mathbf{U} \) from the factorization \( \mathbf{A} = \mathbf{P L U} \); the unit diagonal elements of \( \mathbf{L} \) are not stored. |
| ipiv | Output vector. The pivot indices, i.e. for \( 1 \leq i \leq \min(m,n) \), row \( i \) of the matrix was interchanged with row ipiv(i-1). |
| int nda::lapack::getri | ( | A && | a, |
| IPIV const & | ipiv ) |
#include <nda/lapack/getri.hpp>
Interface to the LAPACK getri routine.
Computes the inverse of an \( n \times n \) matrix \( \mathbf{A} \) matrix using the LU factorization computed by nda::lapack::getrf.
This method inverts \( \mathbf{U} \) and then computes \( \mathrm{inv}(\mathbf{A}) \) by solving the system \( \mathrm{inv}(\mathbf{A}) L = \mathrm{inv}(\mathbf{U}) \) for \( \mathrm{inv}(\mathbf{A}) \).
| A | nda::MemoryMatrix type. |
| IPIV | nda::MemoryVector type. |
| a | Input/output matrix. On entry, the factors \( \mathbf{L} \) and \( \mathbf{U} \) from the factorization \( \mathbf{A} = \mathbf{P L U} \) as computed by nda::lapack::getrf. On exit, if INFO == 0, the inverse of the original matrix \( \mathbf{A} \). |
| ipiv | Input vector. The pivot indices from nda::lapack::getrf, i.e. for \( 1 \leq i \leq n \), row i of the matrix was interchanged with row ipiv(i). |
| int nda::lapack::getrs | ( | A const & | a, |
| B && | b, | ||
| IPIV const & | ipiv ) |
#include <nda/lapack/getrs.hpp>
Interface to the LAPACK getrs routine.
Solves a system of linear equations
with a general \( n \times n \) matrix \( \mathbf{A} \) using the LU factorization computed by nda::lapack::getrf.
| A | nda::Matrix type. |
| B | nda::MemoryArray type. |
| IPIV | nda::MemoryVector type. |
| a | Input matrix. The factors \( \mathbf{L} \) and \( \mathbf{U} \) from the factorization \( \mathbf{A} = \mathbf{P L U} \) as computed by nda::lapack::getrf. |
| b | Input/output matrix. On entry, the right hand side matrix \( \mathbf{B} \). On exit, the solution matrix \( \mathbf{X} \). |
| ipiv | Input vector. The pivot indices from nda::lapack::getrf, i.e. for \( 1 \leq i \leq n \), row i of the matrix was interchanged with row ipiv(i). |
| int nda::lapack::gtsv | ( | DL && | dl, |
| D && | d, | ||
| DU && | du, | ||
| B && | b ) |
#include <nda/lapack/gtsv.hpp>
Interface to the LAPACK gtsv routine.
Solves the equation
\[ \mathbf{A} \mathbf{X} = \mathbf{B}, \]
where \( \mathbf{A} \) is an \( n \times n \) tridiagonal matrix, by Gaussian elimination with partial pivoting.
Note that the equation \( \mathbf{A}^H \mathbf{X} = \mathbf{B} \) may be solved by interchanging the order of the arguments containing the subdiagonal elements.
| DL | nda::MemoryVector type. |
| D | nda::MemoryVector type. |
| DU | nda::MemoryVector type. |
| B | nda::MemoryArray type. |
| dl | Input/Output vector. On entry, it must contain the \( n - 1 \) subdiagonal elements of \( \mathbf{A} \). On exit, it is overwritten by the \( n - 2 \) elements of the second superdiagonal of the upper triangular matrix \( \mathbf{U} \) from the LU factorization of \( \mathbf{A} \). |
| d | Input/Output vector. On entry, it must contain the diagonal elements of \( \mathbf{A} \). On exit, it is overwritten by the \( n \) diagonal elements of \( \mathbf{U} \). |
| du | Input/Output vector. On entry, it must contain the \( n - 1 \) superdiagonal elements of \( \mathbf{A} \). On exit, it is overwritten by the \( n - 1 \) elements of the first superdiagonal of \( \mathbf{U} \) . |
| b | Input/Output array. On entry, the \( n \times n_{\mathrm{rhs}} \) right hand side matrix \( \mathbf{B} \). On exit, if INFO == 0, the \( n \times n_{\mathrm{rhs}} \) solution matrix \( \mathbf{X} \). |
| int nda::lapack::heev | ( | A && | a, |
| W && | w, | ||
| char | jobz = 'V' ) |
#include <nda/lapack/heev.hpp>
Interface to the LAPACK heev routine.
Computes all eigenvalues \( \lambda_i \) and, optionally, eigenvectors \( \mathbf{v}_i \) of a complex hermitian matrix eigenvalue problem of the form
\[ \mathbf{A} \mathbf{v}_i = \lambda_i \mathbf{v}_i \; , \]
for a given complex hermitian matrix \( \mathbf{A} \).
| A | nda::MemoryMatrix with complex value type. |
| a | Input/output matrix. On entry, the hermitian matrix \( \mathbf{A} \). On exit, if jobz = V, \(\mathbf{A} \) contains the orthonormal eigenvectors of the matrix \( \mathbf{A} \). If jobz = N, then on exit \( \mathbf{A} \) is destroyed. |
| w | Output vector. The eigenvalues in ascending order. |
| jobz | Character indicating whether to compute eigenvectors and eigenvalues ('V') or eigenvalues only ('N'). |
| int nda::lapack::hegv | ( | A && | a, |
| B && | b, | ||
| W && | w, | ||
| char | jobz = 'V', | ||
| int | itype = 1 ) |
#include <nda/lapack/hegv.hpp>
Interface to the LAPACK hegv routine.
Computes all eigenvalues \( \lambda_i \) and, optionally, eigenvectors \( \mathbf{v}_i \) of a complex generalized Hermitian-definite eigenvalue problem of the form
Here \( \mathbf{A} \) and \( \mathbf{B} \) are assumed to be Hermitian and \( \mathbf{B} \) is also positive definite.
| A | nda::MemoryMatrix with complex value type. |
| B | nda::MemoryMatrix with complex value type. |
| a | Input/output matrix. On entry, the Hermitian matrix \( \mathbf{A} \). On exit, if jobz = V, \(\mathbf{A} \) contains the matrix \( \mathbf{V} \) of normalized eigenvectors such that \( \mathbf{V}^H \mathbf{B} \mathbf{V} = \mathbf{I} \) (if itype = 1 or itype = 2) or \( \mathbf{V}^H \mathbf{B}^{-1} \mathbf{V} = \mathbf{I} \) (if itype = 3). If jobz = N, then on exit \( \mathbf{A} \) is destroyed. |
| b | Input/output matrix. On entry, the symmetric positive definite matrix \( \mathbf{B} \). On exit, the part of \( \mathbf{B} \) containing the matrix is overwritten by the triangular factor \( \mathbf{U} \) or \( \mathbf{L} \) from a Cholesky factorization. |
| w | Output vector. The eigenvalues in ascending order. |
| jobz | Character indicating whether to compute eigenvectors and eigenvalues ('V') or eigenvalues only ('N'). |
| itype | Specifies the problem to be solved. |
| int nda::lapack::orgqr | ( | A && | a, |
| TAU && | tau ) |
#include <nda/lapack/orgqr.hpp>
Interface to the LAPACK orgqr routine.
Generates an \( m \times n \) real matrix \( \mathbf{Q} \) with orthonormal columns, which is defined as the first \( n \) columns of a product of \( k \) elementary reflectors of order \( m \)
\[ \mathbf{Q} = \mathbf{H}(1) \mathbf{H}(2) \ldots \mathbf{H}(k) \; , \]
as returned by nda::lapack::geqp3.
Each \( \mathbf{H}(i) \) has the form
\[ \mathbf{H}(i) = \mathbf{I} - \tau_i * \mathbf{v}_i \mathbf{v}_i^T \]
where \( \tau_i \) is a real scalar, and \( \mathbf{v}_i \) is a real vector with
| A | nda::MemoryMatrix with double value type. |
| TAU | nda::MemoryVector with double value type. |
| a | Input/output matrix. On entry, the ith column must contain the vector which defines the elementary reflector \( H(i) \; , i = 1,2,...,k \), as returned by nda::lapack::geqp3 in the first \( k \) columns. On exit, the \( m \times n \) matrix \( \mathbf{Q} \). |
| tau | Input vector. \( \tau_i \) must contain the scalar factor of the elementary reflector \(\mathbf{H}(i) \), as returned by nda::lapack::geqp3. |
| int nda::lapack::syev | ( | A && | a, |
| W && | w, | ||
| char | jobz = 'V' ) |
#include <nda/lapack/syev.hpp>
Interface to the LAPACK syev routine.
Computes all eigenvalues \( \lambda_i \) and, optionally, eigenvectors \( \mathbf{v}_i \) of a real symmetric eigenvalue problem of the form
\[ \mathbf{A} \mathbf{v}_i = \lambda_i \mathbf{v}_i \; , \]
for a given real symmetric matrix \( \mathbf{A} \).
| A | nda::MemoryMatrix with double value type. |
| a | Input/output matrix. On entry, the symmetric matrix \( \mathbf{A} \). On exit, if jobz = V, \(\mathbf{A} \) contains the orthonormal eigenvectors of the matrix \( \mathbf{A} \). If jobz = N, then on exit \( \mathbf{A} \) is destroyed. |
| w | Output vector. The eigenvalues in ascending order. |
| jobz | Character indicating whether to compute eigenvectors and eigenvalues ('V') or eigenvalues only ('N'). |
| int nda::lapack::sygv | ( | A && | a, |
| B && | b, | ||
| W && | w, | ||
| char | jobz = 'V', | ||
| int | itype = 1 ) |
#include <nda/lapack/sygv.hpp>
Interface to the LAPACK sygv routine.
Computes all eigenvalues \( \lambda_i \) and, optionally, eigenvectors \( \mathbf{v}_i \) of a real generalized symmetric-definite eigenvalue problem of the form
Here \( \mathbf{A} \) and \( \mathbf{B} \) are assumed to be symmetric and \( \mathbf{B} \) is also positive definite.
| A | nda::MemoryMatrix with double value type. |
| B | nda::MemoryMatrix with double value type. |
| a | Input/output matrix. On entry, the symmetric matrix \( \mathbf{A} \). On exit, if jobz = V, \(\mathbf{A} \) contains the matrix \( \mathbf{V} \) of normalized eigenvectors such that \( \mathbf{V}^T \mathbf{B} \mathbf{V} = \mathbf{I} \) (if itype = 1 or itype = 2) or \( \mathbf{V}^T \mathbf{B}^{-1} \mathbf{V} = \mathbf{I} \) (if itype = 3). If jobz = N, then on exit \( \mathbf{A} \) is destroyed. |
| b | Input/output matrix. On entry, the symmetric positive definite matrix \( \mathbf{B} \). On exit, the part of \( \mathbf{B} \) containing the matrix is overwritten by the triangular factor \( \mathbf{U} \) or \( \mathbf{L} \) from a Cholesky factorization. |
| w | Output vector. The eigenvalues in ascending order. |
| jobz | Character indicating whether to compute eigenvectors and eigenvalues ('V') or eigenvalues only ('N'). |
| itype | Specifies the problem to be solved. |
| int nda::lapack::ungqr | ( | A && | a, |
| TAU && | tau ) |
#include <nda/lapack/ungqr.hpp>
Interface to the LAPACK ungqr routine.
Generates an \( m \times n \) complex matrix \( \mathbf{Q} \) with orthonormal columns, which is defined as the first \( n \) columns of a product of \( k \) elementary reflectors of order \( m \)
\[ \mathbf{Q} = \mathbf{H}(1) \mathbf{H}(2) \ldots \mathbf{H}(k) \; , \]
as returned by nda::lapack::geqp3.
Each \( \mathbf{H}(i) \) has the form
\[ \mathbf{H}(i) = \mathbf{I} - \tau_i * \mathbf{v}_i \mathbf{v}_i^H \]
where \( \tau_i \) is a complex scalar, and \( \mathbf{v}_i \) is a complex vector with
| A | nda::MemoryMatrix with complex value type. |
| TAU | nda::MemoryVector with complex value type. |
| a | Input/output matrix. On entry, the ith column must contain the vector which defines the elementary reflector \( H(i) \; , i = 1,2,...,k \), as returned by nda::lapack::geqp3 in the first \( k \) columns. On exit, the \( m \times n \) matrix \( \mathbf{Q} \). |
| tau | Input vector. \( \tau_i \) must contain the scalar factor of the elementary reflector \(\mathbf{H}(i) \), as returned by nda::lapack::geqp3. |