MODULE Lapack_interfaces

Module Lapack_interfaces contains/exports generic interfaces for selected routines/drivers available in the LAPACK library (LAPACK) for use within the framework of STATPACK. Note, however, that contrary to the BLAS routines, which are used inside the STATPACK library if the cpp macro _BLAS is activated at compilation of STATPACK, LAPACK routines are not presently used inside STATPACK and the cpp macro _LAPACK is not defined in STATPACK and will have no effect at compilation.

However, use of the interface blocks exported by module Lapack_interfaces unsures that calls to LAPACK routines are correct, when used with STATPACK. Generic interfaces are presently provided for the following LAPACK routines and drivers:

  • Tridiagonal reduction of a real symmetric or complex hermitian matrix:

    sytrd()

    Generic interface for SSYTRD, DSYTRD, CHETRD and ZHETRD subroutines (decomposition)

    Examples: ex1_lapack_sytrd.F90 ex2_lapack_sytrd.F90

    orgtr()

    Generic interface for SORGTR, CORGTR, CUNGTR and ZUNGTR subroutines (generation of orthogonal matrix)

    Examples: ex1_lapack_orgtr.F90

    ormtr()

    Generic interface for SORMTR, CORMTR, CUNMTR and ZUNMTR subroutines (multiplication by orthogonal matrix)

    Examples: ex1_lapack_ormtr.F90 ex2_lapack_ormtr.F90

  • Eigenvalues and eigenvectors decomposition of a real symmetric or complex hermitian matrix:

    syev()

    Generic interface for SSYEV, DSYEV, CHEEV and ZHEEV subroutines (implicit QR/QL method)

    Examples: ex1_lapack_syev.F90 ex2_lapack_syev.F90

    syevd()

    Generic interface for SSYEVD, DSYEVD, CHEEVD and ZHEEVD subroutines (divide and conquer method)

    Examples: ex1_lapack_syevd.F90 ex2_lapack_syevd.F90

    syevr()

    Generic interface for SSYEVR, DSYEVR, CHEEVR and ZHEEVR subroutines (MRRR method)

    Examples: ex1_lapack_syevr.F90 ex2_lapack_syevr.F90 ex3_lapack_syevr.F90

    syevx()

    Generic interface for SSYEVX, DSYEVX, CHEEVX and ZHEEVX subroutines (bisection and inverse iteration)

    Examples: ex1_lapack_syevx.F90 ex2_lapack_syevx.F90 ex3_lapack_syevx.F90

  • Eigenvalues and eigenvectors decomposition of a real symmetric or complex hermitian matrix in packed storage:

    spev()

    Generic interface for SSPEV, DSPEV, CHPEV and ZHPEV subroutines (implicit QR/QL method)

    Examples: ex1_lapack_spev.F90

    spevd()

    Generic interface for SSPEVD, DSPEVD, CHPEVD and ZHPEVD subroutines (divide and conquer method)

    Examples: ex1_lapack_spevd.F90

    spevx()

    Generic interface for SSPEVX, DSPEVX, CHPEVX and ZHPEVX subroutines (bisection and inverse iteration)

    Examples: ex1_lapack_spevx.F90 ex3_lapack_spevx.F90

  • Eigenvalues and eigenvectors decomposition of a real symmetric tridiagonal matrix:

    steqr()

    Generic interface for SSTEQR, DSTEQR, CSTEQR and ZSTEQR subroutines (implicit QR/QL method)

    stedc()

    Generic interface for SSTEDC, DSTEDC, CSTEDC and ZSTEDC subroutines (divide and conquer method)

    stemr()

    Generic interface for SSTEMR, DSTEMR, CSTEMR and ZSTEMR subroutines (MRRR method)

    Examples: ex1_lapack_stemr.F90 ex2_lapack_stemr.F90 ex3_lapack_stemr.F90

    stevx()

    Generic interface for SSTEVX and DSTEVX subroutines (partial or full spectrum by bisection and inverse iteration)

    stev()

    Generic interface for SSTEV and DSTEV subroutines (eigenvalues by the Pal-Walker-Kahan variant of the QL or QR algorithm or eigenvalues/eigenvectors by implicit QR/QL method)

    stevd()

    Generic interface for SSTEVD and DSTEVD subroutines (eigenvalues by the Pal-Walker-Kahan variant of the QL or QR algorithm or eigenvalues/eigenvectors by divide and conquer method)

    stevr()

    Generic interface for SSTEVR and DSTEVR subroutines (full spectrum by MRRR method and partial spectrum by bisection and inverse iteration)

  • Eigenvalues and eigenvectors of a real generalized symmetric-definite or complex generalized hermitian-definite problem:

    sygv()

    Generic interface for SSYGV, DSYGV, CHEGV and ZHEGV subroutines (implicit QR/QL method)

    sygvd()

    Generic interface for SSYGVD, DSYGVD, CHEGVD and ZHEGVD subroutines (divide and conquer method)

    sygvx()

    Generic interface for SSYGVX, DSYGVX, CHEGVX and ZHEGVX subroutines (bisection and inverse iteration)

  • Eigenvalues and eigenvectors of a real or complex general matrix:

    geev()

    Generic interface for SGEEV, DGEEV, CGEEV and ZGEEV subroutines (QR method)

    geevx()

    Generic interface for SGEEVX, DGEEVX, CGEEVX and ZGEEVX subroutines (QR method with balancing)

  • Bidiagonal reduction of a real or complex general matrix:

    gebrd()

    Generic interface for SGEBRD, DGEBRD, CGEBRD and ZGEBRD subroutines (decomposition)

    Examples: ex1_lapack_gebrd.F90 ex2_lapack_gebrd.F90

    orgbr()

    Generic interface for SORGBR, DORGBR, CUNGBR and ZUNGBR subroutines (generation of orthogonal matrices)

    Examples: ex1_lapack_orgbr.F90

    ormbr()

    Generic interface for SORMBR, CORMBR, CUNMBR and ZUNMBR subroutines (multiplication by orthogonal matrices)

    Examples: ex1_lapack_ormbr.F90 ex2_lapack_ormbr.F90

  • Singular Value Decomposition (SVD) of a real or complex general matrix:

    gesvd()

    Generic interface for SGESVD, DGESVD, CGESVD and ZGESVD subroutines (implicit QR method)

    Examples: ex1_lapack_gesvd.F90 ex2_lapack_gesvd.F90 ex3_lapack_gesvd.F90

    gesdd()

    Generic interface for SGESDD, DGESDD, CGESDD and ZGESDD subroutines (divide and conquer method)

    Examples: ex1_lapack_gesdd.F90 ex2_lapack_gesdd.F90 ex3_lapack_gesdd.F90

    gesvdx()

    Generic interface for SGESVDX, DGESVDX, CGESVDX and ZGESVDX SGESVD subroutines (bisection/inverse iteration, including partial SVD)

    Examples: ex1_lapack_gesvdx.F90 ex2_lapack_gesvdx.F90 ex3_lapack_gesvdx.F90

  • Singular Value Decomposition (SVD) of a real bidiagonal matrix:

    bdsqr()

    Generic interface for SBDSQR, DBDSQR, CBDSQR and ZBDSQR subroutines (implicit QR method)

    bdsdc()

    Generic interface for SBDSDC and DBDSDC subroutines (divide and conquer method)

    bdsvdx()

    Generic interface for SBDSVDX and DBDSVDX subroutines (bisection/inverse iteration, including partial SVD)

  • Solution of a real or complex system of linear equations with a general matrix and several right hand side vectors:

    gesv()

    Generic interface for SGESV, DGESV, CGESV and ZGESV subroutines (LU decomposition)

    Examples: ex1_lapack_gesv.F90 ex2_lapack_gesv.F90

  • Solution of a real or complex system of linear equations with a symmetric matrix and several right hand side vectors:

    sysv()

    Generic interface for SSYSV, DSYSV, CSYSV and ZSYSV subroutines (diagonal pivoting method)

    Examples: ex1_lapack_sysv.F90 ex2_lapack_sysv.F90

  • Solution of a real or complex system of linear equations with a symmetric or hermitian positive definite matrix and several right hand side vectors:

    posv()

    Generic interface for SPOSV, DPOSV, CPOSV and ZPOSV subroutines (Cholesky decomposition)

    Examples: ex1_lapack_posv.F90 ex2_lapack_posv.F90

  • Minimum-norm solution of a real or complex linear least square problem with several right hand side vectors:

    gelss()

    Generic interface for SGELSS, DGELSS, CGELSS and ZGELSS subroutines (SVD via implicit QR method)

    Examples: ex1_lapack_gelss.F90 ex2_lapack_gelss.F90

    gelsd()

    Generic interface for SGELSD, DGELSD, CGELSD and ZGELSD subroutines (SVD via divide and conquer method)

    Examples: ex1_lapack_gelsd.F90 ex2_lapack_gelsd.F90

    gelsy()

    Generic interface for SGELSY, DGELSY, CGELSY and ZGELSY subroutines (complete orthogonal factorization)

    Examples: ex1_lapack_gelsy.F90 ex2_lapack_gelsy.F90

  • Solution of a real or complex, overdetermined or underdetermined, linear system with a coefficient matrix of full rank and several right hand side vectors:

    gels()

    Generic interface for SGELS, DGELS, CGELS and ZGELS subroutines (QR/QL method)

    Examples: ex1_lapack_gels.F90 ex2_lapack_gels.F90

  • Computation of an unblocked QR factorization of a real or complex matrix:

    geqr2()

    Generic interface for SGEQR2, DGEQR2, CGEQR2 and ZGEQR2 subroutines

  • Computation of a blocked QR factorization of a real or complex matrix:

    geqrf()

    Generic interface for SGEQRF, DGEQRF, CGEQRF and ZGEQRF subroutines

  • Computation of a blocked QR factorization with column pivoting of a real or complex matrix using level 3 BLAS:

    geqp3()

    Generic interface for SGEQP3, DGEQP3, CGEQP3 and ZGEQP3 subroutines

  • Computation of a blocked QR factorization of a “triangular-pentagonal” real or complex matrix, which is composed of a triangular block A and pentagonal block B:

    tpqrt()

    Generic interface for STPQRT, DTPQRT, CTPQRT, ZTPQRT subroutines

  • Application of an orthogonal matrix Q obtained from a “triangular-pentagonal” block reflector to a general real or complex matrix:

    tpmqrt()

    Generic interface for STPMQRT, DTPMQRT, CTPMQRT, ZTPMQRT subroutines

