OpenFAST
Wind turbine multiphysics simulator
Data Types | Functions/Subroutines
nwtc_lapack Module Reference

This code provides a wrapper for the LAPACK routines currently used at the NWTC (mainly codes in the FAST framework). More...

Data Types

interface  lapack_axpy
 
interface  lapack_copy
 straight-up lapack routines (from ExtPtfm_MCKF): More...
 
interface  lapack_gbsv
 Computes the solution to system of linear equations A * X = B for GB matrices. More...
 
interface  lapack_gemm
 Computes scalar1*op( A )*op( B ) + scalar2*C where op(x) = x or op(x) = x**T for matrices A, B, and C. More...
 
interface  lapack_gemv
 
interface  lapack_gesv
 Computes the solution to system of linear equations A * X = B for GE matrices. More...
 
interface  lapack_gesvd
 Compute the SVD for a general matrix A = USV^T. More...
 
interface  lapack_getrf
 Factor matrix into A=PLU. More...
 
interface  lapack_getri
 Compute the inverse of a matrix using the LU factorization. More...
 
interface  lapack_getrs
 Solve system(s) of linear equations Ax=PLUx=b. More...
 
interface  lapack_ggev
 Compute generalized eigenvalues and/or eigenvectors for a pair of N-by-N real nonsymmetric matrices (A,B). More...
 
interface  lapack_posv
 Compute the solution to system of linear equations A * X = B for PO matrices. More...
 
interface  lapack_pptrf
 Compute the Cholesky factorization of a real symmetric positive definite matrix A stored in packed format. More...
 
interface  lapack_tpttr
 Unpack packed (1D) to regular matrix format (2D) More...
 

Functions/Subroutines

subroutine lapack_dgbsv (N, KL, KU, NRHS, AB, IPIV, B, ErrStat, ErrMsg)
 general banded solve: Computes the solution to system of linear equations A * X = B for GB (general, banded) matrices. More...
 
subroutine lapack_sgbsv (N, KL, KU, NRHS, AB, IPIV, B, ErrStat, ErrMsg)
 general banded solve: Computes the solution to system of linear equations A * X = B for GB (general, banded) matrices. More...
 
subroutine lapack_dgemm (TRANSA, TRANSB, ALPHA, A, B, BETA, C, ErrStat, ErrMsg)
 general matrix multiply: computes C = alpha*op( A )*op( B ) + beta*C where op(x) = x or op(x) = x**T for matrices A, B, and C use LAPACK_GEMM (nwtc_lapack::lapack_gemm) instead of this specific function. More...
 
subroutine lapack_sgemm (TRANSA, TRANSB, ALPHA, A, B, BETA, C, ErrStat, ErrMsg)
 general matrix multiply: computes C = alpha*op( A )*op( B ) + beta*C where op(x) = x or op(x) = x**T for matrices A, B, and C use LAPACK_GEMM (nwtc_lapack::lapack_gemm) instead of this specific function. More...
 
subroutine lapack_dgesv (N, A, IPIV, B, ErrStat, ErrMsg)
 general solve: Computes the solution to system of linear equations A * X = B for GE matrices. More...
 
subroutine lapack_sgesv (N, A, IPIV, B, ErrStat, ErrMsg)
 general solve: Computes the solution to system of linear equations A * X = B for GE matrices. More...
 
subroutine lapack_dgetrf (M, N, A, IPIV, ErrStat, ErrMsg)
 general matrix factorization: Factor matrix into A=PLU. More...
 
subroutine lapack_sgetrf (M, N, A, IPIV, ErrStat, ErrMsg)
 general matrix factorization: Factor matrix into A=PLU. More...
 
subroutine lapack_dgetrs (TRANS, N, A, IPIV, B, ErrStat, ErrMsg)
 general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b. More...
 
subroutine lapack_dgetrs1 (TRANS, N, A, IPIV, B, ErrStat, ErrMsg)
 general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b. More...
 
subroutine lapack_sgetrs (TRANS, N, A, IPIV, B, ErrStat, ErrMsg)
 general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b. More...
 
subroutine lapack_sgetrs1 (TRANS, N, A, IPIV, B, ErrStat, ErrMsg)
 general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b. More...
 
subroutine lapack_dgetri (N, A, IPIV, WORK, LWORK, ErrStat, ErrMsg)
 Compute the inverse of a general matrix using the LU factorization. More...
 
subroutine lapack_sgetri (N, A, IPIV, WORK, LWORK, ErrStat, ErrMsg)
 Compute the inverse of a general matrix using the LU factorization. More...
 
subroutine lapack_dggev (JOBVL, JOBVR, N, A, B, ALPHAR, ALPHAI, BETA, VL, VR, WORK, LWORK, ErrStat, ErrMsg)
 Compute generalized eigenvalues and/or eigenvectors for a pair of N-by-N real nonsymmetric matrices (A,B). More...
 
subroutine lapack_sggev (JOBVL, JOBVR, N, A, B, ALPHAR, ALPHAI, BETA, VL, VR, WORK, LWORK, ErrStat, ErrMsg)
 Compute generalized eigenvalues and/or eigenvectors for a pair of N-by-N real nonsymmetric matrices (A,B). More...
 
