LAPACK++ 2024.05.31
LAPACK C++ API
Loading...
Searching...
No Matches
Householder reflectors and plane rotations

Functions

void lapack::larf (lapack::Side side, int64_t m, int64_t n, double const *v, int64_t incv, double tau, double *C, int64_t ldc)
 
void lapack::larf (lapack::Side side, int64_t m, int64_t n, float const *v, int64_t incv, float tau, float *C, int64_t ldc)
 
void lapack::larf (lapack::Side side, int64_t m, int64_t n, std::complex< double > const *v, int64_t incv, std::complex< double > tau, std::complex< double > *C, int64_t ldc)
 Applies a elementary reflector H to a m-by-n matrix C, from either the left or the right.
 
void lapack::larf (lapack::Side side, int64_t m, int64_t n, std::complex< float > const *v, int64_t incv, std::complex< float > tau, std::complex< float > *C, int64_t ldc)
 
void lapack::larfb (lapack::Side side, lapack::Op trans, lapack::Direction direction, lapack::StoreV storev, int64_t m, int64_t n, int64_t k, double const *V, int64_t ldv, double const *T, int64_t ldt, double *C, int64_t ldc)
 
void lapack::larfb (lapack::Side side, lapack::Op trans, lapack::Direction direction, lapack::StoreV storev, int64_t m, int64_t n, int64_t k, float const *V, int64_t ldv, float const *T, int64_t ldt, float *C, int64_t ldc)
 
void lapack::larfb (lapack::Side side, lapack::Op trans, lapack::Direction direction, lapack::StoreV storev, int64_t m, int64_t n, int64_t k, std::complex< double > const *V, int64_t ldv, std::complex< double > const *T, int64_t ldt, std::complex< double > *C, int64_t ldc)
 Applies a block reflector \(H\) or its transpose \(H^H\) to a m-by-n matrix C, from either the left or the right.
 
void lapack::larfb (lapack::Side side, lapack::Op trans, lapack::Direction direction, lapack::StoreV storev, int64_t m, int64_t n, int64_t k, std::complex< float > const *V, int64_t ldv, std::complex< float > const *T, int64_t ldt, std::complex< float > *C, int64_t ldc)
 
void lapack::larfg (int64_t n, double *alpha, double *X, int64_t incx, double *tau)
 
void lapack::larfg (int64_t n, float *alpha, float *X, int64_t incx, float *tau)
 
void lapack::larfg (int64_t n, std::complex< double > *alpha, std::complex< double > *X, int64_t incx, std::complex< double > *tau)
 Generates an elementary reflector H of order n, such that:
 
void lapack::larfg (int64_t n, std::complex< float > *alpha, std::complex< float > *X, int64_t incx, std::complex< float > *tau)
 
void lapack::larfgp (int64_t n, double *alpha, double *X, int64_t incx, double *tau)
 
void lapack::larfgp (int64_t n, float *alpha, float *X, int64_t incx, float *tau)
 
void lapack::larfgp (int64_t n, std::complex< double > *alpha, std::complex< double > *X, int64_t incx, std::complex< double > *tau)
 Generates an elementary reflector H of order n, such that:
 
void lapack::larfgp (int64_t n, std::complex< float > *alpha, std::complex< float > *X, int64_t incx, std::complex< float > *tau)
 
void lapack::larft (lapack::Direction direction, lapack::StoreV storev, int64_t n, int64_t k, double const *V, int64_t ldv, double const *tau, double *T, int64_t ldt)
 
void lapack::larft (lapack::Direction direction, lapack::StoreV storev, int64_t n, int64_t k, float const *V, int64_t ldv, float const *tau, float *T, int64_t ldt)
 
void lapack::larft (lapack::Direction direction, lapack::StoreV storev, int64_t n, int64_t k, std::complex< double > const *V, int64_t ldv, std::complex< double > const *tau, std::complex< double > *T, int64_t ldt)
 Forms the triangular factor T of a complex block reflector H of order n, which is defined as a product of k elementary reflectors.
 
void lapack::larft (lapack::Direction direction, lapack::StoreV storev, int64_t n, int64_t k, std::complex< float > const *V, int64_t ldv, std::complex< float > const *tau, std::complex< float > *T, int64_t ldt)
 
void lapack::larfx (lapack::Side side, int64_t m, int64_t n, double const *v, double tau, double *C, int64_t ldc)
 
void lapack::larfx (lapack::Side side, int64_t m, int64_t n, float const *v, float tau, float *C, int64_t ldc)
 
void lapack::larfx (lapack::Side side, int64_t m, int64_t n, std::complex< double > const *v, std::complex< double > tau, std::complex< double > *C, int64_t ldc)
 Applies a elementary reflector H to a m-by-n matrix C, from either the left or the right.
 
