NcmLapack

NcmLapack — Encapsulated LAPACK functions.

Functions

Object Hierarchy

    GBoxed
    ╰── NcmLapackWS

Description

This object is dedicated to encapsulate functions from LAPACK choosing the most suitable backend.

Priority order: (1) LAPACK and (2) GSL. It no longer tries to use clapack or lapacke, it is faster and simpler to stick to fortran's lapack.

The description of each function follows its respective LAPACK documentation.

Functions

ncm_lapack_ws_new ()

NcmLapackWS *
ncm_lapack_ws_new (void);

Creates a new Lapack workspace object.

Returns

a newly created NcmLapackWS.

[transfer full]


ncm_lapack_ws_dup ()

NcmLapackWS *
ncm_lapack_ws_dup (NcmLapackWS *ws);

Duplicates a Lapack workspace object.

Parameters

ws

a NcmLapackWS

 

Returns

a copy of ws .

[transfer full]


ncm_lapack_ws_free ()

void
ncm_lapack_ws_free (NcmLapackWS *ws);

Frees a Lapack workspace object.

Parameters

ws

a NcmLapackWS

 

ncm_lapack_ws_clear ()

void
ncm_lapack_ws_clear (NcmLapackWS **ws);

Clears a Lapack workspace object.

Parameters

ws

a NcmLapackWS

 

ncm_lapack_dptsv ()

gint
ncm_lapack_dptsv (gdouble *d,
                  gdouble *e,
                  gdouble *b,
                  gdouble *x,
                  gint n);

This function computes the solution to a real system of linear equations $A*X = B$ (B = b ), where $A$ is an N-by-N (N = n ) symmetric positive definite tridiagonal matrix, and $X$ and $B$ are N-by-NRHS (NRHS = 1) matrices.

$A$ is factored as $A = L*D*L^T$, and the factored form of $A$ is then used to solve the system of equations.

Parameters

d

array of doubles with dimension n

 

e

array of doubles with dimension n -1

 

b

array of doubles with dimension n

 

x

array of doubles with dimension n

 

n

The order of the matrix $A$ (>= 0)

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: i, the leading minor of order i is not positive definite, and the solution has not been computed. The factorization has not been completed unless i = N.


ncm_lapack_dpotrf ()

gint
ncm_lapack_dpotrf (gchar uplo,
                   gint n,
                   gdouble *a,
                   gint lda);

This function computes the Cholesky factorization of a real symmetric positive definite matrix a .

The factorization has the form $A = U^T * U$, if uplo = 'U', or $A = L * L^T$, if uplo = 'L', where A = a , $U$ is an upper triangular matrix and $L$ is lower triangular.

Parameters

uplo

'U' upper triangle of a is stored; 'L' lower triangle of a is stored

 

n

The order of the matrix a , n >= 0

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1,n )

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: i, the leading minor of order i is not positive definite, and the factorization could not be completed.


ncm_lapack_dpotri ()

gint
ncm_lapack_dpotri (gchar uplo,
                   gint n,
                   gdouble *a,
                   gint lda);

This function computes the inverse of a real symmetric positive definite matrix a = A using the Cholesky factorization $A = U^T*U$ or $A = L*L^T$ computed by ncm_lapack_dpotrf().

Parameters

uplo

'U' upper triangle of a is stored; 'L' lower triangle of a is stored

 

n

The order of the matrix a , n >= 0

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1,n )

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dpotrs ()

gint
ncm_lapack_dpotrs (gchar uplo,
                   gint n,
                   gint nrhs,
                   gdouble *a,
                   gint lda,
                   gdouble *b,
                   gint ldb);

This function computes the solution of $A X = B$ for a real symmetric positive definite matrix a = A using the Cholesky factorization $A = U^T*U$ or $A = L*L^T$ already performed by ncm_lapack_dpotrf(). On entry b contain the vectors $B$ and on exit b contain the solutions if the return is 0.

Parameters

uplo

'U' upper triangle of a is stored; 'L' lower triangle of a is stored

 

n

