Modelica.Math.Matrices.LAPACK

Interface to LAPACK library (should usually not directly be used but only indirectly via Modelica.Math.Matrices)

Information


This package contains external Modelica functions as interface to the LAPACK library (http://www.netlib.org/lapack) that provides FORTRAN subroutines to solve linear algebra tasks. Usually, these functions are not directly called, but only via the much more convenient interface of Modelica.Math.Matrices. The documentation of the LAPACK functions is a copy of the original FORTRAN code.

The details of LAPACK are described in:

Anderson E., Bai Z., Bischof C., Blackford S., Demmel J., Dongarra J., Du Croz J., Greenbaum A., Hammarling S., McKenney A., and Sorensen D.:
Lapack Users' Guide. Third Edition, SIAM, 1999.

This package contains a direct interface to the LAPACK subroutines

Extends from Modelica.Icons.Library (Icon for library).

Package Content

NameDescription
Modelica.Math.Matrices.LAPACK.dgeev dgeev Compute eigenvalues and (right) eigenvectors for real nonsymmetrix matrix A
Modelica.Math.Matrices.LAPACK.dgeev_eigenValues dgeev_eigenValues Compute eigenvalues for real nonsymmetrix matrix A
Modelica.Math.Matrices.LAPACK.dgegv dgegv Compute generalized eigenvalues and eigenvectors for a (A,B) system
Modelica.Math.Matrices.LAPACK.dgels_vec dgels_vec Solves overdetermined or underdetermined real linear equations A*x=b with a b vector
Modelica.Math.Matrices.LAPACK.dgelsx_vec dgelsx_vec Computes the minimum-norm solution to a real linear least squares problem with rank deficient A
Modelica.Math.Matrices.LAPACK.dgesv dgesv Solve real system of linear equations A*X=B with a B matrix
Modelica.Math.Matrices.LAPACK.dgesv_vec dgesv_vec Solve real system of linear equations A*x=b with a b vector
Modelica.Math.Matrices.LAPACK.dgglse_vec dgglse_vec Solve a linear equality constrained least squares problem
Modelica.Math.Matrices.LAPACK.dgtsv dgtsv Solve real system of linear equations A*X=B with B matrix and tridiagonal A
Modelica.Math.Matrices.LAPACK.dgtsv_vec dgtsv_vec Solve real system of linear equations A*x=b with b vector and tridiagonal A
Modelica.Math.Matrices.LAPACK.dgbsv dgbsv Solve real system of linear equations A*X=B with a B matrix
Modelica.Math.Matrices.LAPACK.dgbsv_vec dgbsv_vec Solve real system of linear equations A*x=b with a b vector
Modelica.Math.Matrices.LAPACK.dgesvd dgesvd Determine singular value decomposition
Modelica.Math.Matrices.LAPACK.dgesvd_sigma dgesvd_sigma Determine singular values
Modelica.Math.Matrices.LAPACK.dgetrf dgetrf Compute LU factorization of square or rectangular matrix A (A = P*L*U)
Modelica.Math.Matrices.LAPACK.dgetrs dgetrs Solves a system of linear equations with the LU decomposition from dgetrf(..)
Modelica.Math.Matrices.LAPACK.dgetrs_vec dgetrs_vec Solves a system of linear equations with the LU decomposition from dgetrf(..)
Modelica.Math.Matrices.LAPACK.dgetri dgetri Computes the inverse of a matrix using the LU factorization from dgetrf(..)
Modelica.Math.Matrices.LAPACK.dgeqpf dgeqpf Compute QR factorization of square or rectangular matrix A with column pivoting (A(:,p) = Q*R)
Modelica.Math.Matrices.LAPACK.dorgqr dorgqr Generates a Real orthogonal matrix Q which is defined as the product of elementary reflectors as returned from dgeqpf


Modelica.Math.Matrices.LAPACK.dgeev Modelica.Math.Matrices.LAPACK.dgeev

Compute eigenvalues and (right) eigenvectors for real nonsymmetrix matrix A

Information

Lapack documentation
    Purpose
    =======
    DGEEV computes for an N-by-N real nonsymmetric matrix A, the
    eigenvalues and, optionally, the left and/or right eigenvectors.
    The right eigenvector v(j) of A satisfies
                     A * v(j) = lambda(j) * v(j)
    where lambda(j) is its eigenvalue.
    The left eigenvector u(j) of A satisfies
                  u(j)**H * A = lambda(j) * u(j)**H
    where u(j)**H denotes the conjugate transpose of u(j).
    The computed eigenvectors are normalized to have Euclidean norm
    equal to 1 and largest component real.
    Arguments
    =========
    JOBVL   (input) CHARACTER*1
            = 'N': left eigenvectors of A are not computed;
            = 'V': left eigenvectors of A are computed.
    JOBVR   (input) CHARACTER*1
            = 'N': right eigenvectors of A are not computed;
            = 'V': right eigenvectors of A are computed.
    N       (input) INTEGER
            The order of the matrix A. N >= 0.
    A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
            On entry, the N-by-N matrix A.
            On exit, A has been overwritten.
    LDA     (input) INTEGER
            The leading dimension of the array A.  LDA >= max(1,N).
    WR      (output) DOUBLE PRECISION array, dimension (N)
    WI      (output) DOUBLE PRECISION array, dimension (N)
            WR and WI contain the real and imaginary parts,
            respectively, of the computed eigenvalues.  Complex
            conjugate pairs of eigenvalues appear consecutively
            with the eigenvalue having the positive imaginary part
            first.
    VL      (output) DOUBLE PRECISION array, dimension (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 JOBVL = 'N', VL is not referenced.
            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)-st 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).
    LDVL    (input) INTEGER
            The leading dimension of the array VL.  LDVL >= 1; if
            JOBVL = 'V', LDVL >= N.
    VR      (output) DOUBLE PRECISION array, dimension (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 JOBVR = 'N', VR is not referenced.
            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)-st 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).
    LDVR    (input) INTEGER
            The leading dimension of the array VR.  LDVR >= 1; if
            JOBVR = 'V', LDVR >= N.
    WORK    (workspace/output) DOUBLE PRECISION array, dimension (LWORK)

            On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
    LWORK   (input) INTEGER
            The dimension of the array WORK.  LWORK >= max(1,3*N), and
            if JOBVL = 'V' or JOBVR = 'V', LWORK >= 4*N.  For good
            performance, LWORK must generally be larger.
    INFO    (output) INTEGER
            = 0:  successful exit
            < 0:  if INFO = -i, the i-th argument had an illegal value.
            > 0:  if INFO = i, the QR algorithm failed to compute all the
                  eigenvalues, and no eigenvectors have been computed;
                  elements i+1:N of WR and WI contain eigenvalues which
                  have converged.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, size(A, 1)]  

Outputs

TypeNameDescription
RealeigenReal[size(A, 1)]Real part of eigen values
RealeigenImag[size(A, 1)]Imaginary part of eigen values
RealeigenVectors[size(A, 1), size(A, 1)]Right eigen vectors
Integerinfo 

Modelica definition

function dgeev 
  "Compute eigenvalues and (right) eigenvectors for real nonsymmetrix matrix A"

  extends Modelica.Icons.Function;
  input Real A[:, size(A, 1)];
  output Real eigenReal[size(A, 1)] "Real part of eigen values";
  output Real eigenImag[size(A, 1)] "Imaginary part of eigen values";
  output Real eigenVectors[size(A, 1), size(A, 1)] "Right eigen vectors";
  output Integer info;
protected 
  Integer n=size(A, 1);
  Integer lwork=12*n;
  Real Awork[n, n]=A;
  Real work[lwork];

external "Fortran 77" dgeev("N", "V", n, Awork, n, eigenReal, eigenImag,
    eigenVectors, n, eigenVectors, n, work, size(work, 1), info);
end dgeev;

Modelica.Math.Matrices.LAPACK.dgeev_eigenValues Modelica.Math.Matrices.LAPACK.dgeev_eigenValues

Compute eigenvalues for real nonsymmetrix matrix A

Information