Consult the official LAPACK site at LAPACK, the hyper-text documentation at LAPACK documentation or the nice summary available at http://www.icl.utk.edu/~mgates3/docs/lapack.html for the definition/documentation of the LAPACK routines.

See the FORTRAN programs ex1_lapack_gesdd.F90 and ex2_lapack_gesdd.F90 for working examples of using the STATPACK Fortran 95 generic interface gesdd() (defined in the Lapack_interfaces module) for subroutines xGESDD() (where x can be S, D, C or Z) available in the LAPACK library for performing a Singular Value Decomposition (SVD) of a real/complex matrix of kind stnd by the divide and conquer method, inside the framework of STATPACK.

Since the LAPACK library provides routines only for single and double precision real/complex data, the interface blocks defined in the module Lapack_interfaces will work obviously only if the real/complex kind type stnd defined in module Select_Parameters is equivalent to single or double precision real/complex data.

Finally, note that you can add at your convenience interface blocks for other LAPACK routines in module Lapack_interfaces, which is in the file Module_Lapack_Interfaces.F90. Here, is an example of the generic interface syevd() for the SSYEVD, DSYEVD, CHEEVD and ZHEEVD subroutines available in LAPACK, which can be used as a model for creating a generic interface for other LAPACK subroutines:

!
!   Generic interface for SSYEVD, DSYEVD, CHEEVD and ZHEEVD subroutines in LAPACK
!
interface syevd
!
    subroutine ssyevd( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK, INFO )