The order of the matrix a , n >= 0

 

nrhs

Number of right-hand-side vectors to solve

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1, n )

 

b

array of doubles with dimension (n , ldb )

 

ldb

The leading dimension of the array b , ldb >= max (1, n )

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dposv ()

gint
ncm_lapack_dposv (gchar uplo,
                  gint n,
                  gint nrhs,
                  gdouble *a,
                  gint lda,
                  gdouble *b,
                  gint ldb);

This function computes the solution of $A X = B$ for a real symmetric positive definite matrix a = A using the Cholesky factorization $A = U^T*U$ or $A = L*L^T$. On entry b contain the vectors $B$ and on exit b contain the solutions if the return is 0.

Parameters

uplo

'U' upper triangle of a is stored; 'L' lower triangle of a is stored

 

n

The order of the matrix a , n >= 0

 

nrhs

Number of right-hand-side vectors to solve

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1, n )

 

b

array of doubles with dimension (n , ldb )

 

ldb

The leading dimension of the array b , ldb >= max (1, n )

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dsytrf ()

gint
ncm_lapack_dsytrf (gchar uplo,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gint *ipiv,
                   NcmLapackWS *ws);

This function computes the factorization of a real symmetric matrix a , using the Bunch-Kaufman diagonal pivoting method.

Parameters

uplo

'U' upper triangle of a is stored; 'L' lower triangle of a is stored

 

n

The order of the matrix a , n >= 0

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1, n )

 

ipiv

Information about decomposition swaps and blocks

 

ws

a NcmLapackWS

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dsytrs ()

gint
ncm_lapack_dsytrs (gchar uplo,
                   gint n,
                   gint nrhs,
                   gdouble *a,
                   gint lda,
                   gint *ipiv,
                   gdouble *b,
                   gint ldb);

This function computes the solution of $A X = B$ for a real symmetric positive definite matrix a = A using the Cholesky factorization $A = U^T*U$ or $A = L*L^T$ already performed by ncm_lapack_dpotrf(). On entry b contain the vectors $B$ and on exit b contain the solutions if the return is 0.

Parameters

uplo

'U' upper triangle of a is stored; 'L' lower triangle of a is stored

 

n

The order of the matrix a , n >= 0

 

nrhs

Number of right-hand-side vectors to solve

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1, n )

 

ipiv

Information about decomposition swaps and blocks

 

b

array of doubles with dimension (n , ldb )

 

ldb

The leading dimension of the array b , ldb >= max (1, n )

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dsytri ()

gint
ncm_lapack_dsytri (gchar uplo,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gint *ipiv,
                   NcmLapackWS *ws);

This function compute the inverse of a real symmetric indefinite matrix a using the factorization a = U*D*U**T or a = L*D*L**T computed by ncm_lapack_dsytrf().

Parameters

uplo

'U' upper triangle of a is stored; 'L' lower triangle of a is stored

 

n

The order of the matrix a , n >= 0

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1, n )

 

ipiv

Information about decomposition swaps and blocks

 

ws

a NcmLapackWS

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dsysv ()

gint
ncm_lapack_dsysv (gchar uplo,
                  gint n,
                  gint nrhs,
                  gdouble *a,
                  gint lda,
                  gint *ipiv,
                  gdouble *b,
                  gint ldb,
                  gdouble *work,
                  gint lwork);

Purpose

DSYSV uses the diagonal pivoting factorization to compute the solution to a real system of linear equations A * X = B, where A is an N-by-N symmetric matrix and X and B are N-by-NRHS matrices.

Parameters

uplo

UPLO is CHARACTER*1

 

n

N is INTEGER

 

nrhs

NRHS is INTEGER

 

a

A is DOUBLE PRECISION array, dimension (LDA,N)

 

lda

LDA is INTEGER

 

ipiv

IPIV is INTEGER array, dimension (N)

 

b

B is DOUBLE PRECISION array, dimension (LDB,NRHS)

 

ldb

LDB is INTEGER

 

work

WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))

 

lwork

LWORK is INTEGER

 

Returns