Lapack documentation
    Purpose
    =======
    DGEEV computes for an N-by-N real nonsymmetric matrix A, the
    eigenvalues and, optionally, the left and/or right eigenvectors.
    The right eigenvector v(j) of A satisfies
                     A * v(j) = lambda(j) * v(j)
    where lambda(j) is its eigenvalue.
    The left eigenvector u(j) of A satisfies
                  u(j)**H * A = lambda(j) * u(j)**H
    where u(j)**H denotes the conjugate transpose of u(j).
    The computed eigenvectors are normalized to have Euclidean norm
    equal to 1 and largest component real.
    Arguments
    =========
    JOBVL   (input) CHARACTER*1
            = 'N': left eigenvectors of A are not computed;
            = 'V': left eigenvectors of A are computed.
    JOBVR   (input) CHARACTER*1
            = 'N': right eigenvectors of A are not computed;
            = 'V': right eigenvectors of A are computed.
    N       (input) INTEGER
            The order of the matrix A. N >= 0.
    A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
            On entry, the N-by-N matrix A.
            On exit, A has been overwritten.
    LDA     (input) INTEGER
            The leading dimension of the array A.  LDA >= max(1,N).
    WR      (output) DOUBLE PRECISION array, dimension (N)
    WI      (output) DOUBLE PRECISION array, dimension (N)
            WR and WI contain the real and imaginary parts,
            respectively, of the computed eigenvalues.  Complex
            conjugate pairs of eigenvalues appear consecutively
            with the eigenvalue having the positive imaginary part
            first.
    VL      (output) DOUBLE PRECISION array, dimension (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 JOBVL = 'N', VL is not referenced.
            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)-st 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).
    LDVL    (input) INTEGER
            The leading dimension of the array VL.  LDVL >= 1; if
            JOBVL = 'V', LDVL >= N.
    VR      (output) DOUBLE PRECISION array, dimension (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 JOBVR = 'N', VR is not referenced.
            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)-st 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).
    LDVR    (input) INTEGER
            The leading dimension of the array VR.  LDVR >= 1; if
            JOBVR = 'V', LDVR >= N.
    WORK    (workspace/output) DOUBLE PRECISION array, dimension (LWORK)

            On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
    LWORK   (input) INTEGER
            The dimension of the array WORK.  LWORK >= max(1,3*N), and
            if JOBVL = 'V' or JOBVR = 'V', LWORK >= 4*N.  For good
            performance, LWORK must generally be larger.
    INFO    (output) INTEGER
            = 0:  successful exit
            < 0:  if INFO = -i, the i-th argument had an illegal value.
            > 0:  if INFO = i, the QR algorithm failed to compute all the
                  eigenvalues, and no eigenvectors have been computed;
                  elements i+1:N of WR and WI contain eigenvalues which
                  have converged.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, size(A, 1)]  

Outputs

TypeNameDescription
RealEigenReal[size(A, 1)] 
RealEigenImag[size(A, 1)] 
Integerinfo 

Modelica definition

function dgeev_eigenValues 
  "Compute eigenvalues for real nonsymmetrix matrix A"

  extends Modelica.Icons.Function;
  input Real A[:, size(A, 1)];
  output Real EigenReal[size(A, 1)];
  output Real EigenImag[size(A, 1)];

  /*
      output Real Eigenvectors[size(A, 1), size(A, 1)]=zeros(size(A, 1), size(
          A, 1)); */
  output Integer info;
protected 
  Integer lwork=8*size(A, 1);
  Real Awork[size(A, 1), size(A, 1)]=A;
  Real work[lwork];
  Real EigenvectorsL[size(A, 1), size(A, 1)]=zeros(size(A, 1), size(A, 1));

  /*
    external "Fortran 77" dgeev("N", "V", size(A, 1), Awork, size(A, 1),
        EigenReal, EigenImag, EigenvectorsL, size(EigenvectorsL, 1),
        Eigenvectors, size(Eigenvectors, 1), work, size(work, 1), info)
*/
external "Fortran 77" dgeev("N", "N", size(A, 1), Awork, size(A, 1),
    EigenReal, EigenImag, EigenvectorsL, size(EigenvectorsL, 1),
    EigenvectorsL, size(EigenvectorsL, 1), work, size(work, 1), info);

end dgeev_eigenValues;

Modelica.Math.Matrices.LAPACK.dgegv Modelica.Math.Matrices.LAPACK.dgegv

Compute generalized eigenvalues and eigenvectors for a (A,B) system

Information

Purpose
=======

For a pair of N-by-N real nonsymmetric matrices A, B:
   compute the generalized eigenvalues (alphar +/- alphai*i, beta)
   compute the left and/or right generalized eigenvectors
           (VL and VR)
The second action is optional -- see the description of JOBVL and
JOBVR below.
A generalized eigenvalue for a pair of matrices (A,B) is, roughly
speaking, a scalar w or a ratio  alpha/beta = w, such that  A - w*B
is singular.  It is usually represented as the pair (alpha,beta),
as there is a reasonable interpretation for beta=0, and even for
both being zero.  A good beginning reference is the book, "Matrix
Computations", by G. Golub & C. van Loan (Johns Hopkins U. Press)
A right generalized eigenvector corresponding to a generalized
eigenvalue  w  for a pair of matrices (A,B) is a vector  r  such
that  (A - w B) r = 0 .  A left generalized eigenvector is a vector
                       H
l  such that  (A - w B) l = 0 .
Note: this routine performs "full balancing" on A and B -- see
"Further Details", below.
Arguments
=========
JOBVL   (input) CHARACTER*1
        = 'N':  do not compute the left generalized eigenvectors;
        = 'V':  compute the left generalized eigenvectors.
JOBVR   (input) CHARACTER*1
        = 'N':  do not compute the right generalized eigenvectors;
        = 'V':  compute the right generalized eigenvectors.
N       (input) INTEGER
        The number of rows and columns in the matrices A, B, VL, and
        VR.  N >= 0.
A       (input/workspace) DOUBLE PRECISION array, dimension (LDA, N)
        On entry, the first of the pair of matrices whose
        generalized eigenvalues and (optionally) generalized
        eigenvectors are to be computed.
        On exit, the contents will have been destroyed.  (For a
        description of the contents of A on exit, see "Further
        Details", below.)
LDA     (input) INTEGER
        The leading dimension of A.  LDA >= max(1,N).
B       (input/workspace) DOUBLE PRECISION array, dimension (LDB, N)
        On entry, the second of the pair of matrices whose
        generalized eigenvalues and (optionally) generalized
        eigenvectors are to be computed.
        On exit, the contents will have been destroyed.  (For a
        description of the contents of B on exit, see "Further
        Details", below.)
LDB     (input) INTEGER
        The leading dimension of B.  LDB >= max(1,N).
ALPHAR  (output) DOUBLE PRECISION array, dimension (N)
ALPHAI  (output) DOUBLE PRECISION array, dimension (N)
BETA    (output) DOUBLE PRECISION array, dimension (N)
        On 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).
VL      (output) DOUBLE PRECISION array, dimension (LDVL,N)
        If JOBVL = 'V', the left generalized eigenvectors.  (See
        "Purpose", above.)  Real eigenvectors take one column,
        complex take two columns, the first for the real part and
        the second for the imaginary part.  Complex eigenvectors
        correspond to an eigenvalue with positive imaginary part.
        Each eigenvector will be scaled so the largest component
        will have abs(real part) + abs(imag. part) = 1, *except*
        that for eigenvalues with alpha=beta=0, a zero vector will
        be returned as the corresponding eigenvector.
        Not referenced if JOBVL = 'N'.
LDVL    (input) INTEGER
        The leading dimension of the matrix VL. LDVL >= 1, and
        if JOBVL = 'V', LDVL >= N.
VR      (output) DOUBLE PRECISION array, dimension (LDVR,N)
        If JOBVL = 'V', the right generalized eigenvectors.  (See
        "Purpose", above.)  Real eigenvectors take one column,
        complex take two columns, the first for the real part and
        the second for the imaginary part.  Complex eigenvectors
        correspond to an eigenvalue with positive imaginary part.
        Each eigenvector will be scaled so the largest component
        will have abs(real part) + abs(imag. part) = 1, *except*
        that for eigenvalues with alpha=beta=0, a zero vector will
        be returned as the corresponding eigenvector.
        Not referenced if JOBVR = 'N'.
LDVR    (input) INTEGER
        The leading dimension of the matrix VR. LDVR >= 1, and
        if JOBVR = 'V', LDVR >= N.
WORK    (workspace/output) DOUBLE PRECISION array, dimension (LWORK)
        On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
LWORK   (input) INTEGER
        The dimension of the array WORK.  LWORK >= max(1,8*N).
        For good performance, LWORK must generally be larger.
        To compute the optimal value of LWORK, call ILAENV to get
        blocksizes (for DGEQRF, DORMQR, and DORGQR.)  Then compute:
        NB  -- MAX of the blocksizes for DGEQRF, DORMQR, and DORGQR;
        The optimal LWORK is:
            2*N + MAX( 6*N, N*(NB+1) ).
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument had an illegal value.
        = 1,...,N:
              The QZ iteration failed.  No eigenvectors have been
              calculated, but ALPHAR(j), ALPHAI(j), and BETA(j)
              should be correct for j=INFO+1,...,N.
        > N:  errors that usually indicate LAPACK problems:
              =N+1: error return from DGGBAL
              =N+2: error return from DGEQRF
              =N+3: error return from DORMQR
              =N+4: error return from DORGQR
              =N+5: error return from DGGHRD
              =N+6: error return from DHGEQZ (other than failed
                                              iteration)
              =N+7: error return from DTGEVC
              =N+8: error return from DGGBAK (computing VL)
              =N+9: error return from DGGBAK (computing VR)
              =N+10: error return from DLASCL (various calls)
Further Details
===============
Balancing
---------
This driver calls DGGBAL to both permute and scale rows and columns
of A and B.  The permutations PL and PR are chosen so that PL*A*PR
and PL*B*R will be upper triangular except for the diagonal blocks
A(i:j,i:j) and B(i:j,i:j), with i and j as close together as
possible.  The diagonal scaling matrices DL and DR are chosen so
that the pair  DL*PL*A*PR*DR, DL*PL*B*PR*DR have entries close to
one (except for the entries that start out zero.)
After the eigenvalues and eigenvectors of the balanced matrices
have been computed, DGGBAK transforms the eigenvectors back to what
they would have been (in perfect arithmetic) if they had not been
balanced.
Contents of A and B on Exit
-------- -- - --- - -- ----
If any eigenvectors are computed (either JOBVL='V' or JOBVR='V' or
both), then on exit the arrays A and B will contain the real Schur
form[*] of the "balanced" versions of A and B.  If no eigenvectors
are computed, then only the diagonal blocks will be correct.
[*] See DHGEQZ, DGEGS, or read the book "Matrix Computations",

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, size(A, 1)]  
RealB[size(A, 1), size(A, 1)]  