void lapack::larfx (lapack::Side side, int64_t m, int64_t n, std::complex< float > const *v, std::complex< float > tau, std::complex< float > *C, int64_t ldc)
 
void lapack::larfy (lapack::Uplo uplo, int64_t n, double const *V, int64_t incv, double tau, double *C, int64_t ldc)
 
void lapack::larfy (lapack::Uplo uplo, int64_t n, float const *V, int64_t incv, float tau, float *C, int64_t ldc)
 
void lapack::larfy (lapack::Uplo uplo, int64_t n, std::complex< double > const *V, int64_t incv, std::complex< double > tau, std::complex< double > *C, int64_t ldc)
 Applies an elementary reflector, or Householder matrix, H, to an n x n Hermitian matrix C, from both the left and the right.
 
void lapack::larfy (lapack::Uplo uplo, int64_t n, std::complex< float > const *V, int64_t incv, std::complex< float > tau, std::complex< float > *C, int64_t ldc)
 

Detailed Description

Function Documentation

◆ larf()

void lapack::larf ( lapack::Side  side,
int64_t  m,
int64_t  n,
std::complex< double > const *  v,
int64_t  incv,
std::complex< double >  tau,
std::complex< double > *  C,
int64_t  ldc 
)

Applies a elementary reflector H to a m-by-n matrix C, from either the left or the right.

H is represented in the form

\[ H = I - \tau v v^H, \]

where \(\tau\) is a scalar and v is a vector.

If \(\tau = 0,\) then H is taken to be the unit matrix.

To apply \(H^H,\) supply \(\text{conj}(\tau)\) instead of \(\tau.\)

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]side
  • lapack::Side::Left: form \(H C\)
  • lapack::Side::Right: form \(C H\)
[in]mThe number of rows of the matrix C.
[in]nThe number of columns of the matrix C.
[in]vThe vector v in the representation of H. v is not used if tau = 0.
  • If side = Left, the vector v of length 1 + (m-1)*abs(incv);
  • if side = Right, the vector v of length 1 + (n-1)*abs(incv).
[in]incvThe increment between elements of v. incv != 0.
[in]tauThe value tau in the representation of H.
[in,out]CThe m-by-n matrix C, stored in an ldc-by-n array. On entry, the m-by-n matrix C. On exit, C is overwritten by the matrix \(H C\) if side = Left, or \(C H\) if side = Right.
[in]ldcThe leading dimension of the array C. ldc >= max(1,m).

◆ larfb()

void lapack::larfb ( lapack::Side  side,
lapack::Op  trans,
lapack::Direction  direction,
lapack::StoreV  storev,
int64_t  m,
int64_t  n,
int64_t  k,
std::complex< double > const *  V,
int64_t  ldv,
std::complex< double > const *  T,
int64_t  ldt,
std::complex< double > *  C,
int64_t  ldc 
)

Applies a block reflector \(H\) or its transpose \(H^H\) to a m-by-n matrix C, from either the left or the right.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]side
  • lapack::Side::Left: apply \(H\) or \(H^H\) from the Left
  • lapack::Side::Right: apply \(H\) or \(H^H\) from the Right
[in]trans
  • lapack::Op::NoTrans: apply \(H \) (No transpose)
  • lapack::Op::ConjTrans: apply \(H^H\) (Conjugate transpose)
[in]directionIndicates how H is formed from a product of elementary reflectors
  • lapack::Direction::Forward: \(H = H(1) H(2) \dots H(k)\)
  • lapack::Direction::Backward: \(H = H(k) \dots H(2) H(1)\)
[in]storevIndicates how the vectors which define the elementary reflectors are stored:
  • lapack::StoreV::Columnwise
  • lapack::StoreV::Rowwise
[in]mThe number of rows of the matrix C.
[in]nThe number of columns of the matrix C.
[in]kThe order of the matrix T (= the number of elementary reflectors whose product defines the block reflector).
  • If side = Left, m >= k >= 0;
  • if side = Right, n >= k >= 0.
[in]V
  • If storev = Columnwise:
    • if side = Left, the m-by-k matrix V, stored in an ldv-by-k array;
    • if side = Right, the n-by-k matrix V, stored in an ldv-by-k array.
  • If storev = Rowwise:
    • if side = Left, the k-by-m matrix V, stored in an ldv-by-m array;
    • if side = Right, the k-by-n matrix V, stored in an ldv-by-n array.
  • See Further Details.
[in]ldvThe leading dimension of the array V.
  • If storev = Columnwise and side = Left, ldv >= max(1,m);
  • if storev = Columnwise and side = Right, ldv >= max(1,n);
  • if storev = Rowwise, ldv >= k.