INFO is INTEGER

  • = 0: successful exit

  • < 0: if INFO = -i, the i-th argument had an illegal value

  • > 0: if INFO = i, and i is

  • <= N: D(i,i) is exactly zero. The factorization has been completed but the factor D is exactly singular, so the solution could not be computed.


ncm_lapack_dsysvx ()

gint
ncm_lapack_dsysvx (gchar fact,
                   gchar uplo,
                   gint n,
                   gint nrhs,
                   gdouble *a,
                   gint lda,
                   gdouble *af,
                   gint ldaf,
                   gint *ipiv,
                   gdouble *b,
                   gint ldb,
                   gdouble *x,
                   gint ldx,
                   gdouble *rcond,
                   gdouble *ferr,
                   gdouble *berr,
                   gdouble *work,
                   gint lwork,
                   gint *iwork);

Purpose

DSYSVX uses the diagonal pivoting factorization to compute the solution to a real system of linear equations A * X = B, where A is an N-by-N symmetric matrix and X and B are N-by-NRHS matrices.

Error bounds on the solution and a condition estimate are also provided.


Description

The following steps are performed:

  1. If FACT = 'N', the diagonal pivoting method is used to factor A. The form of the factorization is

    • A = U * D * U**T, if UPLO = 'U', or

    • A = L * D * L**T, if UPLO = 'L', where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.

  2. If some D(i,i)=0, so that D is exactly singular, then the routine returns with INFO = i. Otherwise, the factored form of A is used to estimate the condition number of the matrix A. If the reciprocal of the condition number is less than machine precision, INFO = N+1 is returned as a warning, but the routine still goes on to solve for X and compute error bounds as described below.

  3. The system of equations is solved for X using the factored form of A.

  4. Iterative refinement is applied to improve the computed solution matrix and calculate error bounds and backward error estimates for it.

Parameters

fact

FACT is CHARACTER*1

 

uplo

UPLO is CHARACTER*1

 

n

N is INTEGER

 

nrhs

NRHS is INTEGER

 

a

A is DOUBLE PRECISION array, dimension (LDA,N)

 

lda

LDA is INTEGER

 

af

AF is DOUBLE PRECISION array, dimension (LDAF,N)

 

ldaf

LDA is INTEGER

 

ipiv

IPIV is INTEGER array, dimension (N)

 

b

B is DOUBLE PRECISION array, dimension (LDB,NRHS)

 

ldb

LDB is INTEGER

 

x

X is DOUBLE PRECISION array, dimension (LDX,NRHS)

 

ldx

LDX is INTEGER

 

rcond

RCOND is DOUBLE PRECISION

 

ferr

FERR is DOUBLE PRECISION array, dimension (NRHS)

 

berr

BERR is DOUBLE PRECISION array, dimension (NRHS)

 

work

WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))

 

lwork

LWORK is INTEGER

 

iwork

IWORK is INTEGER array, dimension (N)

 

Returns

INFO is INTEGER

  • = 0: successful exit

  • < 0: if INFO = -i, the i-th argument had an illegal value

  • > 0: if INFO = i, and i is

  • <= N: D(i,i) is exactly zero. The factorization has been completed but the factor D is exactly singular, so the solution and error bounds could not be computed. RCOND = 0 is returned.

  • = N+1: D is nonsingular, but RCOND is less than machine precision, meaning that the matrix is singular to working precision. Nevertheless, the solution and error bounds are computed because there are a number of situations where the computed solution can be more accurate than the value of RCOND would suggest.


ncm_lapack_dsysvxx ()

gint
ncm_lapack_dsysvxx (gchar fact,
                    gchar uplo,
                    gint n,
                    gint nrhs,
                    gdouble *a,
                    gint lda,
                    gdouble *af,
                    gint ldaf,
                    gint *ipiv,
                    gchar *equed,
                    gdouble *s,
                    gdouble *b,
                    gint ldb,
                    gdouble *x,
                    gint ldx,
                    gdouble *rcond,
                    gdouble *rpvgrw,
                    gdouble *berr,
                    const gint n_err_bnds,
                    gdouble *err_bnds_norm,
                    gdouble *err_bnds_comp,
                    const gint nparams,
                    gdouble *params,
                    gdouble *work,
                    gint *iwork);