Outputs

TypeNameDescription
RealalphaReal[size(A, 1)]Real part of alpha (eigenvalue=(alphaReal+i*alphaImag)/beta)
RealalphaImag[size(A, 1)]Imaginary part of alpha
Realbeta[size(A, 1)]Denominator of eigenvalue
Integerinfo 

Modelica definition

function dgegv 
  "Compute generalized eigenvalues and eigenvectors for a (A,B) system"
  extends Modelica.Icons.Function;
  input Real A[:, size(A, 1)];
  input Real B[size(A,1), size(A, 1)];
  output Real alphaReal[size(A, 1)] 
    "Real part of alpha (eigenvalue=(alphaReal+i*alphaImag)/beta)";
  output Real alphaImag[size(A, 1)] "Imaginary part of alpha";
  output Real beta[size(A,1)] "Denominator of eigenvalue";
  output Integer info;
protected 
  Integer n=size(A, 1);
  Integer lwork=12*n;
  Real Awork[n, n]=A;
  Real Bwork[n, n]=B;
  Real work[lwork];
  Real dummy1[1,1];
  Real dummy2[1,1];

  external "Fortran 77" dgegv("N", "N", n, Awork, n, Bwork, n, alphaReal, alphaImag, beta,
             dummy1, 1, dummy2, 1, work, size(work, 1), info);
end dgegv;

Modelica.Math.Matrices.LAPACK.dgels_vec Modelica.Math.Matrices.LAPACK.dgels_vec

Solves overdetermined or underdetermined real linear equations A*x=b with a b vector

Information

Lapack documentation
  Purpose
  =======

  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 undetermined 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.

  Arguments
  =========

  TRANS   (input) CHARACTER
          = 'N': the linear system involves A;
          = 'T': the linear system involves A**T.

  M       (input) INTEGER
          The number of rows of the matrix A.  M >= 0.

  N       (input) INTEGER
          The number of columns of the matrix A.  N >= 0.

  NRHS    (input) INTEGER
          The number of right hand sides, i.e., the number of
          columns of the matrices B and X. NRHS >=0.

  A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
          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     (input) INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).

  B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
          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, 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     (input) INTEGER
          The leading dimension of the array B. LDB >= MAX(1,M,N).

  WORK    (workspace) DOUBLE PRECISION array, dimension (LWORK)
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.

  LWORK   (input) INTEGER
          The dimension of the array WORK.
          LWORK >= min(M,N) + MAX(1,M,N,NRHS).
          For optimal performance,
          LWORK >= min(M,N) + MAX(1,M,N,NRHS) * NB
          where NB is the optimum block size.

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

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, :]  
Realb[size(A, 1)]  

Outputs

TypeNameDescription
Realx[nx]solution is in first size(A,2) rows
Integerinfo 

Modelica definition

function dgels_vec 
  "Solves overdetermined or underdetermined real linear equations A*x=b with a b vector"

  extends Modelica.Icons.Function;
  input Real A[:, :];
  input Real b[size(A,1)];
  output Real x[nx]= cat(1,b,zeros(nx-nrow)) 
    "solution is in first size(A,2) rows";
  output Integer info;
protected 
  Integer nrow=size(A,1);
  Integer ncol=size(A,2);
  Integer nx=max(nrow,ncol);
  Integer lwork=min(nrow,ncol) + nx;
  Real work[lwork];
  Real Awork[nrow,ncol]=A;
  external "FORTRAN 77" dgels("N", nrow, ncol, 1, Awork, nrow, x,
                              nx, work, lwork, info);


end dgels_vec;

Modelica.Math.Matrices.LAPACK.dgelsx_vec Modelica.Math.Matrices.LAPACK.dgelsx_vec

Computes the minimum-norm solution to a real linear least squares problem with rank deficient A

Information

Lapack documentation
  Purpose
  =======

  DGELSX computes the minimum-norm solution to a real linear least
  squares problem:
      minimize || A * X - B ||
  using a complete orthogonal factorization 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 routine first computes a QR factorization with column pivoting:
      A * P = Q * [ R11 R12 ]
                  [  0  R22 ]
  with R11 defined as the largest leading submatrix whose estimated
  condition number is less than 1/RCOND.  The order of R11, RANK,
  is the effective rank of A.

  Then, R22 is considered to be negligible, and R12 is annihilated
  by orthogonal transformations from the right, arriving at the
  complete orthogonal factorization:
     A * P = Q * [ T11 0 ] * Z
                 [  0  0 ]
  The minimum-norm solution is then
     X = P * Z' [ inv(T11)*Q1'*B ]
                [        0       ]
  where Q1 consists of the first RANK columns of Q.

  Arguments
  =========

  M       (input) INTEGER
          The number of rows of the matrix A.  M >= 0.

  N       (input) INTEGER
          The number of columns of the matrix A.  N >= 0.

  NRHS    (input) INTEGER
          The number of right hand sides, i.e., the number of
          columns of matrices B and X. NRHS >= 0.

  A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the M-by-N matrix A.
          On exit, A has been overwritten by details of its
          complete orthogonal factorization.

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

  B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
          On entry, the M-by-NRHS right hand side matrix B.
          On exit, the N-by-NRHS solution matrix X.
          If m >= n and RANK = n, the residual sum-of-squares for
          the solution in the i-th column is given by the sum of
          squares of elements N+1:M in that column.

  LDB     (input) INTEGER
          The leading dimension of the array B. LDB >= max(1,M,N).

  JPVT    (input/output) INTEGER array, dimension (N)
          On entry, if JPVT(i) .ne. 0, the i-th column of A is an
          initial column, otherwise it is a free column.  Before
          the QR factorization of A, all initial columns are
          permuted to the leading positions; only the remaining
          free columns are moved as a result of column pivoting
          during the factorization.
          On exit, if JPVT(i) = k, then the i-th column of A*P
          was the k-th column of A.

  RCOND   (input) DOUBLE PRECISION
          RCOND is used to determine the effective rank of A, which
          is defined as the order of the largest leading triangular
          submatrix R11 in the QR factorization with pivoting of A,
          whose estimated condition number < 1/RCOND.

  RANK    (output) INTEGER
          The effective rank of A, i.e., the order of the submatrix
          R11.  This is the same as the order of the submatrix T11
          in the complete orthogonal factorization of A.

  WORK    (workspace) DOUBLE PRECISION array, dimension
                      (max( min(M,N)+3*N, 2*min(M,N)+NRHS )),

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

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, :]  
Realb[size(A, 1)]  
Realrcond0.0Reciprocal condition number to estimate rank

Outputs

TypeNameDescription
Realx[max(nrow, ncol)]solution is in first size(A,2) rows
Integerinfo 
IntegerrankEffective rank of A

Modelica definition

function dgelsx_vec 
  "Computes the minimum-norm solution to a real linear least squares problem with rank deficient A"

  extends Modelica.Icons.Function;
  input Real A[:, :];
  input Real b[size(A,1)];
  input Real rcond=0.0 "Reciprocal condition number to estimate rank";
  output Real x[max(nrow,ncol)]= cat(1,b,zeros(max(nrow,ncol)-nrow)) 
    "solution is in first size(A,2) rows";
  output Integer info;
  output Integer rank "Effective rank of A";