[in]TThe k-by-k matrix T, stored in an ldt-by-k array. The triangular k-by-k matrix T in the representation of the block reflector.
[in]ldtThe leading dimension of the array T. ldt >= k.
[in,out]CThe m-by-n matrix C, stored in an ldc-by-n array. On entry, the m-by-n matrix C. On exit, C is overwritten by \(H C\) or \(H^H C\) or \(C H\) or \(C H^H\).
[in]ldcThe leading dimension of the array C. ldc >= max(1,m).
Further Details

The shape of the matrix V and the storage of the vectors which define the H(i) is best illustrated by the following example with n = 5 and k = 3. The elements equal to 1 are not stored. The rest of the array is not used.

direction = Forward and          direction = Forward and
storev = Columnwise:             storev = Rowwise:

V = (  1       )                 V = (  1 v1 v1 v1 v1 )
    ( v1  1    )                     (     1 v2 v2 v2 )
    ( v1 v2  1 )                     (        1 v3 v3 )
    ( v1 v2 v3 )
    ( v1 v2 v3 )

direction = Backward and         direction = Backward and
storev = Columnwise:             storev = Rowwise:

V = ( v1 v2 v3 )                 V = ( v1 v1  1       )
    ( v1 v2 v3 )                     ( v2 v2 v2  1    )
    (  1 v2 v3 )                     ( v3 v3 v3 v3  1 )
    (     1 v3 )
    (        1 )

◆ larfg()

void lapack::larfg ( int64_t  n,
std::complex< double > *  alpha,
std::complex< double > *  X,
int64_t  incx,
std::complex< double > *  tau 
)

Generates an elementary reflector H of order n, such that:

\[ H^H \begin{bmatrix} \alpha \\ x \end{bmatrix} = \begin{bmatrix} \beta \\ 0 \end{bmatrix}; \quad H^H H = I. \]

where \(\alpha\) and \(\beta\) are scalars, with \(\beta\) real, and x is an (n-1)-element vector. H is represented in the form

\[ H = I - \tau \begin{bmatrix} 1 \\ v \end{bmatrix} \begin{bmatrix} 1 & v^H \end{bmatrix}, \]

where \(\tau\) is a scalar and v is a (n-1)-element vector. For complex H, note that H is not hermitian.

If the elements of x are all zero and alpha is real, then \(\tau = 0\) and H is taken to be the unit matrix.

Otherwise \(1 \le \text{real}(\tau) \le 2\) and \(|\tau - 1| \le 1.\)

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]nThe order of the elementary reflector.
[in,out]alphaOn entry, the value alpha. On exit, it is overwritten with the value beta.
[in,out]XThe vector X of length 1+(n-2)*abs(incx). On entry, the vector x. On exit, it is overwritten with the vector v.
[in]incxThe increment between elements of X. incx > 0.
[out]tauThe value tau.

◆ larfgp()

void lapack::larfgp ( int64_t  n,
std::complex< double > *  alpha,
std::complex< double > *  X,
int64_t  incx,
std::complex< double > *  tau 
)

Generates an elementary reflector H of order n, such that:

\[ H^H \begin{bmatrix} \alpha \\ x \end{bmatrix} = \begin{bmatrix} \beta \\ 0 \end{bmatrix}; \quad H^H H = I. \]

where \(\alpha\) and \(\beta\) are scalars, with \(\beta\) real and non-negative, and x is an (n-1)-element vector. H is represented in the form

\[ H = I - \tau \begin{bmatrix} 1 \\ v \end{bmatrix} \begin{bmatrix} 1 & v^H \end{bmatrix}, \]

where \(\tau\) is a scalar and v is a (n-1)-element vector. For complex H, note that H is not hermitian.

If the elements of x are all zero and alpha is real, then \(\tau = 0\) and H is taken to be the unit matrix.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Since
LAPACK 3.2.2
Parameters
[in]nThe order of the elementary reflector.
[in,out]alphaOn entry, the value alpha. On exit, it is overwritten with the value beta.
[in,out]XThe vector X of length 1+(n-2)*abs(incx). On entry, the vector x. On exit, it is overwritten with the vector v.
[in]incxThe increment between elements of X. incx > 0.
[out]tauThe value tau.

◆ larft()

void lapack::larft ( lapack::Direction  direction,
lapack::StoreV  storev,
int64_t  n,
int64_t  k,
std::complex< double > const *  V,
int64_t  ldv,
std::complex< double > const *  tau,
std::complex< double > *  T,
int64_t  ldt 
)

Forms the triangular factor T of a complex block reflector H of order n, which is defined as a product of k elementary reflectors.

If direction = Forward, \(H = H(1) H(2) \dots H(k)\) and T is upper triangular;

If direction = Backward, \(H = H(k) \dots H(2) H(1)\) and T is lower triangular.

If storev = Columnwise, the vector which defines the elementary reflector H(i) is stored in the i-th column of the array V, and

\[ H = I - V T V^H. \]

If storev = Rowwise, the vector which defines the elementary reflector H(i) is stored in the i-th row of the array V, and