subroutine lapack_dposv (UPLO, N, NRHS, A, B, ErrStat, ErrMsg)
 Compute the solution to system of linear equations A * X = B for PO (positive-definite) matrices. More...
 
subroutine lapack_sposv (UPLO, N, NRHS, A, B, ErrStat, ErrMsg)
 Compute the solution to system of linear equations A * X = B for PO (positive-definite) matrices. More...
 
subroutine lapack_dpptrf (UPLO, N, AP, ErrStat, ErrMsg)
 Compute the Cholesky factorization of a real symmetric positive definite matrix A stored in packed format. More...
 
subroutine lapack_spptrf (UPLO, N, AP, ErrStat, ErrMsg)
 Compute the Cholesky factorization of a real symmetric positive definite matrix A stored in packed format. More...
 
subroutine lapack_dgesvd (JOBU, JOBVT, M, N, A, S, U, VT, WORK, LWORK, ErrStat, ErrMsg)
 Compute singular value decomposition (SVD) for a general matrix, A. More...
 
subroutine lapack_sgesvd (JOBU, JOBVT, M, N, A, S, U, VT, WORK, LWORK, ErrStat, ErrMsg)
 Compute singular value decomposition (SVD) for a general matrix, A. More...
 
subroutine lapack_dtpttr (UPLO, N, AP, A, LDA, ErrStat, ErrMsg)
 Unpack a by-column-packed array into a 2D matrix format See documentation in DTPTTR/STPTTR source code. More...
 
subroutine lapack_stpttr (UPLO, N, AP, A, LDA, ErrStat, ErrMsg)
 Unpack a by-column-packed array into a 2D matrix format. More...
 

Detailed Description

This code provides a wrapper for the LAPACK routines currently used at the NWTC (mainly codes in the FAST framework).

This enables us to call generic routines (not single- or double-precision specific ones) so that we don't have to change source code to compile in double vs. single precision.

Function/Subroutine Documentation

◆ lapack_dgbsv()