protected 
  Integer nrow=size(A,1);
  Integer ncol=size(A,2);
  Integer nx=max(nrow,ncol);
  Integer lwork=max( min(nrow,ncol)+3*ncol, 2*min(nrow,ncol)+1);
  Real work[lwork];
  Real Awork[nrow,ncol]=A;
  Integer jpvt[ncol]=zeros(ncol);
  external "FORTRAN 77" dgelsx(nrow, ncol, 1, Awork, nrow, x, nx, jpvt,
                              rcond, rank, work, lwork, info);


end dgelsx_vec;

Modelica.Math.Matrices.LAPACK.dgesv Modelica.Math.Matrices.LAPACK.dgesv

Solve real system of linear equations A*X=B with a B matrix

Information

Lapack documentation:
    Purpose
    =======
    DGESV computes the solution to a real system of linear equations
       A * X = B,
    where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
    The LU decomposition with partial pivoting and row interchanges is
    used to factor A as
       A = P * L * U,
    where P is a permutation matrix, L is unit lower triangular, and U is

    upper triangular.  The factored form of A is then used to solve the
    system of equations A * X = B.
    Arguments
    =========
    N       (input) INTEGER
            The number of linear equations, i.e., the order of the
            matrix A.  N >= 0.
    NRHS    (input) INTEGER
            The number of right hand sides, i.e., the number of columns
            of the matrix B.  NRHS >= 0.
    A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
            On 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.
    LDA     (input) INTEGER
            The leading dimension of the array A.  LDA >= max(1,N).
    IPIV    (output) INTEGER array, dimension (N)
            The pivot indices that define the permutation matrix P;
            row i of the matrix was interchanged with row IPIV(i).
    B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
            On entry, the N-by-NRHS matrix of right hand side matrix B.
            On exit, if INFO = 0, the N-by-NRHS solution matrix X.
    LDB     (input) INTEGER
            The leading dimension of the array B.  LDB >= max(1,N).
    INFO    (output) INTEGER
            = 0:  successful exit
            < 0:  if INFO = -i, the i-th argument had an illegal value
            > 0:  if INFO = i, U(i,i) is exactly zero.  The factorization

                  has been completed, but the factor U is exactly
                  singular, so the solution could not be computed.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, size(A, 1)]  
RealB[size(A, 1), :]  

Outputs

TypeNameDescription
RealX[size(A, 1), size(B, 2)] 
Integerinfo 

Modelica definition

function dgesv 
  "Solve real system of linear equations A*X=B with a B matrix"
  extends Modelica.Icons.Function;
  input Real A[:, size(A, 1)];
  input Real B[size(A, 1), :];
  output Real X[size(A, 1), size(B, 2)]=B;
  output Integer info;
protected 
  Real Awork[size(A, 1), size(A, 1)]=A;
  Integer ipiv[size(A, 1)];

external "FORTRAN 77" dgesv(size(A, 1), size(B, 2), Awork, size(A, 1), ipiv,
     X, size(A, 1), info);
end dgesv;

Modelica.Math.Matrices.LAPACK.dgesv_vec Modelica.Math.Matrices.LAPACK.dgesv_vec

Solve real system of linear equations A*x=b with a b vector

Information

Same as function LAPACK.dgesv, but right hand side is a vector and not a matrix.
For details of the arguments, see documentation of dgesv.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, size(A, 1)]  
Realb[size(A, 1)]  

Outputs

TypeNameDescription
Realx[size(A, 1)] 
Integerinfo 

Modelica definition

function dgesv_vec 
  "Solve real system of linear equations A*x=b with a b vector"
  extends Modelica.Icons.Function;
  input Real A[:, size(A, 1)];
  input Real b[size(A, 1)];
  output Real x[size(A, 1)]=b;
  output Integer info;
protected 
  Real Awork[size(A, 1), size(A, 1)]=A;
  Integer ipiv[size(A, 1)];

external "FORTRAN 77" dgesv(size(A, 1), 1, Awork, size(A, 1), ipiv, x, size(
    A, 1), info);
end dgesv_vec;

Modelica.Math.Matrices.LAPACK.dgglse_vec Modelica.Math.Matrices.LAPACK.dgglse_vec

Solve a linear equality constrained least squares problem

Information

Lapack documentation

  Purpose
  =======

  DGGLSE solves the linear equality constrained least squares (LSE)
  problem:

          minimize || A*x - c ||_2   subject to B*x = d

  using a generalized RQ factorization of matrices A and B, where A is
  M-by-N, B is P-by-N, assume P <= N <= M+P, and ||.||_2 denotes vector
  2-norm. It is assumed that

                       rank(B) = P                                  (1)

  and the null spaces of A and B intersect only trivially, i.e.,

   intersection of Null(A) and Null(B) = {0} <=> rank( ( A ) ) = N  (2)
                                                     ( ( B ) )

  where N(A) denotes the null space of matrix A. Conditions (1) and (2)
  ensure that the problem LSE has a unique solution.

  Arguments
  =========

  M       (input) INTEGER
          The number of rows of the matrix A.  M >= 0.

  N       (input) INTEGER
          The number of columns of the matrices A and B. N >= 0.
          Assume that P <= N <= M+P.

  P       (input) INTEGER
          The number of rows of the matrix B.  P >= 0.

  A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the P-by-M matrix A.
          On exit, A is destroyed.

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

  B       (input/output) DOUBLE PRECISION array, dimension (LDB,N)
          On entry, the P-by-N matrix B.
          On exit, B is destroyed.

  LDB     (input) INTEGER
          The leading dimension of the array B. LDB >= max(1,P).

  C       (input/output) DOUBLE PRECISION array, dimension (M)
          On entry, C contains the right hand side vector for the
          least squares part of the LSE problem.
          On exit, the residual sum of squares for the solution
          is given by the sum of squares of elements N-P+1 to M of
          vector C.

  D       (input/output) DOUBLE PRECISION array, dimension (P)
          On entry, D contains the right hand side vector for the
          constrained equation.
          On exit, D is destroyed.

  X       (output) DOUBLE PRECISION array, dimension (N)
          On exit, X is the solution of the LSE problem.

  WORK    (workspace) DOUBLE PRECISION array, dimension (LWORK)
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.

  LWORK   (input) INTEGER
          The dimension of the array WORK. LWORK >= N+P+max(N,M,P).
          For optimum performance LWORK >=
          N+P+max(M,P,N)*max(NB1,NB2), where NB1 is the optimal
          blocksize for the QR factorization of M-by-N matrix A.
          NB2 is the optimal blocksize for the RQ factorization of
          P-by-N matrix B.

  INFO    (output) INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, :] Minimize |A*x - c|^2
Realc[size(A, 1)]  
RealB[:, size(A, 2)] subject to B*x=d
Reald[size(B, 1)]  

Outputs

TypeNameDescription
Realx[size(A, 2)]solution vector
Integerinfo 

Modelica definition

function dgglse_vec 
  "Solve a linear equality constrained least squares problem"
  extends Modelica.Icons.Function;
  input Real A[:,:] "Minimize |A*x - c|^2";
  input Real c[size(A,1)];
  input Real B[:,size(A,2)] "subject to B*x=d";
  input Real d[size(B,1)];
  output Real x[size(A,2)] "solution vector";
  output Integer info;
protected 
  Integer nrow_A=size(A,1);
  Integer nrow_B=size(B,1);
  Integer ncol_A=size(A,2) "(min=nrow_B,max=nrow_A+nrow_B) required";
  Real Awork[nrow_A,ncol_A]=A;
  Real Bwork[nrow_B,ncol_A]=B;
  Real cwork[nrow_A] = c;
  Real dwork[nrow_B] = d;
  Integer lwork=ncol_A + nrow_B + max(nrow_A, max(ncol_A, nrow_B))*5;
  Real work[lwork];
  external "FORTRAN 77" dgglse(nrow_A, ncol_A, nrow_B, Awork, nrow_A,
                               Bwork, nrow_B, cwork, dwork, x,
                               work, lwork, info);

end dgglse_vec;

Modelica.Math.Matrices.LAPACK.dgtsv Modelica.Math.Matrices.LAPACK.dgtsv

Solve real system of linear equations A*X=B with B matrix and tridiagonal A

Information