!     ..
!     .. Scalar Arguments ..
       CHARACTER        JOBZ, UPLO
       INTEGER          INFO, LDA, LIWORK, LWORK, N
!     ..
!     .. Array Arguments ..
       INTEGER          IWORK( * )
       REAL             A( LDA, * ), W( * ), WORK( * )
    end subroutine
!
    subroutine dsyevd( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK, INFO )
!     ..
!     .. Scalar Arguments ..
       CHARACTER        JOBZ, UPLO
       INTEGER          INFO, LDA, LIWORK, LWORK, N
!     ..
!     .. Array Arguments ..
       INTEGER          IWORK( * )
       DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * )
    end subroutine
!
    subroutine cheevd( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO )
!     ..
!     .. Scalar Arguments ..
       CHARACTER        JOBZ, UPLO
       INTEGER          INFO, LDA, LIWORK, LRWORK, LWORK, N
!     ..
!     .. Array Arguments ..
       INTEGER          IWORK( * )
       REAL             W( * ), RWORK( * )
       COMPLEX          A( LDA, * ), WORK( * )
    end subroutine
!
    subroutine zheevd( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO )
!     ..
!     .. Scalar Arguments ..
       CHARACTER        JOBZ, UPLO
       INTEGER          INFO, LDA, LIWORK, LRWORK, LWORK, N
!     ..
!     .. Array Arguments ..
       INTEGER          IWORK( * )
       DOUBLE PRECISION W( * ), RWORK( * )
       COMPLEX*16       A( LDA, * ), WORK( * )
    end subroutine
!
end interface
Flag Counter