1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140
|
*DECK SGEDI
SUBROUTINE SGEDI (A, LDA, N, IPVT, DET, WORK, JOB)
C***BEGIN PROLOGUE SGEDI
C***PURPOSE Compute the determinant and inverse of a matrix using the
C factors computed by SGECO or SGEFA.
C***LIBRARY SLATEC (LINPACK)
C***CATEGORY D2A1, D3A1
C***TYPE SINGLE PRECISION (SGEDI-S, DGEDI-D, CGEDI-C)
C***KEYWORDS DETERMINANT, INVERSE, LINEAR ALGEBRA, LINPACK, MATRIX
C***AUTHOR Moler, C. B., (U. of New Mexico)
C***DESCRIPTION
C
C SGEDI computes the determinant and inverse of a matrix
C using the factors computed by SGECO or SGEFA.
C
C On Entry
C
C A REAL(LDA, N)
C the output from SGECO or SGEFA.
C
C LDA INTEGER
C the leading dimension of the array A .
C
C N INTEGER
C the order of the matrix A .
C
C IPVT INTEGER(N)
C the pivot vector from SGECO or SGEFA.
C
C WORK REAL(N)
C work vector. Contents destroyed.
C
C JOB INTEGER
C = 11 both determinant and inverse.
C = 01 inverse only.
C = 10 determinant only.
C
C On Return
C
C A inverse of original matrix if requested.
C Otherwise unchanged.
C
C DET REAL(2)
C determinant of original matrix if requested.
C Otherwise not referenced.
C Determinant = DET(1) * 10.0**DET(2)
C with 1.0 .LE. ABS(DET(1)) .LT. 10.0
C or DET(1) .EQ. 0.0 .
C
C Error Condition
C
C A division by zero will occur if the input factor contains
C a zero on the diagonal and the inverse is requested.
C It will not occur if the subroutines are called correctly
C and if SGECO has set RCOND .GT. 0.0 or SGEFA has set
C INFO .EQ. 0 .
C
C***REFERENCES J. J. Dongarra, J. R. Bunch, C. B. Moler, and G. W.
C Stewart, LINPACK Users' Guide, SIAM, 1979.
C***ROUTINES CALLED SAXPY, SSCAL, SSWAP
C***REVISION HISTORY (YYMMDD)
C 780814 DATE WRITTEN
C 890831 Modified array declarations. (WRB)
C 890831 REVISION DATE from Version 3.2
C 891214 Prologue converted to Version 4.0 format. (BAB)
C 900326 Removed duplicate information from DESCRIPTION section.
C (WRB)
C 920501 Reformatted the REFERENCES section. (WRB)
C***END PROLOGUE SGEDI
INTEGER LDA,N,IPVT(*),JOB
REAL A(LDA,*),DET(2),WORK(*)
C
REAL T
REAL TEN
INTEGER I,J,K,KB,KP1,L,NM1
C***FIRST EXECUTABLE STATEMENT SGEDI
C
C COMPUTE DETERMINANT
C
IF (JOB/10 .EQ. 0) GO TO 70
DET(1) = 1.0E0
DET(2) = 0.0E0
TEN = 10.0E0
DO 50 I = 1, N
IF (IPVT(I) .NE. I) DET(1) = -DET(1)
DET(1) = A(I,I)*DET(1)
IF (DET(1) .EQ. 0.0E0) GO TO 60
10 IF (ABS(DET(1)) .GE. 1.0E0) GO TO 20
DET(1) = TEN*DET(1)
DET(2) = DET(2) - 1.0E0
GO TO 10
20 CONTINUE
30 IF (ABS(DET(1)) .LT. TEN) GO TO 40
DET(1) = DET(1)/TEN
DET(2) = DET(2) + 1.0E0
GO TO 30
40 CONTINUE
50 CONTINUE
60 CONTINUE
70 CONTINUE
C
C COMPUTE INVERSE(U)
C
IF (MOD(JOB,10) .EQ. 0) GO TO 150
DO 100 K = 1, N
A(K,K) = 1.0E0/A(K,K)
T = -A(K,K)
CALL SSCAL(K-1,T,A(1,K),1)
KP1 = K + 1
IF (N .LT. KP1) GO TO 90
DO 80 J = KP1, N
T = A(K,J)
A(K,J) = 0.0E0
CALL SAXPY(K,T,A(1,K),1,A(1,J),1)
80 CONTINUE
90 CONTINUE
100 CONTINUE
C
C FORM INVERSE(U)*INVERSE(L)
C
NM1 = N - 1
IF (NM1 .LT. 1) GO TO 140
DO 130 KB = 1, NM1
K = N - KB
KP1 = K + 1
DO 110 I = KP1, N
WORK(I) = A(I,K)
A(I,K) = 0.0E0
110 CONTINUE
DO 120 J = KP1, N
T = WORK(J)
CALL SAXPY(N,T,A(1,J),1,A(1,K),1)
120 CONTINUE
L = IPVT(K)
IF (L .NE. K) CALL SSWAP(N,A(1,K),1,A(1,L),1)
130 CONTINUE
140 CONTINUE
150 CONTINUE
RETURN
END
|