Lapack documentation:
    Purpose
    =======
    DGTSV  solves the equation
       A*X = B,
    where A is an N-by-N tridiagonal matrix, by Gaussian elimination with

    partial pivoting.
    Note that the equation  A'*X = B  may be solved by interchanging the

    order of the arguments DU and DL.
    Arguments
    =========
    N       (input) INTEGER
            The order of the matrix A.  N >= 0.
    NRHS    (input) INTEGER
            The number of right hand sides, i.e., the number of columns
            of the matrix B.  NRHS >= 0.
    DL      (input/output) DOUBLE PRECISION array, dimension (N-1)
            On entry, DL must contain the (n-1) subdiagonal elements of
            A.
            On exit, DL is overwritten by the (n-2) elements of the
            second superdiagonal of the upper triangular matrix U from
            the LU factorization of A, in DL(1), ..., DL(n-2).
    D       (input/output) DOUBLE PRECISION array, dimension (N)
            On entry, D must contain the diagonal elements of A.
            On exit, D is overwritten by the n diagonal elements of U.
    DU      (input/output) DOUBLE PRECISION array, dimension (N-1)
            On entry, DU must contain the (n-1) superdiagonal elements
            of A.
            On exit, DU is overwritten by the (n-1) elements of the first

            superdiagonal of U.
    B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
            On entry, the N-by-NRHS right hand side matrix B.
            On exit, if INFO = 0, the N-by-NRHS solution matrix X.
    LDB     (input) INTEGER
            The leading dimension of the array B.  LDB >= max(1,N).
    INFO    (output) INTEGER
            = 0:  successful exit
            < 0:  if INFO = -i, the i-th argument had an illegal value
            > 0:  if INFO = i, U(i,i) is exactly zero, and the solution
                  has not been computed.  The factorization has not been

                  completed unless i = N.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
Realsuperdiag[:]  
Realdiag[size(superdiag, 1) + 1]  
Realsubdiag[size(superdiag, 1)]  
RealB[size(diag, 1), :]  

Outputs

TypeNameDescription
RealX[size(B, 1), size(B, 2)] 
Integerinfo 

Modelica definition

function dgtsv 
  "Solve real system of linear equations A*X=B with B matrix and tridiagonal A"

  extends Modelica.Icons.Function;
  input Real superdiag[:];
  input Real diag[size(superdiag, 1) + 1];
  input Real subdiag[size(superdiag, 1)];
  input Real B[size(diag, 1), :];
  output Real X[size(B, 1), size(B, 2)]=B;
  output Integer info;
protected 
  Real superdiagwork[size(superdiag, 1)]=superdiag;
  Real diagwork[size(diag, 1)]=diag;
  Real subdiagwork[size(subdiag, 1)]=subdiag;

external "FORTRAN 77" dgtsv(size(diag, 1), size(B, 2), subdiagwork,
    diagwork, superdiagwork, X, size(B, 1), info);
end dgtsv;

Modelica.Math.Matrices.LAPACK.dgtsv_vec Modelica.Math.Matrices.LAPACK.dgtsv_vec

Solve real system of linear equations A*x=b with b vector and tridiagonal A

Information

Same as function LAPACK.dgtsv, but right hand side is a vector and not a matrix.
For details of the arguments, see documentation of dgtsv.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
Realsuperdiag[:]  
Realdiag[size(superdiag, 1) + 1]  
Realsubdiag[size(superdiag, 1)]  
Realb[size(diag, 1)]  

Outputs

TypeNameDescription
Realx[size(b, 1)] 
Integerinfo 

Modelica definition

function dgtsv_vec 
  "Solve real system of linear equations A*x=b with b vector and tridiagonal A"

  extends Modelica.Icons.Function;
  input Real superdiag[:];
  input Real diag[size(superdiag, 1) + 1];
  input Real subdiag[size(superdiag, 1)];
  input Real b[size(diag, 1)];
  output Real x[size(b, 1)]=b;
  output Integer info;
protected 
  Real superdiagwork[size(superdiag, 1)]=superdiag;
  Real diagwork[size(diag, 1)]=diag;
  Real subdiagwork[size(subdiag, 1)]=subdiag;

external "FORTRAN 77" dgtsv(size(diag, 1), 1, subdiagwork, diagwork,
    superdiagwork, x, size(b, 1), info);
end dgtsv_vec;

Modelica.Math.Matrices.LAPACK.dgbsv Modelica.Math.Matrices.LAPACK.dgbsv

Solve real system of linear equations A*X=B with a B matrix

Information

Lapack documentation:
Purpose
=======
DGBSV computes the solution to a real system of linear equations
A * X = B, where A is a band matrix of order N with KL subdiagonals
and KU superdiagonals, and X and B are N-by-NRHS matrices.
The LU decomposition with partial pivoting and row interchanges is
used to factor A as A = L * U, where L is a product of permutation
and unit lower triangular matrices with KL subdiagonals, and U is
upper triangular with KL+KU superdiagonals.  The factored form of A
is then used to solve the system of equations A * X = B.
Arguments
=========
N       (input) INTEGER
        The number of linear equations, i.e., the order of the
        matrix A.  N >= 0.
KL      (input) INTEGER
        The number of subdiagonals within the band of A.  KL >= 0.
KU      (input) INTEGER
        The number of superdiagonals within the band of A.  KU >= 0.
NRHS    (input) INTEGER
        The number of right hand sides, i.e., the number of columns
        of the matrix B.  NRHS >= 0.
AB      (input/output) DOUBLE PRECISION array, dimension (LDAB,N)
        On 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.
        See below for further details.
LDAB    (input) INTEGER
        The leading dimension of the array AB.  LDAB >= 2*KL+KU+1.
IPIV    (output) INTEGER array, dimension (N)
        The pivot indices that define the permutation matrix P;
        row i of the matrix was interchanged with row IPIV(i).
B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
        On entry, the N-by-NRHS right hand side matrix B.
        On exit, if INFO = 0, the N-by-NRHS solution matrix X.
LDB     (input) INTEGER
        The leading dimension of the array B.  LDB >= max(1,N).
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument had an illegal value
        > 0:  if INFO = i, U(i,i) is exactly zero.  The factorization
              has been completed, but the factor U is exactly
              singular, and the solution has not been computed.
Further Details
===============
The band storage scheme is illustrated by the following example, when
M = N = 6, KL = 2, KU = 1:
On entry:                       On exit:
    *    *    *    +    +    +       *    *    *   u14  u25  u36
    *    *    +    +    +    +       *    *   u13  u24  u35  u46
    *   a12  a23  a34  a45  a56      *   u12  u23  u34  u45  u56
   a11  a22  a33  a44  a55  a66     u11  u22  u33  u44  u55  u66
   a21  a32  a43  a54  a65   *      m21  m32  m43  m54  m65   *
   a31  a42  a53  a64   *    *      m31  m42  m53  m64   *    *
Array elements marked * are not used by the routine; elements marked
+ need not be set on entry, but are required by the routine to store
elements of U because of fill-in resulting from the row interchanges.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
Integern Number of equations
IntegerkLower Number of lower bands
IntegerkUpper Number of upper bands
RealA[2*kLower + kUpper + 1, n]  
RealB[n, :]  

Outputs

TypeNameDescription
RealX[n, size(B, 2)] 
Integerinfo 

Modelica definition

function dgbsv 
  "Solve real system of linear equations A*X=B with a B matrix"
  extends Modelica.Icons.Function;
  input Integer n "Number of equations";
  input Integer kLower "Number of lower bands";
  input Integer kUpper "Number of upper bands";
  input Real A[2*kLower + kUpper + 1, n];
  input Real B[n, :];
  output Real X[n, size(B, 2)]=B;
  output Integer info;
protected 
  Real Awork[size(A, 1), size(A, 2)]=A;
  Integer ipiv[n];

external "FORTRAN 77" dgbsv(n, kLower, kUpper, size(B, 2), Awork, size(
    Awork, 1), ipiv, X, n, info);
end dgbsv;

Modelica.Math.Matrices.LAPACK.dgbsv_vec Modelica.Math.Matrices.LAPACK.dgbsv_vec

Solve real system of linear equations A*x=b with a b vector

Information

Lapack documentation:

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
Integern Number of equations
IntegerkLower Number of lower bands
IntegerkUpper Number of upper bands
RealA[2*kLower + kUpper + 1, n]  
Realb[n]  

Outputs

TypeNameDescription
Realx[n] 
Integerinfo 

Modelica definition

function dgbsv_vec 
  "Solve real system of linear equations A*x=b with a b vector"
  extends Modelica.Icons.Function;
  input Integer n "Number of equations";
  input Integer kLower "Number of lower bands";
  input Integer kUpper "Number of upper bands";
  input Real A[2*kLower + kUpper + 1, n];
  input Real b[n];
  output Real x[n]=b;
  output Integer info;
protected 
  Real Awork[size(A, 1), size(A, 2)]=A;
  Integer ipiv[n];