Purpose

DSYSVXX uses the diagonal pivoting factorization to compute the solution to a double precision system of linear equations A * X = B, where A is an N-by-N symmetric matrix and X and B are N-by-NRHS matrices.

If requested, both norm-wise and maximum component-wise error bounds are returned. DSYSVXX will return a solution with a tiny guaranteed error (O(eps) where eps is the working machine precision) unless the matrix is very ill-conditioned, in which case a warning is returned. Relevant condition numbers also are calculated and returned.

DSYSVXX accepts user-provided factorizations and equilibration factors; see the definitions of the FACT and EQUED options. Solving with refinement and using a factorization from a previous DSYSVXX call will also produce a solution with either O(eps) errors or warnings, but we cannot make that claim for general user-provided factorizations and equilibration factors if they differ from what DSYSVXX would itself produce.


Description

The following steps are performed:

  1. If FACT = 'E', double precision scaling factors are computed to equilibrate the system:

    • diag(S)*A*diag(S) *inv(diag(S))*X = diag(S)*B Whether or not the system will be equilibrated depends on the scaling of the matrix A, but if equilibration is used, A is overwritten by diag(S)*A*diag(S) and B by diag(S)*B.

  2. If FACT = 'N' or 'E', the LU decomposition is used to factor the matrix A (after equilibration if FACT = 'E') as

    • A = U * D * U**T, if UPLO = 'U', or

    • A = L * D * L**T, if UPLO = 'L', where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.

  3. If some D(i,i)=0, so that D is exactly singular, then the routine returns with INFO = i. Otherwise, the factored form of A is used to estimate the condition number of the matrix A (see argument RCOND). If the reciprocal of the condition number is less than machine precision, the routine still goes on to solve for X and compute error bounds as described below.

  4. The system of equations is solved for X using the factored form of A.

  5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero), the routine will use iterative refinement to try to get a small error and error bounds. Refinement calculates the residual to at least twice the working precision.

  6. If equilibration was used, the matrix X is premultiplied by diag(R) so that it solves the original system before equilibration.

Some optional parameters are bundled in the PARAMS array. These settings determine how refinement is performed, but often the defaults are acceptable. If the defaults are acceptable, users can pass NPARAMS = 0 which prevents the source code from accessing the PARAMS argument.

Parameters

fact

FACT is CHARACTER*1

 

uplo

UPLO is CHARACTER*1

 

n

N is INTEGER

 

nrhs

NRHS is INTEGER

 

a

A is DOUBLE PRECISION array, dimension (LDA,N)

 

lda

LDA is INTEGER

 

af

AF is DOUBLE PRECISION array, dimension (LDAF,N)

 

ldaf

LDA is INTEGER

 

ipiv

IPIV is INTEGER array, dimension (N)

 

equed

EQUED is CHARACTER*1

 

s

S is DOUBLE PRECISION array, dimension (N)

 

b

B is DOUBLE PRECISION array, dimension (LDB,NRHS)

 

ldb

LDB is INTEGER

 

x

X is DOUBLE PRECISION array, dimension (LDX,NRHS)

 

ldx

LDX is INTEGER

 

rcond

RCOND is DOUBLE PRECISION

 

rpvgrw

RPVGRW is DOUBLE PRECISION

 

berr

BERR is DOUBLE PRECISION array, dimension (NRHS)

 

n_err_bnds

N_ERR_BNDS is INTEGER

 

err_bnds_norm

ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)

 

err_bnds_comp

ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)

 

nparams

NPARAMS is INTEGER

 

params

PARAMS is DOUBLE PRECISION array, dimension (NPARAMS)

 

work

WORK is DOUBLE PRECISION array, dimension (4*N)

 

iwork

IWORK is INTEGER array, dimension (N)

 

Returns