subroutine nwtc_lapack::lapack_dgbsv ( integer, intent(in)  N,
integer, intent(in)  KL,
integer, intent(in)  KU,
integer, intent(in)  NRHS,
real(r8ki), dimension( :, : ), intent(inout)  AB,
integer, dimension( : ), intent(out)  IPIV,
real(r8ki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general banded solve: Computes the solution to system of linear equations A * X = B for GB (general, banded) matrices.

use LAPACK_GBSV (nwtc_lapack::lapack_gbsv) instead of this specific function.

Parameters
[in]klThe number of subdiagonals within the band of A. KL >= 0.
[in]kuThe number of superdiagonals within the band of A. KU >= 0.
[in]nThe number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]abOn entry, the matrix A in band storage, in rows KL+1 to 2*KL+KU+1; rows 1 to KL of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows: AB(KL+KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+KL) On exit, details of the factorization: U is stored as an upper triangular band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and the multipliers used during the factorization are stored in rows KL+KU+2 to 2*KL+KU+1.
[in,out]bOn entry, the N-by-NRHS right hand side matrix B. On exit, if INFO = 0, the N-by-NRHS solution matrix X.
[out]ipivThe pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dgemm()

subroutine nwtc_lapack::lapack_dgemm ( character(1), intent(in)  TRANSA,
character(1), intent(in)  TRANSB,
real(r8ki), intent(in)  ALPHA,
real(r8ki), dimension( :, : ), intent(in)  A,
real(r8ki), dimension( :, : ), intent(in)  B,
real(r8ki), intent(in)  BETA,
real(r8ki), dimension( :, : ), intent(inout)  C,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general matrix multiply: computes C = alpha*op( A )*op( B ) + beta*C where op(x) = x or op(x) = x**T for matrices A, B, and C use LAPACK_GEMM (nwtc_lapack::lapack_gemm) instead of this specific function.

Parameters
[in]transaOn entry, TRANSA specifies the form of op( A ) to be used in the matrix multiplication as follows: TRANSA = 'N' or 'n', op( A ) = A. TRANSA = 'T' or 't', op( A ) = A**T.
[in]transbOn entry, TRANSB specifies the form of op( A ) to be used in the matrix multiplication as follows: TRANSB = 'N' or 'n', op( B ) = B. TRANSB = 'T' or 't', op( B ) = B**T.
[in]alphaOn entry, ALPHA specifies the scalar alpha.
[in]betaOn entry, BETA specifies the scalar beta. When BETA is supplied as zero then C need not be set on input.
[in]aMatrix A
[in]bMatrix B
[in,out]cMatrix C: Before entry, C must contain the matrix C, except when beta is zero, in which case C need not be set on entry. On exit, the array C is overwritten by the m by n matrix ( alpha*op( A )*op( B ) + beta*C ).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dgesv()

subroutine nwtc_lapack::lapack_dgesv ( integer, intent(in)  N,
real(r8ki), dimension( :, : ), intent(inout)  A,
integer, dimension( : ), intent(out)  IPIV,
real(r8ki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general solve: Computes the solution to system of linear equations A * X = B for GE matrices.

use LAPACK_GESV (nwtc_lapack::lapack_gesv) instead of this specific function.

Parameters
[in]nThe number of linear equations, i.e., the order of the matrix A. N >= 0.
[in,out]aOn entry, the N-by-N coefficient matrix A. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in,out]bOn entry, the N-by-NRHS matrix of right hand side matrix B. On exit, if INFO = 0, the N-by-NRHS solution matrix X.
[out]ipivThe pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dgesvd()

subroutine nwtc_lapack::lapack_dgesvd ( character(1), intent(in)  JOBU,
character(1), intent(in)  JOBVT,
integer, intent(in)  M,
integer, intent(in)  N,
real(r8ki), dimension( :, : ), intent(inout)  A,
real(r8ki), dimension( : ), intent(out)  S,
real(r8ki), dimension( :, : ), intent(out)  U,
real(r8ki), dimension( :, : ), intent(out)  VT,
real(r8ki), dimension( : ), intent(inout)  WORK,
integer, intent(in)  LWORK,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute singular value decomposition (SVD) for a general matrix, A.

use LAPACK_DGESVD (nwtc_lapack::lapack_dgesvd) instead of this specific function.

Parameters
[in]jobu'A': all M columns of U are returned in array U; 'S': the first min(m,n) columns of U (the left singular vectors) are returned in the array U; 'O': the first min(m,n) columns of U (the left singular vectors) are overwritten on the array A; 'N': no columns of U (no left singular vectors) are computed.
[in]jobvt'A': all N rows of V^T are returned in the array VT; 'S': the first min(m,n) rows of V^T (the right singular vectors) are returned in the array VT; 'O': the first min(m,n) rows of V^T (the right singular vectors) are overwritten on the array A; 'N': no rows of V**T (no right singular vectors) are computed.
[in]mThe number of rows of the input matrix A. M >= 0.
[in]nThe number of columns of the input matrix A. N >= 0.
[in,out]aA is DOUBLE PRECISION array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, if JOBU = 'O', A is overwritten with the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBVT = 'O', A is overwritten with the first min(m,n) rows of V**T (the right singular vectors, stored rowwise); if JOBU .ne. 'O' and JOBVT .ne. 'O', the contents of A are destroyed.
[out]sS is DOUBLE PRECISION array, dimension (min(M,N)) The singular values of A, sorted so that S(i) >= S(i+1).
[out]uU is DOUBLE PRECISION array, dimension (LDU,UCOL) (LDU,M) if JOBU = 'A' or (LDU,min(M,N)) if JOBU = 'S'. If JOBU = 'A', U contains the M-by-M orthogonal matrix U; if JOBU = 'S', U contains the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBU = 'N' or 'O', U is not referenced.
[out]vtVT is DOUBLE PRECISION array, dimension (LDVT,N) If JOBVT = 'A', VT contains the N-by-N orthogonal matrix V**T; if JOBVT = 'S', VT contains the first min(m,n) rows of V**T (the right singular vectors, stored rowwise); if JOBVT = 'N' or 'O', VT is not referenced.
[in,out]workdimension (MAX(1,LWORK)). On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]lworkThe dimension of the array WORK. LWORK >= MAX(1,5*MIN(M,N)) for the paths (see comments inside code):
  • PATH 1 (M much larger than N, JOBU='N')
  • PATH 1t (N much larger than M, JOBVT='N') LWORK >= MAX(1,3*MIN(M,N) + MAX(M,N),5*MIN(M,N)) for the other paths For good performance, LWORK should generally be larger. If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dgetrf()

subroutine nwtc_lapack::lapack_dgetrf ( integer, intent(in)  M,
integer, intent(in)  N,
real(r8ki), dimension( :, : ), intent(inout)  A,
integer, dimension( : ), intent(out)  IPIV,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general matrix factorization: Factor matrix into A=PLU.

use LAPACK_GETRF (nwtc_lapack::lapack_getrf) instead of this specific function.

Parameters
[in]mThe number of rows of the matrix A. M >= 0.
[in]nThe number of columns of the matrix A. N >= 0.
[in,out]aOn entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[out]ipivThe pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dgetri()

subroutine nwtc_lapack::lapack_dgetri ( integer, intent(in)  N,
real(r8ki), dimension( :, : ), intent(inout)  A,
integer, dimension( : ), intent(in)  IPIV,
real(r8ki), dimension( : ), intent(inout)  WORK,
integer, intent(in)  LWORK,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute the inverse of a general matrix using the LU factorization.

use LAPACK_GETRI (nwtc_lapack::lapack_getri) instead of this specific function.

Parameters
[in]nThe order of the matrix A. N >= 0.
[in]lworkThe dimension of the array WORK. LWORK >= max(1,N). For optimal performance LWORK >= N*NB, where NB is the optimal blocksize returned by ILAENV. If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.
[in]ipivdimension (N). The pivot indices from DGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).
[in,out]aOn entry, the factors L and U from the factorization A = P*L*U as computed by DGETRF. On exit, if INFO = 0, the inverse of the original matrix A.
[in,out]workOn exit, if INFO=0, then WORK(1) returns the optimal LWORK.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dgetrs()

subroutine nwtc_lapack::lapack_dgetrs ( character(1), intent(in)  TRANS,
integer, intent(in)  N,
real(r8ki), dimension( :, : ), intent(in)  A,
integer, dimension( : ), intent(in)  IPIV,
real(r8ki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b.

use LAPACK_GETRS (nwtc_lapack::lapack_getrs) instead of this specific function.

Parameters
[in]transSpecifies the form of the system of equations: = 'N': A * X = B (No transpose) = 'T': A**T* X = B (Transpose) = 'C': A**T* X = B (Conjugate transpose = Transpose)
[in]nThe order of the matrix A. N >= 0.
[in]ipivThe pivot indices from DGETRF (nwtc_lapack::lapack_getrf); for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).
[in]aThe factors L and U from the factorization A = P*L*U as computed by DGETRF.
[in,out]bOn entry, the right hand side matrix B. On exit, the solution matrix X.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dgetrs1()

subroutine nwtc_lapack::lapack_dgetrs1 ( character(1), intent(in)  TRANS,
integer, intent(in)  N,
real(r8ki), dimension( :, : ), intent(in)  A,
integer, dimension( : ), intent(in)  IPIV,
real(r8ki), dimension( : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b.

use LAPACK_GETRS (nwtc_lapack::lapack_getrs) instead of this specific function.

Parameters
[in]transSpecifies the form of the system of equations: = 'N': A * X = B (No transpose) = 'T': A**T* X = B (Transpose) = 'C': A**T* X = B (Conjugate transpose = Transpose)
[in]nThe order of the matrix A. N >= 0.
[in]ipivThe pivot indices from DGETRF (nwtc_lapack::lapack_getrf); for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).
[in]aThe factors L and U from the factorization A = P*L*U as computed by DGETRF.
[in,out]bOn entry, the right hand side matrix B. On exit, the solution matrix X.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dggev()

subroutine nwtc_lapack::lapack_dggev ( character(1), intent(in)  JOBVL,
character(1), intent(in)  JOBVR,
integer, intent(in)  N,
real(r8ki), dimension( :, : ), intent(inout)  A,
real(r8ki), dimension( :, : ), intent(inout)  B,
real(r8ki), dimension( : ), intent(out)  ALPHAR,
real(r8ki), dimension( : ), intent(out)  ALPHAI,
real(r8ki), dimension( : ), intent(out)  BETA,
real(r8ki), dimension( :, : ), intent(out)  VL,
real(r8ki), dimension( :, : ), intent(out)  VR,
real(r8ki), dimension( : ), intent(inout)  WORK,
integer, intent(in)  LWORK,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute generalized eigenvalues and/or eigenvectors for a pair of N-by-N real nonsymmetric matrices (A,B).

use LAPACK_GGEV (nwtc_lapack::lapack_ggev) instead of this specific function.

Parameters
[in]jobvl= 'N': do not compute the left generalized eigenvectors; = 'V': compute the left generalized eigenvectors.
[in]jobvr= 'N': do not compute the right generalized eigenvectors; = 'V': compute the right generalized eigenvectors.
[in]nThe order of the matrices A, B, VL, and VR. N >= 0.
[in]lworkThe dimension of the array WORK. LWORK >= max(1,8*N). For good performance, LWORK must generally be larger. If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.
[in,out]adimension (LDA, N). On entry, the matrix A in the pair (A,B). On exit, A has been overwritten.
[in,out]bdimension (LDB, N). On entry, the matrix B in the pair (A,B). On exit, B has been overwritten.
[out]alphardimension (N). See comments for variable "Beta"
[out]alphaidimension (N). See comments for variable "Beta".
[out]betaOn exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will be the generalized eigenvalues. If ALPHAI(j) is zero, then the j-th eigenvalue is real; if positive, then the j-th and (j+1)-st eigenvalues are a complex conjugate pair, with ALPHAI(j+1) negative.

Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j) may easily over- or underflow, and BETA(j) may even be zero. Thus, the user should avoid naively computing the ratio alpha/beta. However, ALPHAR and ALPHAI will be always less than and usually comparable with norm(A) in magnitude, and BETA always less than and usually comparable with norm(B).

Parameters
[out]vldimension (LDVL,N). If JOBVL = 'V', the left eigenvectors u(j) are stored one after another in the columns of VL, in the same order as their eigenvalues. If the j-th eigenvalue is real, then u(j) = VL(:,j), the j-th column of VL. If the j-th and (j+1)-th eigenvalues form a complex conjugate pair, then u(j) = VL(:,j)+i*VL(:,j+1) and u(j+1) = VL(:,j)-i*VL(:,j+1). Each eigenvector is scaled so the largest component has abs(real part)+abs(imag. part)=1. Not referenced if JOBVL = 'N'.
[out]vrdimension (LDVR,N). If JOBVR = 'V', the right eigenvectors v(j) are stored one after another in the columns of VR, in the same order as their eigenvalues. If the j-th eigenvalue is real, then v(j) = VR(:,j), the j-th column of VR. If the j-th and (j+1)-th eigenvalues form a complex conjugate pair, then v(j) = VR(:,j)+i*VR(:,j+1) and v(j+1) = VR(:,j)-i*VR(:,j+1). Each eigenvector is scaled so the largest component has abs(real part)+abs(imag. part)=1. Not referenced if JOBVR = 'N'.
[in,out]workdimension (MAX(1,LWORK)). On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_dposv()

subroutine nwtc_lapack::lapack_dposv ( character(1), intent(in)  UPLO,
integer, intent(in)  N,
integer, intent(in)  NRHS,
real(r8ki), dimension( :, : ), intent(inout)  A,
real(r8ki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute the solution to system of linear equations A * X = B for PO (positive-definite) matrices.

use LAPACK_POSV (nwtc_lapack::lapack_posv) instead of this specific function.

Parameters
[in]nThe number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]aOn entry, the symmetric matrix A. If UPLO = 'U', the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, if INFO = 0, the factor U or L from the Cholesky factorization A = U**T*U or A = L*L**T.
[in,out]bOn entry, the N-by-NRHS matrix of right hand side matrix B. On exit, if INFO = 0, the N-by-NRHS solution matrix X.
[out]errstatError level
[out]errmsgMessage describing error
[in]uplo'U': Upper triangle of A is stored; 'L': Lower triangle of A is stored.

◆ lapack_dpptrf()

subroutine nwtc_lapack::lapack_dpptrf ( character(1), intent(in)  UPLO,
integer, intent(in)  N,
real(r8ki), dimension( : ), intent(inout)  AP,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute the Cholesky factorization of a real symmetric positive definite matrix A stored in packed format.

use LAPACK_PPTRF (nwtc_lapack::lapack_pptrf) instead of this specific function.

Parameters
[in]nThe order of the matrix A. N >= 0.
[in,out]apAP is REAL array, dimension (N*(N+1)/2) On entry, the upper or lower triangle of the symmetric matrix A, packed columnwise in a linear array. The j-th column of A is stored in the array AP as follows: if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j; if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n. See below for further details. On exit, if INFO = 0, the triangular factor U or L from the Cholesky factorization A = U**T*U or A = L*L**T, in the same storage format as A.

Further details: The packed storage scheme is illustrated by the following example when N = 4, UPLO = 'U':

Two-dimensional storage of the symmetric matrix A:

a11 a12 a13 a14 a22 a23 a24 a33 a34 (aij = aji) a44

Packed storage of the upper triangle of A:

AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]

Parameters
[out]errstatError level
[out]errmsgMessage describing error
[in]uplo'U': Upper triangle of A is stored; 'L': Lower triangle of A is stored.

◆ lapack_dtpttr()

subroutine nwtc_lapack::lapack_dtpttr ( character(1), intent(in)  UPLO,
integer, intent(in)  N,
real(r8ki), dimension( : ), intent(in)  AP,
real(r8ki), dimension( :,: ), intent(out)  A,
integer, intent(in)  LDA,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Unpack a by-column-packed array into a 2D matrix format See documentation in DTPTTR/STPTTR source code.

Parameters
[in]uplo= 'U': A is an upper triangular matrix; 'L': A is a lower triangular matrix
[in]nThe order of matrix A and AP.
[in]ldaThe leading dimension of the matrix A. LDA ? max(1,N)
[out]errstatError level
[out]errmsgMessage describing error
[in]apPacked array
[out]aUnpacked array : Note AP(1)=A(1,1); AP(2)=A(1,2); AP(3)=A(2,2); AP(4)=A(1,3) etc. by column, upper triang

◆ lapack_sgbsv()

subroutine nwtc_lapack::lapack_sgbsv ( integer, intent(in)  N,
integer, intent(in)  KL,
integer, intent(in)  KU,
integer, intent(in)  NRHS,
real(siki), dimension( :, : ), intent(inout)  AB,
integer, dimension( : ), intent(out)  IPIV,
real(siki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general banded solve: Computes the solution to system of linear equations A * X = B for GB (general, banded) matrices.

use LAPACK_GBSV (nwtc_lapack::lapack_gbsv) instead of this specific function.

Parameters
[in]klThe number of subdiagonals within the band of A. KL >= 0.
[in]kuThe number of superdiagonals within the band of A. KU >= 0.
[in]nThe number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]abOn entry, the matrix A in band storage, in rows KL+1 to 2*KL+KU+1; rows 1 to KL of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows: AB(KL+KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+KL) On exit, details of the factorization: U is stored as an upper triangular band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and the multipliers used during the factorization are stored in rows KL+KU+2 to 2*KL+KU+1.
[in,out]bOn entry, the N-by-NRHS right hand side matrix B. On exit, if INFO = 0, the N-by-NRHS solution matrix X.
[out]ipivThe pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sgemm()

subroutine nwtc_lapack::lapack_sgemm ( character(1), intent(in)  TRANSA,
character(1), intent(in)  TRANSB,
real(siki), intent(in)  ALPHA,
real(siki), dimension( :, : ), intent(in)  A,
real(siki), dimension( :, : ), intent(in)  B,
real(siki), intent(in)  BETA,
real(siki), dimension( :, : ), intent(inout)  C,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general matrix multiply: computes C = alpha*op( A )*op( B ) + beta*C where op(x) = x or op(x) = x**T for matrices A, B, and C use LAPACK_GEMM (nwtc_lapack::lapack_gemm) instead of this specific function.

Parameters
[in]transaOn entry, TRANSA specifies the form of op( A ) to be used in the matrix multiplication as follows: TRANSA = 'N' or 'n', op( A ) = A. TRANSA = 'T' or 't', op( A ) = A**T.
[in]transbOn entry, TRANSB specifies the form of op( A ) to be used in the matrix multiplication as follows: TRANSB = 'N' or 'n', op( B ) = B. TRANSB = 'T' or 't', op( B ) = B**T.
[in]alphaOn entry, ALPHA specifies the scalar alpha.
[in]betaOn entry, BETA specifies the scalar beta. When BETA is supplied as zero then C need not be set on input.
[in]aMatrix A
[in]bMatrix B
[in,out]cMatrix C: Before entry, C must contain the matrix C, except when beta is zero, in which case C need not be set on entry. On exit, the array C is overwritten by the m by n matrix ( alpha*op( A )*op( B ) + beta*C ).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sgesv()

subroutine nwtc_lapack::lapack_sgesv ( integer, intent(in)  N,
real(siki), dimension( :, : ), intent(inout)  A,
integer, dimension( : ), intent(out)  IPIV,
real(siki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general solve: Computes the solution to system of linear equations A * X = B for GE matrices.

use LAPACK_GESV (nwtc_lapack::lapack_gesv) instead of this specific function.

Parameters
[in]nThe number of linear equations, i.e., the order of the matrix A. N >= 0.
[in,out]aOn entry, the N-by-N coefficient matrix A. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[in,out]bOn entry, the N-by-NRHS matrix of right hand side matrix B. On exit, if INFO = 0, the N-by-NRHS solution matrix X.
[out]ipivThe pivot indices that define the permutation matrix P; row i of the matrix was interchanged with row IPIV(i).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sgesvd()

subroutine nwtc_lapack::lapack_sgesvd ( character(1), intent(in)  JOBU,
character(1), intent(in)  JOBVT,
integer, intent(in)  M,
integer, intent(in)  N,
real(siki), dimension( :, : ), intent(inout)  A,
real(siki), dimension( : ), intent(out)  S,
real(siki), dimension( :, : ), intent(out)  U,
real(siki), dimension( :, : ), intent(out)  VT,
real(siki), dimension( : ), intent(inout)  WORK,
integer, intent(in)  LWORK,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute singular value decomposition (SVD) for a general matrix, A.

use LAPACK_SGESVD (nwtc_lapack::lapack_sgesvd) instead of this specific function.

Parameters
[in]jobu'A': all M columns of U are returned in array U; 'S': the first min(m,n) columns of U (the left singular vectors) are returned in the array U; 'O': the first min(m,n) columns of U (the left singular vectors) are overwritten on the array A; 'N': no columns of U (no left singular vectors) are computed.
[in]jobvt'A': all N rows of V^T are returned in the array VT; 'S': the first min(m,n) rows of V^T (the right singular vectors) are returned in the array VT; 'O': the first min(m,n) rows of V^T (the right singular vectors) are overwritten on the array A; 'N': no rows of V**T (no right singular vectors) are computed.
[in]mThe number of rows of the input matrix A. M >= 0.
[in]nThe number of columns of the input matrix A. N >= 0.
[in,out]aA is SINGLE PRECISION array, dimension (LDA,N) On entry, the M-by-N matrix A. On exit, if JOBU = 'O', A is overwritten with the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBVT = 'O', A is overwritten with the first min(m,n) rows of V**T (the right singular vectors, stored rowwise); if JOBU .ne. 'O' and JOBVT .ne. 'O', the contents of A are destroyed.
[out]sS is SINGLE PRECISION array, dimension (min(M,N)) The singular values of A, sorted so that S(i) >= S(i+1).
[out]uU is SINGLE PRECISION array, dimension (LDU,UCOL) (LDU,M) if JOBU = 'A' or (LDU,min(M,N)) if JOBU = 'S'. If JOBU = 'A', U contains the M-by-M orthogonal matrix U; if JOBU = 'S', U contains the first min(m,n) columns of U (the left singular vectors, stored columnwise); if JOBU = 'N' or 'O', U is not referenced.
[out]vtVT is SINGLE PRECISION array, dimension (LDVT,N) If JOBVT = 'A', VT contains the N-by-N orthogonal matrix V**T; if JOBVT = 'S', VT contains the first min(m,n) rows of V**T (the right singular vectors, stored rowwise); if JOBVT = 'N' or 'O', VT is not referenced.
[in,out]workdimension (MAX(1,LWORK)). On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]lworkThe dimension of the array WORK. LWORK >= MAX(1,5*MIN(M,N)) for the paths (see comments inside code):
  • PATH 1 (M much larger than N, JOBU='N')
  • PATH 1t (N much larger than M, JOBVT='N') LWORK >= MAX(1,3*MIN(M,N) + MAX(M,N),5*MIN(M,N)) for the other paths For good performance, LWORK should generally be larger. If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sgetrf()

subroutine nwtc_lapack::lapack_sgetrf ( integer, intent(in)  M,
integer, intent(in)  N,
real(siki), dimension( :, : ), intent(inout)  A,
integer, dimension( : ), intent(out)  IPIV,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general matrix factorization: Factor matrix into A=PLU.

use LAPACK_GETRF (nwtc_lapack::lapack_getrf) instead of this specific function.

Parameters
[in]mThe number of rows of the matrix A. M >= 0.
[in]nThe number of columns of the matrix A. N >= 0.
[in,out]aOn entry, the M-by-N matrix to be factored. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
[out]ipivThe pivot indices; for 1 <= i <= min(M,N), row i of the matrix was interchanged with row IPIV(i).
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sgetri()

subroutine nwtc_lapack::lapack_sgetri ( integer, intent(in)  N,
real(siki), dimension( :, : ), intent(inout)  A,
integer, dimension( : ), intent(in)  IPIV,
real(siki), dimension( : ), intent(inout)  WORK,
integer, intent(in)  LWORK,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute the inverse of a general matrix using the LU factorization.

use LAPACK_GETRI (nwtc_lapack::lapack_getri) instead of this specific function.

Parameters
[in]nThe order of the matrix A. N >= 0.
[in]lworkThe dimension of the array WORK. LWORK >= max(1,N). For optimal performance LWORK >= N*NB, where NB is the optimal blocksize returned by ILAENV. If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.
[in]ipivdimension (N). The pivot indices from DGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).
[in,out]aOn entry, the factors L and U from the factorization A = P*L*U as computed by SGETRF. On exit, if INFO = 0, the inverse of the original matrix A.
[in,out]workOn exit, if INFO=0, then WORK(1) returns the optimal LWORK.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sgetrs()

subroutine nwtc_lapack::lapack_sgetrs ( character(1), intent(in)  TRANS,
integer, intent(in)  N,
real(siki), dimension( :, : ), intent(in)  A,
integer, dimension( : ), intent(in)  IPIV,
real(siki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b.

use LAPACK_GETRS (nwtc_lapack::lapack_getrs) instead of this specific function.

Parameters
[in]transSpecifies the form of the system of equations: = 'N': A * X = B (No transpose) = 'T': A**T* X = B (Transpose) = 'C': A**T* X = B (Conjugate transpose = Transpose)
[in]nThe order of the matrix A. N >= 0.
[in]ipivThe pivot indices from DGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).
[in]aThe factors L and U from the factorization A = P*L*U as computed by SGETRF.
[in,out]bOn entry, the right hand side matrix B. On exit, the solution matrix X.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sgetrs1()

subroutine nwtc_lapack::lapack_sgetrs1 ( character(1), intent(in)  TRANS,
integer, intent(in)  N,
real(siki), dimension( :, : ), intent(in)  A,
integer, dimension( : ), intent(in)  IPIV,
real(siki), dimension( : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

general solve of factorized matrix: Solve system of linear equations Ax=PLUx=b.

use LAPACK_GETRS (nwtc_lapack::lapack_getrs) instead of this specific function.

Parameters
[in]transSpecifies the form of the system of equations: = 'N': A * X = B (No transpose) = 'T': A**T* X = B (Transpose) = 'C': A**T* X = B (Conjugate transpose = Transpose)
[in]nThe order of the matrix A. N >= 0.
[in]ipivThe pivot indices from DGETRF; for 1<=i<=N, row i of the matrix was interchanged with row IPIV(i).
[in]aThe factors L and U from the factorization A = P*L*U as computed by SGETRF.
[in,out]bOn entry, the right hand side matrix B. On exit, the solution matrix X.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sggev()

subroutine nwtc_lapack::lapack_sggev ( character(1), intent(in)  JOBVL,
character(1), intent(in)  JOBVR,
integer, intent(in)  N,
real(siki), dimension( :, : ), intent(inout)  A,
real(siki), dimension( :, : ), intent(inout)  B,
real(siki), dimension( : ), intent(out)  ALPHAR,
real(siki), dimension( : ), intent(out)  ALPHAI,
real(siki), dimension( : ), intent(out)  BETA,
real(siki), dimension( :, : ), intent(out)  VL,
real(siki), dimension( :, : ), intent(out)  VR,
real(siki), dimension( : ), intent(inout)  WORK,
integer, intent(in)  LWORK,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute generalized eigenvalues and/or eigenvectors for a pair of N-by-N real nonsymmetric matrices (A,B).

use LAPACK_GGEV (nwtc_lapack::lapack_ggev) instead of this specific function.

Parameters
[in]jobvl= 'N': do not compute the left generalized eigenvectors; = 'V': compute the left generalized eigenvectors.
[in]jobvr= 'N': do not compute the right generalized eigenvectors; = 'V': compute the right generalized eigenvectors.
[in]nThe order of the matrices A, B, VL, and VR. N >= 0.
[in]lworkThe dimension of the array WORK. LWORK >= max(1,8*N). For good performance, LWORK must generally be larger. If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.
[in,out]adimension (LDA, N). On entry, the matrix A in the pair (A,B). On exit, A has been overwritten.
[in,out]bdimension (LDB, N). On entry, the matrix B in the pair (A,B). On exit, B has been overwritten.
[out]alphardimension (N). See comments for variable "Beta"
[out]alphaidimension (N). See comments for variable "Beta".
[out]betaOn exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will be the generalized eigenvalues. If ALPHAI(j) is zero, then the j-th eigenvalue is real; if positive, then the j-th and (j+1)-st eigenvalues are a complex conjugate pair, with ALPHAI(j+1) negative.

Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j) may easily over- or underflow, and BETA(j) may even be zero. Thus, the user should avoid naively computing the ratio alpha/beta. However, ALPHAR and ALPHAI will be always less than and usually comparable with norm(A) in magnitude, and BETA always less than and usually comparable with norm(B).

Parameters
[out]vldimension (LDVL,N). If JOBVL = 'V', the left eigenvectors u(j) are stored one after another in the columns of VL, in the same order as their eigenvalues. If the j-th eigenvalue is real, then u(j) = VL(:,j), the j-th column of VL. If the j-th and (j+1)-th eigenvalues form a complex conjugate pair, then u(j) = VL(:,j)+i*VL(:,j+1) and u(j+1) = VL(:,j)-i*VL(:,j+1). Each eigenvector is scaled so the largest component has abs(real part)+abs(imag. part)=1. Not referenced if JOBVL = 'N'.
[out]vrdimension (LDVR,N). If JOBVR = 'V', the right eigenvectors v(j) are stored one after another in the columns of VR, in the same order as their eigenvalues. If the j-th eigenvalue is real, then v(j) = VR(:,j), the j-th column of VR. If the j-th and (j+1)-th eigenvalues form a complex conjugate pair, then v(j) = VR(:,j)+i*VR(:,j+1) and v(j+1) = VR(:,j)-i*VR(:,j+1). Each eigenvector is scaled so the largest component has abs(real part)+abs(imag. part)=1. Not referenced if JOBVR = 'N'.
[in,out]workdimension (MAX(1,LWORK)). On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[out]errstatError level
[out]errmsgMessage describing error

◆ lapack_sposv()

subroutine nwtc_lapack::lapack_sposv ( character(1), intent(in)  UPLO,
integer, intent(in)  N,
integer, intent(in)  NRHS,
real(siki), dimension( :, : ), intent(inout)  A,
real(siki), dimension( :, : ), intent(inout)  B,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute the solution to system of linear equations A * X = B for PO (positive-definite) matrices.

use LAPACK_POSV (nwtc_lapack::lapack_posv) instead of this specific function.

Parameters
[in]nThe number of linear equations, i.e., the order of the matrix A. N >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
[in,out]aOn entry, the symmetric matrix A. If UPLO = 'U', the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, if INFO = 0, the factor U or L from the Cholesky factorization A = U**T*U or A = L*L**T.
[in,out]bOn entry, the N-by-NRHS matrix of right hand side matrix B. On exit, if INFO = 0, the N-by-NRHS solution matrix X.
[out]errstatError level
[out]errmsgMessage describing error
[in]uplo'U': Upper triangle of A is stored; 'L': Lower triangle of A is stored.

◆ lapack_spptrf()

subroutine nwtc_lapack::lapack_spptrf ( character(1), intent(in)  UPLO,
integer, intent(in)  N,
real(siki), dimension( : ), intent(inout)  AP,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Compute the Cholesky factorization of a real symmetric positive definite matrix A stored in packed format.

use LAPACK_PPTRF (nwtc_lapack::lapack_pptrf) instead of this specific function.

Parameters
[in]nThe order of the matrix A. N >= 0.
[in,out]apAP is REAL array, dimension (N*(N+1)/2) On entry, the upper or lower triangle of the symmetric matrix A, packed columnwise in a linear array. The j-th column of A is stored in the array AP as follows: if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j; if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n. See LAPACK_DPPTRF for further details. On exit, if INFO = 0, the triangular factor U or L from the Cholesky factorization A = U**T*U or A = L*L**T, in the same storage format as A.
[out]errstatError level
[out]errmsgMessage describing error
[in]uplo'U': Upper triangle of A is stored; 'L': Lower triangle of A is stored.

◆ lapack_stpttr()

subroutine nwtc_lapack::lapack_stpttr ( character(1), intent(in)  UPLO,
integer, intent(in)  N,
real(siki), dimension( : ), intent(in)  AP,
real(siki), dimension( :,: ), intent(out)  A,
integer, intent(in)  LDA,
integer(intki), intent(out)  ErrStat,
character(*), intent(out)  ErrMsg 
)

Unpack a by-column-packed array into a 2D matrix format.

Parameters
[in]uplo= 'U': A is an upper triangular matrix; 'L': A is a lower triangular matrix
[in]nThe order of matrix A and AP.
[in]ldaThe leading dimension of the matrix A. LDA ? max(1,N)
[out]errstatError level
[out]errmsgMessage describing error
[in]apPacked array
[out]aUnpacked array : Note AP(1)=A(1,1); AP(2)=A(1,2); AP(3)=A(2,2); AP(4)=A(1,3) etc. by column, upper triang