external "FORTRAN 77" dgbsv(n, kLower, kUpper, 1, Awork, size(Awork, 1),
    ipiv, x, n, info);
end dgbsv_vec;

Modelica.Math.Matrices.LAPACK.dgesvd Modelica.Math.Matrices.LAPACK.dgesvd

Determine singular value decomposition

Information

Lapack documentation:
    Purpose
    =======
    DGESVD computes the singular value decomposition (SVD) of a real
    M-by-N matrix A, optionally computing the left and/or right singular

    vectors. The SVD is written
         A = U * SIGMA * transpose(V)
    where SIGMA is an M-by-N matrix which is zero except for its
    min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and
    V is an N-by-N orthogonal matrix.  The diagonal elements of SIGMA
    are the singular values of A; they are real and non-negative, and
    are returned in descending order.  The first min(m,n) columns of
    U and V are the left and right singular vectors of A.
    Note that the routine returns V**T, not V.
    Arguments
    =========
    JOBU    (input) CHARACTER*1
            Specifies options for computing all or part of the matrix U:

            = '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.
    JOBVT   (input) CHARACTER*1
            Specifies options for computing all or part of the matrix
            V**T:
            = '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.
            JOBVT and JOBU cannot both be 'O'.
    M       (input) INTEGER
            The number of rows of the input matrix A.  M >= 0.
    N       (input) INTEGER
            The number of columns of the input matrix A.  N >= 0.
    A       (input/output) 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.
    LDA     (input) INTEGER
            The leading dimension of the array A.  LDA >= max(1,M).
    S       (output) DOUBLE PRECISION array, dimension (min(M,N))
            The singular values of A, sorted so that S(i) >= S(i+1).
    U       (output) 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.
    LDU     (input) INTEGER
            The leading dimension of the array U.  LDU >= 1; if
            JOBU = 'S' or 'A', LDU >= M.
    VT      (output) 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.
    LDVT    (input) INTEGER
            The leading dimension of the array VT.  LDVT >= 1; if
            JOBVT = 'A', LDVT >= N; if JOBVT = 'S', LDVT >= min(M,N).
    WORK    (workspace/output) DOUBLE PRECISION array, dimension (LWORK)

            On exit, if INFO = 0, WORK(1) returns the optimal LWORK;
            if INFO > 0, WORK(2:MIN(M,N)) contains the unconverged
            superdiagonal elements of an upper bidiagonal matrix B
            whose diagonal is in S (not necessarily sorted). B
            satisfies A = U * B * VT, so it has the same singular values

            as A, and singular vectors related by U and VT.
    LWORK   (input) INTEGER
            The dimension of the array WORK. LWORK >= 1.
            LWORK >= MAX(3*MIN(M,N)+MAX(M,N),5*MIN(M,N)-4).
            For good performance, LWORK should generally be larger.
    INFO    (output) INTEGER
            = 0:  successful exit.
            < 0:  if INFO = -i, the i-th argument had an illegal value.
            > 0:  if DBDSQR did not converge, INFO specifies how many
                  superdiagonals of an intermediate bidiagonal form B
                  did not converge to zero. See the description of WORK
                  above for details.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, :]  

Outputs

TypeNameDescription
Realsigma[min(size(A, 1), size(A, 2))] 
RealU[size(A, 1), size(A, 1)] 
RealVT[size(A, 2), size(A, 2)] 
Integerinfo 

Modelica definition

function dgesvd "Determine singular value decomposition"
  extends Modelica.Icons.Function;
  input Real A[:, :];
  output Real sigma[min(size(A, 1), size(A, 2))];
  output Real U[size(A, 1), size(A, 1)]=zeros(size(A, 1), size(A, 1));
  output Real VT[size(A, 2), size(A, 2)]=zeros(size(A, 2), size(A, 2));
  output Integer info;
protected 
  Real Awork[size(A, 1), size(A, 2)]=A;
  Integer lwork=5*size(A, 1) + 5*size(A, 2);
  Real work[lwork];

external "Fortran 77" dgesvd("A", "A", size(A, 1), size(A, 2), Awork, size(
    A, 1), sigma, U, size(A, 1), VT, size(A, 2), work, lwork, info);
end dgesvd;

Modelica.Math.Matrices.LAPACK.dgesvd_sigma Modelica.Math.Matrices.LAPACK.dgesvd_sigma

Determine singular values

Information

Lapack documentation:
    Purpose
    =======
    DGESVD computes the singular value decomposition (SVD) of a real
    M-by-N matrix A, optionally computing the left and/or right singular

    vectors. The SVD is written
         A = U * SIGMA * transpose(V)
    where SIGMA is an M-by-N matrix which is zero except for its
    min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and
    V is an N-by-N orthogonal matrix.  The diagonal elements of SIGMA
    are the singular values of A; they are real and non-negative, and
    are returned in descending order.  The first min(m,n) columns of
    U and V are the left and right singular vectors of A.
    Note that the routine returns V**T, not V.
    Arguments
    =========
    JOBU    (input) CHARACTER*1
            Specifies options for computing all or part of the matrix U:

            = '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.
    JOBVT   (input) CHARACTER*1
            Specifies options for computing all or part of the matrix
            V**T:
            = '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.
            JOBVT and JOBU cannot both be 'O'.
    M       (input) INTEGER
            The number of rows of the input matrix A.  M >= 0.
    N       (input) INTEGER
            The number of columns of the input matrix A.  N >= 0.
    A       (input/output) 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.
    LDA     (input) INTEGER
            The leading dimension of the array A.  LDA >= max(1,M).
    S       (output) DOUBLE PRECISION array, dimension (min(M,N))
            The singular values of A, sorted so that S(i) >= S(i+1).
    U       (output) 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.
    LDU     (input) INTEGER
            The leading dimension of the array U.  LDU >= 1; if
            JOBU = 'S' or 'A', LDU >= M.
    VT      (output) 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.
    LDVT    (input) INTEGER
            The leading dimension of the array VT.  LDVT >= 1; if
            JOBVT = 'A', LDVT >= N; if JOBVT = 'S', LDVT >= min(M,N).
    WORK    (workspace/output) DOUBLE PRECISION array, dimension (LWORK)

            On exit, if INFO = 0, WORK(1) returns the optimal LWORK;
            if INFO > 0, WORK(2:MIN(M,N)) contains the unconverged
            superdiagonal elements of an upper bidiagonal matrix B
            whose diagonal is in S (not necessarily sorted). B
            satisfies A = U * B * VT, so it has the same singular values

            as A, and singular vectors related by U and VT.
    LWORK   (input) INTEGER
            The dimension of the array WORK. LWORK >= 1.
            LWORK >= MAX(3*MIN(M,N)+MAX(M,N),5*MIN(M,N)-4).
            For good performance, LWORK should generally be larger.
    INFO    (output) INTEGER
            = 0:  successful exit.
            < 0:  if INFO = -i, the i-th argument had an illegal value.
            > 0:  if DBDSQR did not converge, INFO specifies how many
                  superdiagonals of an intermediate bidiagonal form B
                  did not converge to zero. See the description of WORK
                  above for details.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, :]  

Outputs

TypeNameDescription
Realsigma[min(size(A, 1), size(A, 2))] 
Integerinfo 

Modelica definition

function dgesvd_sigma "Determine singular values"
  extends Modelica.Icons.Function;
  input Real A[:, :];
  output Real sigma[min(size(A, 1), size(A, 2))];
  output Integer info;
protected 
  Real Awork[size(A, 1), size(A, 2)]=A;
  Real U[size(A, 1), size(A, 1)];
  Real VT[size(A, 2), size(A, 2)];
  Integer lwork=5*size(A, 1) + 5*size(A, 2);
  Real work[lwork];

external "Fortran 77" dgesvd("N", "N", size(A, 1), size(A, 2), Awork, size(
    A, 1), sigma, U, size(A, 1), VT, size(A, 2), work, lwork, info);
end dgesvd_sigma;

Modelica.Math.Matrices.LAPACK.dgetrf Modelica.Math.Matrices.LAPACK.dgetrf

Compute LU factorization of square or rectangular matrix A (A = P*L*U)

Information

Lapack documentation:
  SUBROUTINE DGETRF( M, N, A, LDA, IPIV, INFO )
-- LAPACK routine (version 1.1) --
   Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
   Courant Institute, Argonne National Lab, and Rice University
   March 31, 1993
   .. Scalar Arguments ..
   INTEGER            INFO, LDA, M, N
   ..
   .. Array Arguments ..
   INTEGER            IPIV( * )
   DOUBLE PRECISION   A( LDA, * )
   ..
Purpose
=======
DGETRF computes an LU factorization of a general M-by-N matrix A
using partial pivoting with row interchanges.
The factorization has the form
   A = P * L * U