INFO is INTEGER

  • = 0: Successful exit. The solution to every right-hand side is guaranteed.

  • < 0: If INFO = -i, the i-th argument had an illegal value

  • > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution and error bounds could not be computed. RCOND = 0 is returned.

  • = N+J: The solution corresponding to the Jth right-hand side is not guaranteed. The solutions corresponding to other right- hand sides K with K > J may not be guaranteed as well, but only the first such right-hand side is reported. If a small componentwise error is not requested (PARAMS(3) = 0.0) then the Jth right-hand side is the first with a normwise error bound that is not guaranteed (the smallest J such that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0) the Jth right-hand side is the first with either a normwise or componentwise error bound that is not guaranteed (the smallest J such that either ERR_BNDS_NORM(J,1) = 0.0 or ERR_BNDS_COMP(J,1) = 0.0). See the definition of ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information about all of the right-hand sides check ERR_BNDS_NORM or ERR_BNDS_COMP.


ncm_lapack_dsyevr ()

gint
ncm_lapack_dsyevr (gchar jobz,
                   gchar range,
                   gchar uplo,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gdouble vl,
                   gdouble vu,
                   gint il,
                   gint iu,
                   gdouble abstol,
                   gint *m,
                   gdouble *w,
                   gdouble *z,
                   gint ldz,
                   gint *isuppz,
                   NcmLapackWS *ws);

Computes selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix a .

Parameters

jobz

a char with value 'N', 'V' or 'I'

 

range

a char with value 'A', 'V' or 'I'

 

uplo

a char with value 'U' or 'L'

 

n

an integer with the order of the matrix a

 

a

a double precision array with dimension (n , lda )

 

lda

an integer with the leading dimension of the array a , lda >= max (1, n )

 

vl

a double precision with the lower bound of the interval to be searched for eigenvalues

 

vu

a double precision with the upper bound of the interval to be searched for eigenvalues

 

il

an integer with the index of the smallest eigenvalue to be returned

 

iu

an integer with the index of the largest eigenvalue to be returned

 

abstol

a double precision with the absolute error tolerance for the eigenvalues

 

m

an integer with the total number of eigenvalues found

 

w

a double precision array with dimension n

 

z

a double precision array with dimension (ldz , n )

 

ldz

an integer with the leading dimension of the array z , ldz >= 1, and if jobz = 'V' or 'I', ldz >= n

 

isuppz

an integer array with dimension (2, n )

 

ws

a NcmLapackWS

 

Returns

an integer with the error code


ncm_lapack_dsyevd ()

gint
ncm_lapack_dsyevd (gchar jobz,
                   gchar uplo,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gdouble *w,
                   NcmLapackWS *ws);

Computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix a .

Parameters

jobz

a char with value 'N', 'V' or 'I'

 

uplo

a char with value 'U' or 'L'

 

n

an integer with the order of the matrix a

 

a

a double precision array with dimension (n , lda )

 

lda

an integer with the leading dimension of the array a , lda >= max (1, n )

 

w

a double precision array with dimension n

 

ws

a NcmLapackWS

 

Returns

an integer with the error code


ncm_lapack_dgeev ()

gint
ncm_lapack_dgeev (gchar jobvl,
                  gchar jobvr,
                  gint n,
                  gdouble *a,
                  gint lda,
                  gdouble *wr,
                  gdouble *wi,
                  gdouble *vl,
                  gint ldvl,
                  gdouble *vr,
                  gint ldvr,
                  gdouble *work,
                  gint lwork);

This function computes the eigensystem for a real matrix a = A.

Calling this function with lwork == -1 computed the ideal lwork in work [0].

Parameters

jobvl

n left eigenvectors of a are not computed, 'V' left eigenvectors of a are computed

 

jobvr

n right eigenvectors of a are not computed, 'V' right eigenvectors of a are computed

 

n

The order of the matrix a , n >= 0

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1, n )

 

wr

contain the real part of the computed eigenvalues

 

wi

contain the imaginary part of the computed eigenvalues

 

vl