\[ H = I - V^H T V. \]

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]directionSpecifies the order in which the elementary reflectors are multiplied to form the block reflector:
  • lapack::Direction::Forward: \(H = H(1) H(2) \dots H(k)\)
  • lapack::Direction::Backward: \(H = H(k) \dots H(2) H(1)\)
[in]storevSpecifies how the vectors which define the elementary reflectors are stored (see also Further Details):
  • lapack::StoreV::Columnwise
  • lapack::StoreV::Rowwise
[in]nThe order of the block reflector H. n >= 0.
[in]kThe order of the triangular factor T (= the number of elementary reflectors). k >= 1.
[in]V
  • If storev = Columnwise, the n-by-k matrix V, stored in an ldv-by-k array;
  • if storev = Rowwise, the k-by-n matrix V, stored in an ldv-by-n array.
    See further details.
[in]ldvThe leading dimension of the array V.
  • If storev = Columnwise, ldv >= max(1,n);
  • if storev = Rowwise, ldv >= k.
[in]tauThe vector tau of length k. tau(i) must contain the scalar factor of the elementary reflector H(i).
[out]TThe k-by-k matrix T, stored in an ldt-by-k array. The k-by-k triangular factor T of the block reflector.
  • If direction = Forward, T is upper triangular;
  • if direction = Backward, T is lower triangular.
    The rest of the array is not used.
[in]ldtThe leading dimension of the array T. ldt >= k.
Further Details

The shape of the matrix V and the storage of the vectors which define the H(i) is best illustrated by the following example with n = 5 and k = 3. The elements equal to 1 are not stored.

direction = Forward and          direction = Forward and
storev = Columnwise:             storev = Rowwise:

V = (  1       )                 V = (  1 v1 v1 v1 v1 )
    ( v1  1    )                     (     1 v2 v2 v2 )
    ( v1 v2  1 )                     (        1 v3 v3 )
    ( v1 v2 v3 )
    ( v1 v2 v3 )

direction = Backward and         direction = Backward and
storev = Columnwise:             storev = Rowwise:

V = ( v1 v2 v3 )                 V = ( v1 v1  1       )
    ( v1 v2 v3 )                     ( v2 v2 v2  1    )
    (  1 v2 v3 )                     ( v3 v3 v3 v3  1 )
    (     1 v3 )
    (        1 )

◆ larfx()

void lapack::larfx ( lapack::Side  side,
int64_t  m,
int64_t  n,
std::complex< double > const *  v,
std::complex< double >  tau,
std::complex< double > *  C,
int64_t  ldc 
)

Applies a elementary reflector H to a m-by-n matrix C, from either the left or the right.

H is represented in the form

\[ H = I - \tau v v^H \]

where \(\tau\) is a scalar and v is a vector.

If \(\tau = 0,\) then H is taken to be the unit matrix

This version uses inline code if H has order < 11.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]side
  • lapack::Side::Left: form \(H C\)
  • lapack::Side::Right: form \(C H\)
[in]mThe number of rows of the matrix C.
[in]nThe number of columns of the matrix C.
[in]v
  • If side = Left, the vector v of length m;
  • if side = right, the vector v of length n.
[in]tauThe value tau in the representation of H.
[in,out]CThe m-by-n matrix C, stored in an ldc-by-n array. On entry, the m-by-n matrix C. On exit, C is overwritten by the matrix \(H C\) if side = Left, or \(C H\) if side = Right.
[in]ldcThe leading dimension of the array C. ldc >= max(1,m).

◆ larfy()

void lapack::larfy ( lapack::Uplo  uplo,
int64_t  n,
std::complex< double > const *  V,
int64_t  incv,
std::complex< double >  tau,
std::complex< double > *  C,
int64_t  ldc 
)

Applies an elementary reflector, or Householder matrix, H, to an n x n Hermitian matrix C, from both the left and the right.

H is represented in the form

\[ H = I - \tau v v^H \]

where \(\tau\) is a scalar and \(v\) is a vector.

If \(tau\) is zero, then \(H\) is taken to be the unit matrix.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Since
LAPACK 3.7.0
Parameters
[in]uploWhether the upper or lower triangular part of the Hermitian matrix C is stored.
  • lapack::Uplo::Upper: Upper triangle
  • lapack::Uplo::Lower: Lower triangle
[in]nThe number of rows and columns of the matrix C. n >= 0.
[in]VThe vector V of length 1 + (n-1)*abs(incv).
[in]incvThe increment between successive elements of v. incv must not be zero.
[in]tauThe value tau as described above.
[in,out]CThe n-by-n matrix C, stored in an ldc-by-n array. On entry, the matrix C. On exit, C is overwritten by \(H C H^H\).
[in]ldcThe leading dimension of the array C. ldc >= max( 1, n ).