where P is a permutation matrix, L is lower triangular with unit
diagonal elements (lower trapezoidal if m > n), and U is upper
triangular (upper trapezoidal if m < n).
This is the right-looking Level 3 BLAS version of the algorithm.
Arguments
=========
M       (input) INTEGER
        The number of rows of the matrix A.  M >= 0.
N       (input) INTEGER
        The number of columns of the matrix A.  N >= 0.
A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
        On 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.
LDA     (input) INTEGER
        The leading dimension of the array A.  LDA >= max(1,M).
IPIV    (output) INTEGER array, dimension (min(M,N))
        The pivot indices; for 1 <= i <= min(M,N), row i of the
        matrix was interchanged with row IPIV(i).
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument had an illegal value
        > 0:  if INFO = i, U(i,i) is exactly zero. The factorization
              has been completed, but the factor U is exactly
              singular, and division by zero will occur if it is used
              to solve a system of equations.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, :] Square or rectangular matrix

Outputs

TypeNameDescription
RealLU[size(A, 1), size(A, 2)] 
Integerpivots[min(size(A, 1), size(A, 2))]Pivot vector
IntegerinfoInformation

Modelica definition

function dgetrf 
  "Compute LU factorization of square or rectangular matrix A (A = P*L*U)"

  extends Modelica.Icons.Function;
  input Real A[:, :] "Square or rectangular matrix";
  output Real LU[size(A, 1), size(A, 2)]=A;
  output Integer pivots[min(size(A, 1), size(A, 2))] "Pivot vector";
  output Integer info "Information";

external "FORTRAN 77" dgetrf(size(A, 1), size(A, 2), LU, size(A, 1), pivots,
     info);
end dgetrf;

Modelica.Math.Matrices.LAPACK.dgetrs Modelica.Math.Matrices.LAPACK.dgetrs

Solves a system of linear equations with the LU decomposition from dgetrf(..)

Information

Lapack documentation:
  SUBROUTINE DGETRS( TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO )
-- LAPACK routine (version 1.1) --
   Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
   Courant Institute, Argonne National Lab, and Rice University
   March 31, 1993
   .. Scalar Arguments ..
   CHARACTER          TRANS
   INTEGER            INFO, LDA, LDB, N, NRHS
   ..
   .. Array Arguments ..
   INTEGER            IPIV( * )
   DOUBLE PRECISION   A( LDA, * ), B( LDB, * )
   ..
Purpose
=======
DGETRS solves a system of linear equations
   A * X = B  or  A' * X = B
with a general N-by-N matrix A using the LU factorization computed
by DGETRF.
Arguments
=========
TRANS   (input) CHARACTER*1
        Specifies the form of the system of equations:
        = 'N':  A * X = B  (No transpose)
        = 'T':  A'* X = B  (Transpose)
        = 'C':  A'* X = B  (Conjugate transpose = Transpose)
N       (input) INTEGER
        The order of the matrix A.  N >= 0.
NRHS    (input) INTEGER
        The number of right hand sides, i.e., the number of columns
        of the matrix B.  NRHS >= 0.
A       (input) DOUBLE PRECISION array, dimension (LDA,N)
        The factors L and U from the factorization A = P*L*U
        as computed by DGETRF.
LDA     (input) INTEGER
        The leading dimension of the array A.  LDA >= max(1,N).
IPIV    (input) INTEGER array, dimension (N)
        The pivot indices from DGETRF; for 1<=i<=N, row i of the
        matrix was interchanged with row IPIV(i).
B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
        On entry, the right hand side matrix B.
        On exit, the solution matrix X.
LDB     (input) INTEGER
        The leading dimension of the array B.  LDB >= max(1,N).
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument had an illegal value

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealLU[:, size(LU, 1)] LU factorization of dgetrf of a square matrix
Integerpivots[size(LU, 1)] Pivot vector of dgetrf
RealB[size(LU, 1), :] Right hand side matrix B

Outputs

TypeNameDescription
RealX[size(B, 1), size(B, 2)]Solution matrix X

Modelica definition

function dgetrs 
  "Solves a system of linear equations with the LU decomposition from dgetrf(..)"

  extends Modelica.Icons.Function;
  input Real LU[:, size(LU, 1)] "LU factorization of dgetrf of a square matrix";
  input Integer pivots[size(LU, 1)] "Pivot vector of dgetrf";
  input Real B[size(LU, 1),:] "Right hand side matrix B";
  output Real X[size(B, 1), size(B,2)]=B "Solution matrix X";

protected 
  Real work[size(LU, 1), size(LU, 1)]=LU;
  Integer info;
external "FORTRAN 77" dgetrs("N", size(LU, 1), size(B,2), work, size(LU, 1), pivots,
     X, size(B, 1), info);
end dgetrs;

Modelica.Math.Matrices.LAPACK.dgetrs_vec Modelica.Math.Matrices.LAPACK.dgetrs_vec

Solves a system of linear equations with the LU decomposition from dgetrf(..)

Information

Lapack documentation:
  SUBROUTINE DGETRS( TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO )
-- LAPACK routine (version 1.1) --
   Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
   Courant Institute, Argonne National Lab, and Rice University
   March 31, 1993
   .. Scalar Arguments ..
   CHARACTER          TRANS
   INTEGER            INFO, LDA, LDB, N, NRHS
   ..
   .. Array Arguments ..
   INTEGER            IPIV( * )
   DOUBLE PRECISION   A( LDA, * ), B( LDB, * )
   ..
Purpose
=======
DGETRS solves a system of linear equations
   A * X = B  or  A' * X = B
with a general N-by-N matrix A using the LU factorization computed
by DGETRF.
Arguments
=========
TRANS   (input) CHARACTER*1
        Specifies the form of the system of equations:
        = 'N':  A * X = B  (No transpose)
        = 'T':  A'* X = B  (Transpose)
        = 'C':  A'* X = B  (Conjugate transpose = Transpose)
N       (input) INTEGER
        The order of the matrix A.  N >= 0.
NRHS    (input) INTEGER
        The number of right hand sides, i.e., the number of columns
        of the matrix B.  NRHS >= 0.
A       (input) DOUBLE PRECISION array, dimension (LDA,N)
        The factors L and U from the factorization A = P*L*U
        as computed by DGETRF.
LDA     (input) INTEGER
        The leading dimension of the array A.  LDA >= max(1,N).
IPIV    (input) INTEGER array, dimension (N)
        The pivot indices from DGETRF; for 1<=i<=N, row i of the
        matrix was interchanged with row IPIV(i).
B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
        On entry, the right hand side matrix B.
        On exit, the solution matrix X.
LDB     (input) INTEGER
        The leading dimension of the array B.  LDB >= max(1,N).
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument had an illegal value

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealLU[:, size(LU, 1)] LU factorization of dgetrf of a square matrix
Integerpivots[size(LU, 1)] Pivot vector of dgetrf
Realb[size(LU, 1)] Right hand side vector b

Outputs

TypeNameDescription
Realx[size(b, 1)] 

Modelica definition

function dgetrs_vec 
  "Solves a system of linear equations with the LU decomposition from dgetrf(..)"

  extends Modelica.Icons.Function;
  input Real LU[:, size(LU, 1)] "LU factorization of dgetrf of a square matrix";
  input Integer pivots[size(LU, 1)] "Pivot vector of dgetrf";
  input Real b[size(LU, 1)] "Right hand side vector b";
  output Real x[size(b, 1)]=b;

protected 
  Real work[size(LU, 1), size(LU, 1)]=LU;
  Integer info;
external "FORTRAN 77" dgetrs("N", size(LU, 1), 1, work, size(LU, 1), pivots,
     x, size(b, 1), info);
end dgetrs_vec;

Modelica.Math.Matrices.LAPACK.dgetri Modelica.Math.Matrices.LAPACK.dgetri

Computes the inverse of a matrix using the LU factorization from dgetrf(..)

Information

Lapack documentation:
   SUBROUTINE DGETRI( N, A, LDA, IPIV, WORK, LWORK, INFO )
-- LAPACK routine (version 1.1) --
   Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
   Courant Institute, Argonne National Lab, and Rice University
   March 31, 1993
   .. Scalar Arguments ..
   INTEGER            INFO, LDA, LWORK, N
   ..
   .. Array Arguments ..
   INTEGER            IPIV( * )
   DOUBLE PRECISION   A( LDA, * ), WORK( LWORK )
   ..
Purpose
=======
DGETRI computes the inverse of a matrix using the LU factorization
computed by DGETRF.
This method inverts U and then computes inv(A) by solving the system
inv(A)*L = inv(U) for inv(A).
Arguments
=========
N       (input) INTEGER
        The order of the matrix A.  N >= 0.