if jobvl = 'V', the left eigenvectors $u(j)$ are stored one after another in the rows of vl , in the same order as their eigenvalues

 

ldvl

the leading dimension of the array vl

 

vr

if jobvr = 'V', the left eigenvectors $v(j)$ are stored one after another in the rows of vr , in the same order as their eigenvalues

 

ldvr

the leading dimension of the array vr

 

work

work area, must have lwork allocated doubles

 

lwork

work area size

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dgeevx ()

gint
ncm_lapack_dgeevx (gchar balanc,
                   gchar jobvl,
                   gchar jobvr,
                   gchar sense,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gdouble *wr,
                   gdouble *wi,
                   gdouble *vl,
                   gint ldvl,
                   gdouble *vr,
                   gint ldvr,
                   gint *ilo,
                   gint *ihi,
                   gdouble *scale,
                   gdouble *abnrm,
                   gdouble *rconde,
                   gdouble *rcondv,
                   gdouble *work,
                   gint lwork,
                   gint *iwork);

This function computes the eigensystem for a real matrix a = A.

Calling this function with lwork == -1 computed the ideal lwork in work [0].

Parameters

balanc

a char with value 'N', 'P', 'S' or 'B'

 

jobvl

n left eigenvectors of a are not computed, 'V' left eigenvectors of a are computed

 

jobvr

n right eigenvectors of a are not computed, 'V' right eigenvectors of a are computed

 

sense

a char with value 'N', 'E', 'V' or 'B'

 

n

The order of the matrix a , n >= 0

 

a

array of doubles with dimension (n , lda )

 

lda

The leading dimension of the array a , lda >= max (1, n )

 

wr

contain the real part of the computed eigenvalues

 

wi

contain the imaginary part of the computed eigenvalues

 

vl

if jobvl = 'V', the left eigenvectors $u(j)$ are stored one after another in the rows of vl , in the same order as their eigenvalues

 

ldvl

the leading dimension of the array vl

 

vr

if jobvr = 'V', the left eigenvectors $v(j)$ are stored one after another in the rows of vr , in the same order as their eigenvalues

 

ldvr

the leading dimension of the array vr

 

ilo

an integer with the index of the first eigenvalue to be returned

 

ihi

an integer with the index of the last eigenvalue to be returned

 

scale

an array of doubles with dimension n

 

abnrm

an array of doubles with dimension n

 

rconde

an array of doubles with dimension n

 

rcondv

an array of doubles with dimension n

 

work

work area, must have lwork allocated doubles

 

lwork

work area size

 

iwork

work area, must have liwork allocated integers

 

Returns

i = 0: successful exit

< 0: -i, the i-th argument had an illegal value

> 0: the (i,i) element of the factor U or L is zero, and the inverse could not be computed.


ncm_lapack_dgeqrf ()

gint
ncm_lapack_dgeqrf (gint m,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gdouble *tau,
                   NcmLapackWS *ws);

DGEQRF computes a QR factorization of a real M-by-N matrix A: A = Q * R.

Parameters

m

M is INTEGER The number of rows of the matrix A. M >= 0.

 

n

N is INTEGER The number of columns of the matrix A. N >= 0.

 

a

A is DOUBLE PRECISION array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors (see Further Details).

 

lda

LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M).

 

tau

TAU is DOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).

 

ws

a NcmLapackWS

 

Returns

= 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value


ncm_lapack_dgerqf ()

gint
ncm_lapack_dgerqf (gint m,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gdouble *tau,
                   NcmLapackWS *ws);

DGERQF computes a RQ factorization of a real M-by-N matrix A: A = R * Q.

Parameters

m

M is INTEGER The number of rows of the matrix A. M >= 0.

 

n

N is INTEGER The number of columns of the matrix A. N >= 0.

 

a

A is DOUBLE PRECISION array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors (see Further Details).

 

lda

LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M).

 

tau

TAU is DOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).

 

ws

a NcmLapackWS

 

Returns

= 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value


ncm_lapack_dgeqlf ()

gint
ncm_lapack_dgeqlf (gint m,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gdouble *tau,
                   NcmLapackWS *ws);

