Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

dsgesv (3p)

Name

dsgesv - compute the solution to a real system of linear equations A*X = B

Synopsis

SUBROUTINE DSGESV(N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

INTEGER N, NRHS, LDA, LDB, LDX, ITER, INFO
INTEGER IPIV(*)
REAL SWORK(*)
DOUBLE PRECISION A(LDA,*), B(LDB,*), WORK(N,*), X(LDX,*)

SUBROUTINE DSGESV_64(N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

INTEGER*8 N, NRHS, LDA, LDB, LDX, ITER, INFO
INTEGER*8 IPIV(*)
REAL SWORK(*)
DOUBLE PRECISION A(LDA,*), B(LDB,*), WORK(N,*), X(LDX,*)




F95 INTERFACE
SUBROUTINE SGESV(N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

INTEGER :: M, N, LDA, INFO
INTEGER, DIMENSION(:) :: IPIV
DOUBLE PRECISION, DIMENSION(:,:) :: A, B, X, WORK

SUBROUTINE SGESV_64((N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

INTEGER(8) :: N, NRHS, LDA, LDB, LDX, ITER, INFO
INTEGER(8), DIMENSION(:) :: IPIV
DOUBLE PRECISION, DIMENSION(:,:) :: A, B, X, WORK




C INTERFACE
#include <sunperf.h>

void dsgesv(int n, int nrhs, double *a, int lda, int *ipiv, double  *b,
int  ldb, double *x, int ldx, double *work, float *swork, int
iter, int *info);

void dsgesv_64(long n, long nrhs, double *a, long lda, long *ipiv, dou-
ble  *b,  long  ldb, double *x, long ldx, double *work, float
*swork, long iter, long *info);

Description

Oracle Solaris Studio Performance Library                           dsgesv(3P)



NAME
       dsgesv  - compute the solution to a real system of linear equations A*X
       = B


SYNOPSIS
       SUBROUTINE DSGESV(N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

       INTEGER N, NRHS, LDA, LDB, LDX, ITER, INFO
       INTEGER IPIV(*)
       REAL SWORK(*)
       DOUBLE PRECISION A(LDA,*), B(LDB,*), WORK(N,*), X(LDX,*)

       SUBROUTINE DSGESV_64(N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

       INTEGER*8 N, NRHS, LDA, LDB, LDX, ITER, INFO
       INTEGER*8 IPIV(*)
       REAL SWORK(*)
       DOUBLE PRECISION A(LDA,*), B(LDB,*), WORK(N,*), X(LDX,*)




   F95 INTERFACE
       SUBROUTINE SGESV(N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

       INTEGER :: M, N, LDA, INFO
       INTEGER, DIMENSION(:) :: IPIV
       DOUBLE PRECISION, DIMENSION(:,:) :: A, B, X, WORK

       SUBROUTINE SGESV_64((N, NRHS, A, LDA, IPIV, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)

       INTEGER(8) :: N, NRHS, LDA, LDB, LDX, ITER, INFO
       INTEGER(8), DIMENSION(:) :: IPIV
       DOUBLE PRECISION, DIMENSION(:,:) :: A, B, X, WORK




   C INTERFACE
       #include <sunperf.h>

       void dsgesv(int n, int nrhs, double *a, int lda, int *ipiv, double  *b,
                 int  ldb, double *x, int ldx, double *work, float *swork, int
                 iter, int *info);

       void dsgesv_64(long n, long nrhs, double *a, long lda, long *ipiv, dou-
                 ble  *b,  long  ldb, double *x, long ldx, double *work, float
                 *swork, long iter, long *info);



PURPOSE
       dsgesv 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.
       dsgesv first attempts to factorize the matrix in SINGLE  PRECISION  and
       use this factorization within an iterative refinement procedure to pro-
       duce a solution with DOUBLE PRECISION normwise backward  error  quality
       (see below). If the approach fails the method switches to a DOUBLE PRE-
       CISION factorization and solve.

       The iterative refinement is not going to be a winning strategy  if  the
       ratio SINGLE PRECISION performance over DOUBLE PRECISION performance is
       too small. A reasonable strategy should take the number  of  right-hand
       sides  and the size of the matrix into account. This might be done with
       a call to ILAENV in the future. Up to  now,  we  always  try  iterative
       refinement.  The iterative refinement process is stopped if
           ITER > ITERMAX
       or for all the RHS we have:
           RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX
       where
           o ITER is the number of the current iteration in the iterative
             refinement process
           o RNRM is the infinity-norm of the residual
           o XNRM is the infinity-norm of the solution
           o ANRM is the infinity-operator-norm of the matrix A
           o  EPS  is  the  machine  epsilon returned by DLAMCH('Epsilon') The
       value ITERMAX and BWDMAX are fixed to 30 and 1.0D+00 respectively.


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 or input/output) DOUBLE PRECISION array
                 On entry, the N-by-N coefficient matrix A.  On exit, if iter-
                 ative  refinement  has  been successfully used (INFO.EQ.0 and
                 ITER.GE.0, see description below), then A  is  unchanged,  if
                 double  precision  factorization has been used (INFO.EQ.0 and
                 ITER.LT.0, see description below), then the array A  contains
                 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).  Corresponds
                 either to the single precision  factorization  (if  INFO.EQ.0
                 and  ITER.GE.0)  or  the  double  precision factorization (if
                 INFO.EQ.0 and ITER.LT.0).


       B (input)  DOUBLE PRECISION array, dimension (LDB,NRHS)
                 The N-by-NRHS matrix of right hand side matrix B.


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


       X (output) DOUBLE PRECISION array, dimension (LDX,NRHS)
                 If INFO = 0, the N-by-NRHS solution matrix X.


       LDX (input) INTEGER
                 The leading dimension of the array X.  LDX >= max(1,N).


       WORK    (workspace) DOUBLE PRECISION array, dimension (N,NRHS)
                 This array is used to hold the residual vectors.


       SWORK   (workspace) REAL array, dimension (N*(N+NRHS))
                 This array is used to use the single precision matrix and the
                 right-hand sides or solutions in single precision.


       ITER    (output) INTEGER
                 <  0:  iterative refinement has failed, double precision fac-
                 torization has been performed
                 -1 : taking into account machine parameters, N, NRHS, it is a
                 priori not worth working in SINGLE PRECISION
                 -2  :  overflow of an entry when moving from double to SINGLE
                 PRECISION
                 -3 : failure of SGETRF
                 -31: stop the iterative refinement after the 30th iterations
                 > 0: iterative refinement has been sucessfully used.  Returns
                 the number of iterations


       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) computed in DOUBLE PRECISION is
                 exactly zero.  The factorization has been completed, but  the
                 factor  U  is  exactly singular, so the solution could not be
                 computed.



                                  7 Nov 2015                        dsgesv(3P)