A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
        On 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.
LDA     (input) INTEGER
        The leading dimension of the array A.  LDA >= max(1,N).
IPIV    (input) INTEGER array, dimension (N)
        The pivot indices from DGETRF; for 1<=i<=N, row i of the
        matrix was interchanged with row IPIV(i).
WORK    (workspace) DOUBLE PRECISION array, dimension (LWORK)
        On exit, if INFO=0, then WORK(1) returns the optimal LWORK.
LWORK   (input) INTEGER
        The dimension of the array WORK.  LWORK >= max(1,N).
        For optimal performance LWORK >= N*NB, where NB is
        the optimal blocksize returned by ILAENV.
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument had an illegal value
        > 0:  if INFO = i, U(i,i) is exactly zero; the matrix is
              singular and its inverse could not be computed.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealLU[:, size(LU, 1)] LU factorization of dgetrf of a square matrix
Integerpivots[size(LU, 1)] Pivot vector of dgetrf

Outputs

TypeNameDescription
Realinv[size(LU, 1), size(LU, 2)]Inverse of matrix P*L*U

Modelica definition

function dgetri 
  "Computes the inverse of a matrix using the LU factorization from dgetrf(..)"

  extends Modelica.Icons.Function;
  input Real LU[:, size(LU, 1)] "LU factorization of dgetrf of a square matrix";
  input Integer pivots[size(LU, 1)] "Pivot vector of dgetrf";
  output Real inv[size(LU, 1), size(LU, 2)]=LU "Inverse of matrix P*L*U";

protected 
  Integer lwork=min(10, size(LU, 1))*size(LU, 1) "Length of work array";
  Real work[lwork];
  Integer info;
external "FORTRAN 77" dgetri(size(LU, 1), inv, size(LU, 1), pivots, work,
    lwork, info);
end dgetri;

Modelica.Math.Matrices.LAPACK.dgeqpf Modelica.Math.Matrices.LAPACK.dgeqpf

Compute QR factorization of square or rectangular matrix A with column pivoting (A(:,p) = Q*R)

Information

Lapack documentation:
   SUBROUTINE DGEQPF( M, N, A, LDA, JPVT, TAU, WORK, INFO )
-- LAPACK test routine (version 1.1) --
   Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
   Courant Institute, Argonne National Lab, and Rice University
   March 31, 1993
   .. Scalar Arguments ..
   INTEGER            INFO, LDA, M, N
   ..
   .. Array Arguments ..
   INTEGER            JPVT( * )
   DOUBLE PRECISION   A( LDA, * ), TAU( * ), WORK( * )
   ..
Purpose
=======
DGEQPF computes a QR factorization with column pivoting of a
real M-by-N matrix A: A*P = Q*R.
Arguments
=========
M       (input) INTEGER
        The number of rows of the matrix A. M >= 0.
N       (input) INTEGER
        The number of columns of the matrix A. N >= 0
A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
        On entry, the M-by-N matrix A.
        On exit, the upper triangle of the array contains the
        min(M,N)-by-N upper triangular matrix R; the elements
        below the diagonal, together with the array TAU,
        represent the orthogonal matrix Q as a product of
        min(m,n) elementary reflectors.
LDA     (input) INTEGER
        The leading dimension of the array A. LDA >= max(1,M).
JPVT    (input/output) INTEGER array, dimension (N)
        On entry, if JPVT(i) .ne. 0, the i-th column of A is permuted
        to the front of A*P (a leading column); if JPVT(i) = 0,
        the i-th column of A is a free column.
        On exit, if JPVT(i) = k, then the i-th column of A*P
        was the k-th column of A.
TAU     (output) DOUBLE PRECISION array, dimension (min(M,N))
        The scalar factors of the elementary reflectors.
WORK    (workspace) DOUBLE PRECISION array, dimension (3*N)
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument had an illegal value
Further Details
===============
The matrix Q is represented as a product of elementary reflectors
   Q = H(1) H(2) . . . H(n)
Each H(i) has the form
   H = I - tau * v * v'
where tau is a real scalar, and v is a real vector with
v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i).
The matrix P is represented in jpvt as follows: If
   jpvt(j) = i
then the jth column of P is the ith canonical unit vector.

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealA[:, :] Square or rectangular matrix

Outputs

TypeNameDescription
RealQR[size(A, 1), size(A, 2)]QR factorization in packed format
Realtau[min(size(A, 1), size(A, 2))]The scalar factors of the elementary reflectors of Q
Integerp[size(A, 2)]Pivot vector

Modelica definition

function dgeqpf 
  "Compute QR factorization of square or rectangular matrix A with column pivoting (A(:,p) = Q*R)"

  extends Modelica.Icons.Function;
  input Real A[:, :] "Square or rectangular matrix";
  output Real QR[size(A, 1), size(A, 2)]=A "QR factorization in packed format";
  output Real tau[min(size(A, 1), size(A, 2))] 
    "The scalar factors of the elementary reflectors of Q";
  output Integer p[size(A, 2)]=zeros(size(A, 2)) "Pivot vector";

protected 
  Integer info;
  Integer ncol=size(A, 2) "Column dimension of A";
  Real work[3*ncol] "work array";
external "FORTRAN 77" dgeqpf(size(A, 1), ncol, QR, size(A, 1), p, tau, work,
     info);
end dgeqpf;

Modelica.Math.Matrices.LAPACK.dorgqr Modelica.Math.Matrices.LAPACK.dorgqr

Generates a Real orthogonal matrix Q which is defined as the product of elementary reflectors as returned from dgeqpf

Information

Lapack documentation:
   SUBROUTINE DORGQR( M, N, K, A, LDA, TAU, WORK, LWORK, INFO )
-- LAPACK routine (version 1.1) --
   Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
   Courant Institute, Argonne National Lab, and Rice University
   March 31, 1993
   .. Scalar Arguments ..
   INTEGER            INFO, K, LDA, LWORK, M, N
   ..
   .. Array Arguments ..
   DOUBLE PRECISION   A( LDA, * ), TAU( * ), WORK( LWORK )
   ..
Purpose
=======
DORGQR generates an M-by-N real matrix Q with orthonormal columns,
which is defined as the first N columns of a product of K elementary
reflectors of order M
      Q  =  H(1) H(2) . . . H(k)
as returned by DGEQRF.
Arguments
=========
M       (input) INTEGER
        The number of rows of the matrix Q. M >= 0.
N       (input) INTEGER
        The number of columns of the matrix Q. M >= N >= 0.
K       (input) INTEGER
        The number of elementary reflectors whose product defines the
        matrix Q. N >= K >= 0.
A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
        On entry, the i-th column must contain the vector which
        defines the elementary reflector H(i), for i = 1,2,...,k, as
        returned by DGEQRF in the first k columns of its array
        argument A.
        On exit, the M-by-N matrix Q.
LDA     (input) INTEGER
        The first dimension of the array A. LDA >= max(1,M).
TAU     (input) DOUBLE PRECISION array, dimension (K)
        TAU(i) must contain the scalar factor of the elementary
        reflector H(i), as returned by DGEQRF.
WORK    (workspace) DOUBLE PRECISION array, dimension (LWORK)
        On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
LWORK   (input) INTEGER
        The dimension of the array WORK. LWORK >= max(1,N).
        For optimum performance LWORK >= N*NB, where NB is the
        optimal blocksize.
INFO    (output) INTEGER
        = 0:  successful exit
        < 0:  if INFO = -i, the i-th argument has an illegal value

Extends from Modelica.Icons.Function (Icon for a function).

Inputs

TypeNameDefaultDescription
RealQR[:, :] QR from dgeqpf
Realtau[min(size(QR, 1), size(QR, 2))] The scalar factors of the elementary reflectors of Q

Outputs

TypeNameDescription
RealQ[size(QR, 1), size(QR, 2)]Orthogonal matrix Q

Modelica definition

function dorgqr 
  "Generates a Real orthogonal matrix Q which is defined as the product of elementary reflectors as returned from dgeqpf"

  extends Modelica.Icons.Function;
  input Real QR[:, :] "QR from dgeqpf";
  input Real tau[min(size(QR, 1), size(QR, 2))] 
    "The scalar factors of the elementary reflectors of Q";
  output Real Q[size(QR, 1), size(QR, 2)]=QR "Orthogonal matrix Q";

protected 
  Integer info;
  Integer lwork=min(10, size(QR, 2))*size(QR, 2) "Length of work array";
  Real work[lwork];
external "FORTRAN 77" dorgqr(size(QR, 1), size(QR, 2), size(tau, 1), Q,
    size(Q, 1), tau, work, lwork, info);
end dorgqr;

HTML-documentation generated by Dymola Sun Jan 17 21:12:46 2010.