DGEQLF computes a QL factorization of a real M-by-N matrix A: A = Q * L.

Parameters

m

M is INTEGER The number of rows of the matrix A. M >= 0.

 

n

N is INTEGER The number of columns of the matrix A. N >= 0.

 

a

A is DOUBLE PRECISION array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors (see Further Details).

 

lda

LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M).

 

tau

TAU is DOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).

 

ws

a NcmLapackWS

 

Returns

= 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value


ncm_lapack_dgelqf ()

gint
ncm_lapack_dgelqf (gint m,
                   gint n,
                   gdouble *a,
                   gint lda,
                   gdouble *tau,
                   NcmLapackWS *ws);

DGELQF computes a LQ factorization of a real M-by-N matrix A: A = L * Q.

Parameters

m

M is INTEGER The number of rows of the matrix A. M >= 0.

 

n

N is INTEGER The number of columns of the matrix A. N >= 0.

 

a

A is DOUBLE PRECISION array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array TAU, represent the orthogonal matrix Q as a product of min(m,n) elementary reflectors (see Further Details).

 

lda

LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M).

 

tau

TAU is DOUBLE PRECISION array, dimension (min(M,N)) The scalar factors of the elementary reflectors (see Further Details).

 

ws

a NcmLapackWS

 

Returns

= 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value


ncm_lapack_dggglm_alloc ()

GArray *
ncm_lapack_dggglm_alloc (NcmMatrix *L,
                         NcmMatrix *X,
                         NcmVector *p,
                         NcmVector *d,
                         NcmVector *y);

Calculates and allocs memory to solve the system determined by the parameters.

This function is expect the matrix X and L to be row-major.

Parameters

L

a NcmMatrix

 

X

a NcmMatrix

 

p

a NcmVector

 

d

a NcmVector

 

y

a NcmVector

 

Returns

the newly allocated workspace.

[transfer full][array][element-type double]


ncm_lapack_dggglm_run ()

gint
ncm_lapack_dggglm_run (GArray *ws,
                       NcmMatrix *L,
                       NcmMatrix *X,
                       NcmVector *p,
                       NcmVector *d,
                       NcmVector *y);

Runs the dggglm function using the workspace ws .

This function is expect the matrix X and L to be row-major.

Parameters

ws

a workspace.

[in][array][element-type double]

L

a NcmMatrix

 

X

a NcmMatrix

 

p

a NcmVector

 

d

a NcmVector

 

y

a NcmVector

 

ncm_lapack_dgels ()

gint
ncm_lapack_dgels (gchar trans,
                  const gint m,
                  const gint n,
                  const gint nrhs,
                  gdouble *a,
                  const gint lda,
                  gdouble *b,
                  const gint ldb,
                  double *work,
                  const gint lwork);

DGELS solves overdetermined or underdetermined real linear systems involving an M-by-N matrix A, or its transpose, using a QR or LQ factorization of A. It is assumed that A has full rank.

The following options are provided:

  1. If TRANS = 'N' and m >= n: find the least squares solution of an overdetermined system, i.e., solve the least squares problem minimize || B - A*X ||.

  2. If TRANS = 'N' and m < n: find the minimum norm solution of an underdetermined system A * X = B.

  3. If TRANS = 'T' and m >= n: find the minimum norm solution of an underdetermined system A**T * X = B.

  4. If TRANS = 'T' and m < n: find the least squares solution of an overdetermined system, i.e., solve the least squares problem minimize || B - A**T * X ||.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix X.

Parameters

trans

is a char = 'N': the linear system involves A = 'T': the linear system involves A**T

 

m

is an integer. The number of rows of the matrix A. M >= 0

 

n

is an integer. The number of columns of the matrix A. N >= 0

 

nrhs

is an integer. The number of right hand sides, i.e., the number of columns of the matrices B and X. NRHS >=0

 

a

array of doubles with dimension (n , lda ) On entry, the M-by-N matrix A. On exit, if M >= N, A is overwritten by details of its QR factorization as returned by DGEQRF if M < N, A is overwritten by details of its LQ factorization as returned by DGELQF

 

lda

The leading dimension of the array a , lda >= max (1,n )

 

b

array of doubles with dimension (n , ldb ) On entry, the matrix B of right hand side vectors, stored columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS if TRANS = 'T'. On exit, if INFO = 0, B is overwritten by the solution vectors, stored columnwise: if TRANS = 'N' and m >= n, rows 1 to n of B contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of elements N+1 to M in that column; if TRANS = 'N' and m < n, rows 1 to N of B contain the minimum norm solution vectors; if TRANS = 'T' and m >= n, rows 1 to M of B contain the minimum norm solution vectors; if TRANS = 'T' and m < n, rows 1 to M of B contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of elements M+1 to N in that column.

 

ldb

The leading dimension of the array b , ldb >= max (1, n )

 

work

WORK is DOUBLE PRECISION array, dimension (4*N)

 

lwork

LWORK is INTEGER

 

ncm_lapack_dgelsd ()

gint
ncm_lapack_dgelsd (const gint m,
                   const gint n,
                   const gint nrhs,
                   gdouble *a,
                   const gint lda,
                   gdouble *b,
                   const gint ldb,
                   gdouble *s,
                   gdouble *rcond,
                   gint *rank,
                   NcmLapackWS *ws);

DGELSD computes the minimum-norm solution to a real linear least squares problem: minimize 2-norm(| b - A*x |) using the singular value decomposition (SVD) of A. A is an M-by-N matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix X.

The problem is solved in three steps: (1) Reduce the coefficient matrix A to bidiagonal form with Householder transformations, reducing the original problem into a "bidiagonal least squares problem" (BLS) (2) Solve the BLS using a divide and conquer approach. (3) Apply back all the Householder transformations to solve the original least squares problem.

The effective rank of A is determined by treating as zero those singular values which are less than RCOND times the largest singular value.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters

m

is an integer. The number of rows of the matrix A. M >= 0

 

n

is an integer. The number of columns of the matrix A. N >= 0

 

nrhs

is an integer. The number of right hand sides, i.e., the number of columns of the matrices B and X. NRHS >=0

 

a

array of doubles with dimension (n , lda ) On entry, the M-by-N matrix A. On exit, if M >= N, A is overwritten by details of its QR factorization as returned by DGEQRF if M < N, A is overwritten by details of its LQ factorization as returned by DGELQF

 

lda

The leading dimension of the array a , lda >= max (1,n )

 

b

array of doubles with dimension (n , ldb ) On entry, the matrix B of right hand side vectors, stored columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS if TRANS = 'T'. On exit, if INFO = 0, B is overwritten by the solution vectors, stored columnwise: if TRANS = 'N' and m >= n, rows 1 to n of B contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of elements N+1 to M in that column; if TRANS = 'N' and m < n, rows 1 to N of B contain the minimum norm solution vectors; if TRANS = 'T' and m >= n, rows 1 to M of B contain the minimum norm solution vectors; if TRANS = 'T' and m < n, rows 1 to M of B contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of elements M+1 to N in that column.

 

ldb

The leading dimension of the array b , ldb >= max (1, n )

 

s

S is DOUBLE PRECISION array, dimension (min(M,N)) The singular values of A in decreasing order. The condition number of A in the 2-norm = S(1)/S(min(m,n)).

 

rcond

RCOND is DOUBLE PRECISION RCOND is used to determine the effective rank of A. Singular values S(i) <= RCOND*S(1) are treated as zero. If RCOND < 0, machine precision is used instead.

 

rank

RANK is INTEGER The effective rank of A, i.e., the number of singular values which are greater than RCOND*S(1).

 

ws

a Lapack workspace object NcmLapackWS

 

NCM_LAPACK_CHECK_INFO()

#define NCM_LAPACK_CHECK_INFO(func,info) G_STMT_START { if ((info) != 0) g_error ("# NcmLapack[%s] error %4d", func, (info)); } G_STMT_END