SUBROUTINE DSPSV( UPLO, N, NRHS, AP, IPIV, B, LDB, INFO ) * * -- LAPACK driver routine (version 2.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * March 31, 1993 * * .. Scalar Arguments .. CHARACTER UPLO INTEGER INFO, LDB, N, NRHS * .. * .. Array Arguments .. INTEGER IPIV( * ) DOUBLE PRECISION AP( * ), B( LDB, * ) * .. * * Purpose * ======= * * DSPSV computes the solution to a real system of linear equations * A * X = B, * where A is an N-by-N symmetric matrix stored in packed format and X * and B are N-by-NRHS matrices. * * The diagonal pivoting method is used to factor A as * A = U * D * U**T, if UPLO = 'U', or * A = L * D * L**T, if UPLO = 'L', * where U (or L) is a product of permutation and unit upper (lower) * triangular matrices, D is symmetric and block diagonal with 1-by-1 * and 2-by-2 diagonal blocks. The factored form of A is then used to * solve the system of equations A * X = B. * * Arguments * ========= * * UPLO (input) CHARACTER*1 * = 'U': Upper triangle of A is stored; * = 'L': Lower triangle of A is stored. * * 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. * * AP (input/output) DOUBLE PRECISION array, dimension (N*(N+1)/2) * On entry, the upper or lower triangle of the symmetric matrix * A, packed columnwise in a linear array. The j-th column of A * is stored in the array AP as follows: * if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j; * if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n. * See below for further details. * * On exit, the block diagonal matrix D and the multipliers used * to obtain the factor U or L from the factorization * A = U*D*U**T or A = L*D*L**T as computed by DSPTRF, stored as * a packed triangular matrix in the same storage format as A. * * IPIV (output) INTEGER array, dimension (N) * Details of the interchanges and the block structure of D, as * determined by DSPTRF. If IPIV(k) > 0, then rows and columns * k and IPIV(k) were interchanged, and D(k,k) is a 1-by-1 * diagonal block. If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, * then rows and columns k-1 and -IPIV(k) were interchanged and * D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = 'L' and * IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and * -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 * diagonal block. * * 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, D(i,i) is exactly zero. The factorization * has been completed, but the block diagonal matrix D is * exactly singular, so the solution could not be * computed. * * Further Details * =============== * * The packed storage scheme is illustrated by the following example * when N = 4, UPLO = 'U': * * Two-dimensional storage of the symmetric matrix A: * * a11 a12 a13 a14 * a22 a23 a24 * a33 a34 (aij = aji) * a44 * * Packed storage of the upper triangle of A: * * AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ] * * ===================================================================== * * .. External Functions .. LOGICAL LSAME EXTERNAL LSAME * .. * .. External Subroutines .. EXTERNAL DSPTRF, DSPTRS, XERBLA * .. * .. Intrinsic Functions .. INTRINSIC MAX * .. * .. Executable Statements .. * * Test the input parameters. * INFO = 0 IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN INFO = -1 ELSE IF( N.LT.0 ) THEN INFO = -2 ELSE IF( NRHS.LT.0 ) THEN INFO = -3 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN INFO = -7 END IF IF( INFO.NE.0 ) THEN CALL XERBLA( 'DSPSV ', -INFO ) RETURN END IF * * Compute the factorization A = U*D*U' or A = L*D*L'. * CALL DSPTRF( UPLO, N, AP, IPIV, INFO ) IF( INFO.EQ.0 ) THEN * * Solve the system A*X = B, overwriting B with X. * CALL DSPTRS( UPLO, N, NRHS, AP, IPIV, B, LDB, INFO ) * END IF RETURN * * End of DSPSV * END LOGICAL FUNCTION LSAME( CA, CB ) * * -- LAPACK auxiliary routine (version 2.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * September 30, 1994 * * .. Scalar Arguments .. CHARACTER CA, CB * .. * * Purpose * ======= * * LSAME returns .TRUE. if CA is the same letter as CB regardless of * case. * * Arguments * ========= * * CA (input) CHARACTER*1 * CB (input) CHARACTER*1 * CA and CB specify the single characters to be compared. * * ===================================================================== * * .. Intrinsic Functions .. INTRINSIC ICHAR * .. * .. Local Scalars .. INTEGER INTA, INTB, ZCODE * .. * .. Executable Statements .. * * Test if the characters are equal * LSAME = CA.EQ.CB IF( LSAME ) $ RETURN * * Now test for equivalence if both characters are alphabetic. * ZCODE = ICHAR( 'Z' ) * * Use 'Z' rather than 'A' so that ASCII can be detected on Prime * machines, on which ICHAR returns a value with bit 8 set. * ICHAR('A') on Prime machines returns 193 which is the same as * ICHAR('A') on an EBCDIC machine. * INTA = ICHAR( CA ) INTB = ICHAR( CB ) * IF( ZCODE.EQ.90 .OR. ZCODE.EQ.122 ) THEN * * ASCII is assumed - ZCODE is the ASCII code of either lower or * upper case 'Z'. * IF( INTA.GE.97 .AND. INTA.LE.122 ) INTA = INTA - 32 IF( INTB.GE.97 .AND. INTB.LE.122 ) INTB = INTB - 32 * ELSE IF( ZCODE.EQ.233 .OR. ZCODE.EQ.169 ) THEN * * EBCDIC is assumed - ZCODE is the EBCDIC code of either lower or * upper case 'Z'. * IF( INTA.GE.129 .AND. INTA.LE.137 .OR. $ INTA.GE.145 .AND. INTA.LE.153 .OR. $ INTA.GE.162 .AND. INTA.LE.169 ) INTA = INTA + 64 IF( INTB.GE.129 .AND. INTB.LE.137 .OR. $ INTB.GE.145 .AND. INTB.LE.153 .OR. $ INTB.GE.162 .AND. INTB.LE.169 ) INTB = INTB + 64 * ELSE IF( ZCODE.EQ.218 .OR. ZCODE.EQ.250 ) THEN * * ASCII is assumed, on Prime machines - ZCODE is the ASCII code * plus 128 of either lower or upper case 'Z'. * IF( INTA.GE.225 .AND. INTA.LE.250 ) INTA = INTA - 32 IF( INTB.GE.225 .AND. INTB.LE.250 ) INTB = INTB - 32 END IF LSAME = INTA.EQ.INTB * * RETURN * * End of LSAME * END SUBROUTINE DSPTRF( UPLO, N, AP, IPIV, INFO ) * * -- LAPACK routine (version 2.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * March 31, 1993 * * .. Scalar Arguments .. CHARACTER UPLO INTEGER INFO, N * .. * .. Array Arguments .. INTEGER IPIV( * ) DOUBLE PRECISION AP( * ) * .. * * Purpose * ======= * * DSPTRF computes the factorization of a real symmetric matrix A stored * in packed format using the Bunch-Kaufman diagonal pivoting method: * * A = U*D*U**T or A = L*D*L**T * * where U (or L) is a product of permutation and unit upper (lower) * triangular matrices, and D is symmetric and block diagonal with * 1-by-1 and 2-by-2 diagonal blocks. * * Arguments * ========= * * UPLO (input) CHARACTER*1 * = 'U': Upper triangle of A is stored; * = 'L': Lower triangle of A is stored. * * N (input) INTEGER * The order of the matrix A. N >= 0. * * AP (input/output) DOUBLE PRECISION array, dimension (N*(N+1)/2) * On entry, the upper or lower triangle of the symmetric matrix * A, packed columnwise in a linear array. The j-th column of A * is stored in the array AP as follows: * if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j; * if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n. * * On exit, the block diagonal matrix D and the multipliers used * to obtain the factor U or L, stored as a packed triangular * matrix overwriting A (see below for further details). * * IPIV (output) INTEGER array, dimension (N) * Details of the interchanges and the block structure of D. * If IPIV(k) > 0, then rows and columns k and IPIV(k) were * interchanged and D(k,k) is a 1-by-1 diagonal block. * If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and * columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) * is a 2-by-2 diagonal block. If UPLO = 'L' and IPIV(k) = * IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were * interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. * * INFO (output) INTEGER * = 0: successful exit * < 0: if INFO = -i, the i-th argument had an illegal value * > 0: if INFO = i, D(i,i) is exactly zero. The factorization * has been completed, but the block diagonal matrix D is * exactly singular, and division by zero will occur if it * is used to solve a system of equations. * * Further Details * =============== * * If UPLO = 'U', then A = U*D*U', where * U = P(n)*U(n)* ... *P(k)U(k)* ..., * i.e., U is a product of terms P(k)*U(k), where k decreases from n to * 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 * and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as * defined by IPIV(k), and U(k) is a unit upper triangular matrix, such * that if the diagonal block D(k) is of order s (s = 1 or 2), then * * ( I v 0 ) k-s * U(k) = ( 0 I 0 ) s * ( 0 0 I ) n-k * k-s s n-k * * If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). * If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), * and A(k,k), and v overwrites A(1:k-2,k-1:k). * * If UPLO = 'L', then A = L*D*L', where * L = P(1)*L(1)* ... *P(k)*L(k)* ..., * i.e., L is a product of terms P(k)*L(k), where k increases from 1 to * n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 * and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as * defined by IPIV(k), and L(k) is a unit lower triangular matrix, such * that if the diagonal block D(k) is of order s (s = 1 or 2), then * * ( I 0 0 ) k-1 * L(k) = ( 0 I 0 ) s * ( 0 v I ) n-k-s+1 * k-1 s n-k-s+1 * * If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). * If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), * and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1). * * ===================================================================== * * .. Parameters .. DOUBLE PRECISION ZERO, ONE PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) DOUBLE PRECISION EIGHT, SEVTEN PARAMETER ( EIGHT = 8.0D+0, SEVTEN = 17.0D+0 ) * .. * .. Local Scalars .. LOGICAL UPPER INTEGER IMAX, J, JMAX, K, KC, KK, KNC, KP, KPC, KSTEP, $ KX, NPP DOUBLE PRECISION ABSAKK, ALPHA, C, COLMAX, R1, R2, ROWMAX, S, T * .. * .. External Functions .. LOGICAL LSAME INTEGER IDAMAX EXTERNAL LSAME, IDAMAX * .. * .. External Subroutines .. EXTERNAL DLAEV2, DROT, DSCAL, DSPR, DSWAP, XERBLA * .. * .. Intrinsic Functions .. INTRINSIC ABS, MAX, SQRT * .. * .. Executable Statements .. * * Test the input parameters. * INFO = 0 UPPER = LSAME( UPLO, 'U' ) IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN INFO = -1 ELSE IF( N.LT.0 ) THEN INFO = -2 END IF IF( INFO.NE.0 ) THEN CALL XERBLA( 'DSPTRF', -INFO ) RETURN END IF * * Initialize ALPHA for use in choosing pivot block size. * ALPHA = ( ONE+SQRT( SEVTEN ) ) / EIGHT * IF( UPPER ) THEN * * Factorize A as U*D*U' using the upper triangle of A * * K is the main loop index, decreasing from N to 1 in steps of * 1 or 2 * K = N KC = ( N-1 )*N / 2 + 1 10 CONTINUE KNC = KC * * If K < 1, exit from loop * IF( K.LT.1 ) $ GO TO 70 KSTEP = 1 * * Determine rows and columns to be interchanged and whether * a 1-by-1 or 2-by-2 pivot block will be used * ABSAKK = ABS( AP( KC+K-1 ) ) * * IMAX is the row-index of the largest off-diagonal element in * column K, and COLMAX is its absolute value * IF( K.GT.1 ) THEN IMAX = IDAMAX( K-1, AP( KC ), 1 ) COLMAX = ABS( AP( KC+IMAX-1 ) ) ELSE COLMAX = ZERO END IF * IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN * * Column K is zero: set INFO and continue * IF( INFO.EQ.0 ) $ INFO = K KP = K ELSE IF( ABSAKK.GE.ALPHA*COLMAX ) THEN * * no interchange, use 1-by-1 pivot block * KP = K ELSE * * JMAX is the column-index of the largest off-diagonal * element in row IMAX, and ROWMAX is its absolute value * ROWMAX = ZERO JMAX = IMAX KX = IMAX*( IMAX+1 ) / 2 + IMAX DO 20 J = IMAX + 1, K IF( ABS( AP( KX ) ).GT.ROWMAX ) THEN ROWMAX = ABS( AP( KX ) ) JMAX = J END IF KX = KX + J 20 CONTINUE KPC = ( IMAX-1 )*IMAX / 2 + 1 IF( IMAX.GT.1 ) THEN JMAX = IDAMAX( IMAX-1, AP( KPC ), 1 ) ROWMAX = MAX( ROWMAX, ABS( AP( KPC+JMAX-1 ) ) ) END IF * IF( ABSAKK.GE.ALPHA*COLMAX*( COLMAX / ROWMAX ) ) THEN * * no interchange, use 1-by-1 pivot block * KP = K ELSE IF( ABS( AP( KPC+IMAX-1 ) ).GE.ALPHA*ROWMAX ) THEN * * interchange rows and columns K and IMAX, use 1-by-1 * pivot block * KP = IMAX ELSE * * interchange rows and columns K-1 and IMAX, use 2-by-2 * pivot block * KP = IMAX KSTEP = 2 END IF END IF * KK = K - KSTEP + 1 IF( KSTEP.EQ.2 ) $ KNC = KNC - K + 1 IF( KP.NE.KK ) THEN * * Interchange rows and columns KK and KP in the leading * submatrix A(1:k,1:k) * CALL DSWAP( KP-1, AP( KNC ), 1, AP( KPC ), 1 ) KX = KPC + KP - 1 DO 30 J = KP + 1, KK - 1 KX = KX + J - 1 T = AP( KNC+J-1 ) AP( KNC+J-1 ) = AP( KX ) AP( KX ) = T 30 CONTINUE T = AP( KNC+KK-1 ) AP( KNC+KK-1 ) = AP( KPC+KP-1 ) AP( KPC+KP-1 ) = T IF( KSTEP.EQ.2 ) THEN T = AP( KC+K-2 ) AP( KC+K-2 ) = AP( KC+KP-1 ) AP( KC+KP-1 ) = T END IF END IF * * Update the leading submatrix * IF( KSTEP.EQ.1 ) THEN * * 1-by-1 pivot block D(k): column k now holds * * W(k) = U(k)*D(k) * * where U(k) is the k-th column of U * * Perform a rank-1 update of A(1:k-1,1:k-1) as * * A := A - U(k)*D(k)*U(k)' = A - W(k)*1/D(k)*W(k)' * R1 = ONE / AP( KC+K-1 ) CALL DSPR( UPLO, K-1, -R1, AP( KC ), 1, AP ) * * Store U(k) in column k * CALL DSCAL( K-1, R1, AP( KC ), 1 ) ELSE * * 2-by-2 pivot block D(k): columns k and k-1 now hold * * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k) * * where U(k) and U(k-1) are the k-th and (k-1)-th columns * of U * * Perform a rank-2 update of A(1:k-2,1:k-2) as * * A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )' * = A - ( W(k-1) W(k) )*inv(D(k))*( W(k-1) W(k) )' * * Convert this to two rank-1 updates by using the eigen- * decomposition of D(k) * CALL DLAEV2( AP( KC-1 ), AP( KC+K-2 ), AP( KC+K-1 ), R1, $ R2, C, S ) R1 = ONE / R1 R2 = ONE / R2 CALL DROT( K-2, AP( KNC ), 1, AP( KC ), 1, C, S ) CALL DSPR( UPLO, K-2, -R1, AP( KNC ), 1, AP ) CALL DSPR( UPLO, K-2, -R2, AP( KC ), 1, AP ) * * Store U(k) and U(k-1) in columns k and k-1 * CALL DSCAL( K-2, R1, AP( KNC ), 1 ) CALL DSCAL( K-2, R2, AP( KC ), 1 ) CALL DROT( K-2, AP( KNC ), 1, AP( KC ), 1, C, -S ) END IF END IF * * Store details of the interchanges in IPIV * IF( KSTEP.EQ.1 ) THEN IPIV( K ) = KP ELSE IPIV( K ) = -KP IPIV( K-1 ) = -KP END IF * * Decrease K and return to the start of the main loop * K = K - KSTEP KC = KNC - K GO TO 10 * ELSE * * Factorize A as L*D*L' using the lower triangle of A * * K is the main loop index, increasing from 1 to N in steps of * 1 or 2 * K = 1 KC = 1 NPP = N*( N+1 ) / 2 40 CONTINUE KNC = KC * * If K > N, exit from loop * IF( K.GT.N ) $ GO TO 70 KSTEP = 1 * * Determine rows and columns to be interchanged and whether * a 1-by-1 or 2-by-2 pivot block will be used * ABSAKK = ABS( AP( KC ) ) * * IMAX is the row-index of the largest off-diagonal element in * column K, and COLMAX is its absolute value * IF( K.LT.N ) THEN IMAX = K + IDAMAX( N-K, AP( KC+1 ), 1 ) COLMAX = ABS( AP( KC+IMAX-K ) ) ELSE COLMAX = ZERO END IF * IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN * * Column K is zero: set INFO and continue * IF( INFO.EQ.0 ) $ INFO = K KP = K ELSE IF( ABSAKK.GE.ALPHA*COLMAX ) THEN * * no interchange, use 1-by-1 pivot block * KP = K ELSE * * JMAX is the column-index of the largest off-diagonal * element in row IMAX, and ROWMAX is its absolute value * ROWMAX = ZERO KX = KC + IMAX - K DO 50 J = K, IMAX - 1 IF( ABS( AP( KX ) ).GT.ROWMAX ) THEN ROWMAX = ABS( AP( KX ) ) JMAX = J END IF KX = KX + N - J 50 CONTINUE KPC = NPP - ( N-IMAX+1 )*( N-IMAX+2 ) / 2 + 1 IF( IMAX.LT.N ) THEN JMAX = IMAX + IDAMAX( N-IMAX, AP( KPC+1 ), 1 ) ROWMAX = MAX( ROWMAX, ABS( AP( KPC+JMAX-IMAX ) ) ) END IF * IF( ABSAKK.GE.ALPHA*COLMAX*( COLMAX / ROWMAX ) ) THEN * * no interchange, use 1-by-1 pivot block * KP = K ELSE IF( ABS( AP( KPC ) ).GE.ALPHA*ROWMAX ) THEN * * interchange rows and columns K and IMAX, use 1-by-1 * pivot block * KP = IMAX ELSE * * interchange rows and columns K+1 and IMAX, use 2-by-2 * pivot block * KP = IMAX KSTEP = 2 END IF END IF * KK = K + KSTEP - 1 IF( KSTEP.EQ.2 ) $ KNC = KNC + N - K + 1 IF( KP.NE.KK ) THEN * * Interchange rows and columns KK and KP in the trailing * submatrix A(k:n,k:n) * IF( KP.LT.N ) $ CALL DSWAP( N-KP, AP( KNC+KP-KK+1 ), 1, AP( KPC+1 ), $ 1 ) KX = KNC + KP - KK DO 60 J = KK + 1, KP - 1 KX = KX + N - J + 1 T = AP( KNC+J-KK ) AP( KNC+J-KK ) = AP( KX ) AP( KX ) = T 60 CONTINUE T = AP( KNC ) AP( KNC ) = AP( KPC ) AP( KPC ) = T IF( KSTEP.EQ.2 ) THEN T = AP( KC+1 ) AP( KC+1 ) = AP( KC+KP-K ) AP( KC+KP-K ) = T END IF END IF * * Update the trailing submatrix * IF( KSTEP.EQ.1 ) THEN * * 1-by-1 pivot block D(k): column k now holds * * W(k) = L(k)*D(k) * * where L(k) is the k-th column of L * IF( K.LT.N ) THEN * * Perform a rank-1 update of A(k+1:n,k+1:n) as * * A := A - L(k)*D(k)*L(k)' = A - W(k)*(1/D(k))*W(k)' * R1 = ONE / AP( KC ) CALL DSPR( UPLO, N-K, -R1, AP( KC+1 ), 1, $ AP( KC+N-K+1 ) ) * * Store L(k) in column K * CALL DSCAL( N-K, R1, AP( KC+1 ), 1 ) END IF ELSE * * 2-by-2 pivot block D(k): columns K and K+1 now hold * * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k) * * where L(k) and L(k+1) are the k-th and (k+1)-th columns * of L * IF( K.LT.N-1 ) THEN * * Perform a rank-2 update of A(k+2:n,k+2:n) as * * A := A - ( L(k) L(k+1) )*D(k)*( L(k) L(k+1) )' * = A - ( W(k) W(k+1) )*inv(D(k))*( W(k) W(k+1) )' * * Convert this to two rank-1 updates by using the eigen- * decomposition of D(k) * CALL DLAEV2( AP( KC ), AP( KC+1 ), AP( KNC ), R1, R2, $ C, S ) R1 = ONE / R1 R2 = ONE / R2 CALL DROT( N-K-1, AP( KC+2 ), 1, AP( KNC+1 ), 1, C, $ S ) CALL DSPR( UPLO, N-K-1, -R1, AP( KC+2 ), 1, $ AP( KNC+N-K ) ) CALL DSPR( UPLO, N-K-1, -R2, AP( KNC+1 ), 1, $ AP( KNC+N-K ) ) * * Store L(k) and L(k+1) in columns k and k+1 * CALL DSCAL( N-K-1, R1, AP( KC+2 ), 1 ) CALL DSCAL( N-K-1, R2, AP( KNC+1 ), 1 ) CALL DROT( N-K-1, AP( KC+2 ), 1, AP( KNC+1 ), 1, C, $ -S ) END IF END IF END IF * * Store details of the interchanges in IPIV * IF( KSTEP.EQ.1 ) THEN IPIV( K ) = KP ELSE IPIV( K ) = -KP IPIV( K+1 ) = -KP END IF * * Increase K and return to the start of the main loop * K = K + KSTEP KC = KNC + N - K + 2 GO TO 40 * END IF * 70 CONTINUE RETURN * * End of DSPTRF * END SUBROUTINE DSPTRS( UPLO, N, NRHS, AP, IPIV, B, LDB, INFO ) * * -- LAPACK routine (version 2.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * March 31, 1993 * * .. Scalar Arguments .. CHARACTER UPLO INTEGER INFO, LDB, N, NRHS * .. * .. Array Arguments .. INTEGER IPIV( * ) DOUBLE PRECISION AP( * ), B( LDB, * ) * .. * * Purpose * ======= * * DSPTRS solves a system of linear equations A*X = B with a real * symmetric matrix A stored in packed format using the factorization * A = U*D*U**T or A = L*D*L**T computed by DSPTRF. * * Arguments * ========= * * UPLO (input) CHARACTER*1 * Specifies whether the details of the factorization are stored * as an upper or lower triangular matrix. * = 'U': Upper triangular, form is A = U*D*U**T; * = 'L': Lower triangular, form is A = L*D*L**T. * * 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. * * AP (input) DOUBLE PRECISION array, dimension (N*(N+1)/2) * The block diagonal matrix D and the multipliers used to * obtain the factor U or L as computed by DSPTRF, stored as a * packed triangular matrix. * * IPIV (input) INTEGER array, dimension (N) * Details of the interchanges and the block structure of D * as determined by DSPTRF. * * 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 * * ===================================================================== * * .. Parameters .. DOUBLE PRECISION ONE PARAMETER ( ONE = 1.0D+0 ) * .. * .. Local Scalars .. LOGICAL UPPER INTEGER J, K, KC, KP DOUBLE PRECISION AK, AKM1, AKM1K, BK, BKM1, DENOM * .. * .. External Functions .. LOGICAL LSAME EXTERNAL LSAME * .. * .. External Subroutines .. EXTERNAL DGEMV, DGER, DSCAL, DSWAP, XERBLA * .. * .. Intrinsic Functions .. INTRINSIC MAX * .. * .. Executable Statements .. * INFO = 0 UPPER = LSAME( UPLO, 'U' ) IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN INFO = -1 ELSE IF( N.LT.0 ) THEN INFO = -2 ELSE IF( NRHS.LT.0 ) THEN INFO = -3 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN INFO = -7 END IF IF( INFO.NE.0 ) THEN CALL XERBLA( 'DSPTRS', -INFO ) RETURN END IF * * Quick return if possible * IF( N.EQ.0 .OR. NRHS.EQ.0 ) $ RETURN * IF( UPPER ) THEN * * Solve A*X = B, where A = U*D*U'. * * First solve U*D*X = B, overwriting B with X. * * K is the main loop index, decreasing from N to 1 in steps of * 1 or 2, depending on the size of the diagonal blocks. * K = N KC = N*( N+1 ) / 2 + 1 10 CONTINUE * * If K < 1, exit from loop. * IF( K.LT.1 ) $ GO TO 30 * KC = KC - K IF( IPIV( K ).GT.0 ) THEN * * 1 x 1 diagonal block * * Interchange rows K and IPIV(K). * KP = IPIV( K ) IF( KP.NE.K ) $ CALL DSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB ) * * Multiply by inv(U(K)), where U(K) is the transformation * stored in column K of A. * CALL DGER( K-1, NRHS, -ONE, AP( KC ), 1, B( K, 1 ), LDB, $ B( 1, 1 ), LDB ) * * Multiply by the inverse of the diagonal block. * CALL DSCAL( NRHS, ONE / AP( KC+K-1 ), B( K, 1 ), LDB ) K = K - 1 ELSE * * 2 x 2 diagonal block * * Interchange rows K-1 and -IPIV(K). * KP = -IPIV( K ) IF( KP.NE.K-1 ) $ CALL DSWAP( NRHS, B( K-1, 1 ), LDB, B( KP, 1 ), LDB ) * * Multiply by inv(U(K)), where U(K) is the transformation * stored in columns K-1 and K of A. * CALL DGER( K-2, NRHS, -ONE, AP( KC ), 1, B( K, 1 ), LDB, $ B( 1, 1 ), LDB ) CALL DGER( K-2, NRHS, -ONE, AP( KC-( K-1 ) ), 1, $ B( K-1, 1 ), LDB, B( 1, 1 ), LDB ) * * Multiply by the inverse of the diagonal block. * AKM1K = AP( KC+K-2 ) AKM1 = AP( KC-1 ) / AKM1K AK = AP( KC+K-1 ) / AKM1K DENOM = AKM1*AK - ONE DO 20 J = 1, NRHS BKM1 = B( K-1, J ) / AKM1K BK = B( K, J ) / AKM1K B( K-1, J ) = ( AK*BKM1-BK ) / DENOM B( K, J ) = ( AKM1*BK-BKM1 ) / DENOM 20 CONTINUE KC = KC - K + 1 K = K - 2 END IF * GO TO 10 30 CONTINUE * * Next solve U'*X = B, overwriting B with X. * * K is the main loop index, increasing from 1 to N in steps of * 1 or 2, depending on the size of the diagonal blocks. * K = 1 KC = 1 40 CONTINUE * * If K > N, exit from loop. * IF( K.GT.N ) $ GO TO 50 * IF( IPIV( K ).GT.0 ) THEN * * 1 x 1 diagonal block * * Multiply by inv(U'(K)), where U(K) is the transformation * stored in column K of A. * CALL DGEMV( 'Transpose', K-1, NRHS, -ONE, B, LDB, AP( KC ), $ 1, ONE, B( K, 1 ), LDB ) * * Interchange rows K and IPIV(K). * KP = IPIV( K ) IF( KP.NE.K ) $ CALL DSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB ) KC = KC + K K = K + 1 ELSE * * 2 x 2 diagonal block * * Multiply by inv(U'(K+1)), where U(K+1) is the transformation * stored in columns K and K+1 of A. * CALL DGEMV( 'Transpose', K-1, NRHS, -ONE, B, LDB, AP( KC ), $ 1, ONE, B( K, 1 ), LDB ) CALL DGEMV( 'Transpose', K-1, NRHS, -ONE, B, LDB, $ AP( KC+K ), 1, ONE, B( K+1, 1 ), LDB ) * * Interchange rows K and -IPIV(K). * KP = -IPIV( K ) IF( KP.NE.K ) $ CALL DSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB ) KC = KC + 2*K + 1 K = K + 2 END IF * GO TO 40 50 CONTINUE * ELSE * * Solve A*X = B, where A = L*D*L'. * * First solve L*D*X = B, overwriting B with X. * * K is the main loop index, increasing from 1 to N in steps of * 1 or 2, depending on the size of the diagonal blocks. * K = 1 KC = 1 60 CONTINUE * * If K > N, exit from loop. * IF( K.GT.N ) $ GO TO 80 * IF( IPIV( K ).GT.0 ) THEN * * 1 x 1 diagonal block * * Interchange rows K and IPIV(K). * KP = IPIV( K ) IF( KP.NE.K ) $ CALL DSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB ) * * Multiply by inv(L(K)), where L(K) is the transformation * stored in column K of A. * IF( K.LT.N ) $ CALL DGER( N-K, NRHS, -ONE, AP( KC+1 ), 1, B( K, 1 ), $ LDB, B( K+1, 1 ), LDB ) * * Multiply by the inverse of the diagonal block. * CALL DSCAL( NRHS, ONE / AP( KC ), B( K, 1 ), LDB ) KC = KC + N - K + 1 K = K + 1 ELSE * * 2 x 2 diagonal block * * Interchange rows K+1 and -IPIV(K). * KP = -IPIV( K ) IF( KP.NE.K+1 ) $ CALL DSWAP( NRHS, B( K+1, 1 ), LDB, B( KP, 1 ), LDB ) * * Multiply by inv(L(K)), where L(K) is the transformation * stored in columns K and K+1 of A. * IF( K.LT.N-1 ) THEN CALL DGER( N-K-1, NRHS, -ONE, AP( KC+2 ), 1, B( K, 1 ), $ LDB, B( K+2, 1 ), LDB ) CALL DGER( N-K-1, NRHS, -ONE, AP( KC+N-K+2 ), 1, $ B( K+1, 1 ), LDB, B( K+2, 1 ), LDB ) END IF * * Multiply by the inverse of the diagonal block. * AKM1K = AP( KC+1 ) AKM1 = AP( KC ) / AKM1K AK = AP( KC+N-K+1 ) / AKM1K DENOM = AKM1*AK - ONE DO 70 J = 1, NRHS BKM1 = B( K, J ) / AKM1K BK = B( K+1, J ) / AKM1K B( K, J ) = ( AK*BKM1-BK ) / DENOM B( K+1, J ) = ( AKM1*BK-BKM1 ) / DENOM 70 CONTINUE KC = KC + 2*( N-K ) + 1 K = K + 2 END IF * GO TO 60 80 CONTINUE * * Next solve L'*X = B, overwriting B with X. * * K is the main loop index, decreasing from N to 1 in steps of * 1 or 2, depending on the size of the diagonal blocks. * K = N KC = N*( N+1 ) / 2 + 1 90 CONTINUE * * If K < 1, exit from loop. * IF( K.LT.1 ) $ GO TO 100 * KC = KC - ( N-K+1 ) IF( IPIV( K ).GT.0 ) THEN * * 1 x 1 diagonal block * * Multiply by inv(L'(K)), where L(K) is the transformation * stored in column K of A. * IF( K.LT.N ) $ CALL DGEMV( 'Transpose', N-K, NRHS, -ONE, B( K+1, 1 ), $ LDB, AP( KC+1 ), 1, ONE, B( K, 1 ), LDB ) * * Interchange rows K and IPIV(K). * KP = IPIV( K ) IF( KP.NE.K ) $ CALL DSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB ) K = K - 1 ELSE * * 2 x 2 diagonal block * * Multiply by inv(L'(K-1)), where L(K-1) is the transformation * stored in columns K-1 and K of A. * IF( K.LT.N ) THEN CALL DGEMV( 'Transpose', N-K, NRHS, -ONE, B( K+1, 1 ), $ LDB, AP( KC+1 ), 1, ONE, B( K, 1 ), LDB ) CALL DGEMV( 'Transpose', N-K, NRHS, -ONE, B( K+1, 1 ), $ LDB, AP( KC-( N-K ) ), 1, ONE, B( K-1, 1 ), $ LDB ) END IF * * Interchange rows K and -IPIV(K). * KP = -IPIV( K ) IF( KP.NE.K ) $ CALL DSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB ) KC = KC - ( N-K+2 ) K = K - 2 END IF * GO TO 90 100 CONTINUE END IF * RETURN * * End of DSPTRS * END SUBROUTINE XERBLA( SRNAME, INFO ) * * -- LAPACK auxiliary routine (version 2.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * September 30, 1994 * * .. Scalar Arguments .. CHARACTER*6 SRNAME INTEGER INFO * .. * * Purpose * ======= * * XERBLA is an error handler for the LAPACK routines. * It is called by an LAPACK routine if an input parameter has an * invalid value. A message is printed and execution stops. * * Installers may consider modifying the STOP statement in order to * call system-specific exception-handling facilities. * * Arguments * ========= * * SRNAME (input) CHARACTER*6 * The name of the routine which called XERBLA. * * INFO (input) INTEGER * The position of the invalid parameter in the parameter list * of the calling routine. * * ===================================================================== * * .. Executable Statements .. * WRITE( *, FMT = 9999 )SRNAME, INFO * STOP * 9999 FORMAT( ' ** On entry to ', A6, ' parameter number ', I2, ' had ', $ 'an illegal value' ) * * End of XERBLA * END SUBROUTINE DLAEV2( A, B, C, RT1, RT2, CS1, SN1 ) * * -- LAPACK auxiliary routine (version 2.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * October 31, 1992 * * .. Scalar Arguments .. DOUBLE PRECISION A, B, C, CS1, RT1, RT2, SN1 * .. * * Purpose * ======= * * DLAEV2 computes the eigendecomposition of a 2-by-2 symmetric matrix * [ A B ] * [ B C ]. * On return, RT1 is the eigenvalue of larger absolute value, RT2 is the * eigenvalue of smaller absolute value, and (CS1,SN1) is the unit right * eigenvector for RT1, giving the decomposition * * [ CS1 SN1 ] [ A B ] [ CS1 -SN1 ] = [ RT1 0 ] * [-SN1 CS1 ] [ B C ] [ SN1 CS1 ] [ 0 RT2 ]. * * Arguments * ========= * * A (input) DOUBLE PRECISION * The (1,1) element of the 2-by-2 matrix. * * B (input) DOUBLE PRECISION * The (1,2) element and the conjugate of the (2,1) element of * the 2-by-2 matrix. * * C (input) DOUBLE PRECISION * The (2,2) element of the 2-by-2 matrix. * * RT1 (output) DOUBLE PRECISION * The eigenvalue of larger absolute value. * * RT2 (output) DOUBLE PRECISION * The eigenvalue of smaller absolute value. * * CS1 (output) DOUBLE PRECISION * SN1 (output) DOUBLE PRECISION * The vector (CS1, SN1) is a unit right eigenvector for RT1. * * Further Details * =============== * * RT1 is accurate to a few ulps barring over/underflow. * * RT2 may be inaccurate if there is massive cancellation in the * determinant A*C-B*B; higher precision or correctly rounded or * correctly truncated arithmetic would be needed to compute RT2 * accurately in all cases. * * CS1 and SN1 are accurate to a few ulps barring over/underflow. * * Overflow is possible only if RT1 is within a factor of 5 of overflow. * Underflow is harmless if the input data is 0 or exceeds * underflow_threshold / macheps. * * ===================================================================== * * .. Parameters .. DOUBLE PRECISION ONE PARAMETER ( ONE = 1.0D0 ) DOUBLE PRECISION TWO PARAMETER ( TWO = 2.0D0 ) DOUBLE PRECISION ZERO PARAMETER ( ZERO = 0.0D0 ) DOUBLE PRECISION HALF PARAMETER ( HALF = 0.5D0 ) * .. * .. Local Scalars .. INTEGER SGN1, SGN2 DOUBLE PRECISION AB, ACMN, ACMX, ACS, ADF, CS, CT, DF, RT, SM, $ TB, TN * .. * .. Intrinsic Functions .. INTRINSIC ABS, SQRT * .. * .. Executable Statements .. * * Compute the eigenvalues * SM = A + C DF = A - C ADF = ABS( DF ) TB = B + B AB = ABS( TB ) IF( ABS( A ).GT.ABS( C ) ) THEN ACMX = A ACMN = C ELSE ACMX = C ACMN = A END IF IF( ADF.GT.AB ) THEN RT = ADF*SQRT( ONE+( AB / ADF )**2 ) ELSE IF( ADF.LT.AB ) THEN RT = AB*SQRT( ONE+( ADF / AB )**2 ) ELSE * * Includes case AB=ADF=0 * RT = AB*SQRT( TWO ) END IF IF( SM.LT.ZERO ) THEN RT1 = HALF*( SM-RT ) SGN1 = -1 * * Order of execution important. * To get fully accurate smaller eigenvalue, * next line needs to be executed in higher precision. * RT2 = ( ACMX / RT1 )*ACMN - ( B / RT1 )*B ELSE IF( SM.GT.ZERO ) THEN RT1 = HALF*( SM+RT ) SGN1 = 1 * * Order of execution important. * To get fully accurate smaller eigenvalue, * next line needs to be executed in higher precision. * RT2 = ( ACMX / RT1 )*ACMN - ( B / RT1 )*B ELSE * * Includes case RT1 = RT2 = 0 * RT1 = HALF*RT RT2 = -HALF*RT SGN1 = 1 END IF * * Compute the eigenvector * IF( DF.GE.ZERO ) THEN CS = DF + RT SGN2 = 1 ELSE CS = DF - RT SGN2 = -1 END IF ACS = ABS( CS ) IF( ACS.GT.AB ) THEN CT = -TB / CS SN1 = ONE / SQRT( ONE+CT*CT ) CS1 = CT*SN1 ELSE IF( AB.EQ.ZERO ) THEN CS1 = ONE SN1 = ZERO ELSE TN = -CS / TB CS1 = ONE / SQRT( ONE+TN*TN ) SN1 = TN*CS1 END IF END IF IF( SGN1.EQ.SGN2 ) THEN TN = CS1 CS1 = -SN1 SN1 = TN END IF RETURN * * End of DLAEV2 * END subroutine drot (n,dx,incx,dy,incy,c,s) c c applies a plane rotation. c jack dongarra, linpack, 3/11/78. c modified 12/3/93, array(1) declarations changed to array(*) c double precision dx(*),dy(*),dtemp,c,s integer i,incx,incy,ix,iy,n c if(n.le.0)return if(incx.eq.1.and.incy.eq.1)go to 20 c c code for unequal increments or equal increments not equal c to 1 c ix = 1 iy = 1 if(incx.lt.0)ix = (-n+1)*incx + 1 if(incy.lt.0)iy = (-n+1)*incy + 1 do 10 i = 1,n dtemp = c*dx(ix) + s*dy(iy) dy(iy) = c*dy(iy) - s*dx(ix) dx(ix) = dtemp ix = ix + incx iy = iy + incy 10 continue return c c code for both increments equal to 1 c 20 do 30 i = 1,n dtemp = c*dx(i) + s*dy(i) dy(i) = c*dy(i) - s*dx(i) dx(i) = dtemp 30 continue return end subroutine dscal(n,da,dx,incx) c c scales a vector by a constant. c uses unrolled loops for increment equal to one. c jack dongarra, linpack, 3/11/78. c modified 3/93 to return if incx .le. 0. c modified 12/3/93, array(1) declarations changed to array(*) c double precision da,dx(*) integer i,incx,m,mp1,n,nincx c if( n.le.0 .or. incx.le.0 )return if(incx.eq.1)go to 20 c c code for increment not equal to 1 c nincx = n*incx do 10 i = 1,nincx,incx dx(i) = da*dx(i) 10 continue return c c code for increment equal to 1 c c c clean-up loop c 20 m = mod(n,5) if( m .eq. 0 ) go to 40 do 30 i = 1,m dx(i) = da*dx(i) 30 continue if( n .lt. 5 ) return 40 mp1 = m + 1 do 50 i = mp1,n,5 dx(i) = da*dx(i) dx(i + 1) = da*dx(i + 1) dx(i + 2) = da*dx(i + 2) dx(i + 3) = da*dx(i + 3) dx(i + 4) = da*dx(i + 4) 50 continue return end SUBROUTINE DSPR ( UPLO, N, ALPHA, X, INCX, AP ) * .. Scalar Arguments .. DOUBLE PRECISION ALPHA INTEGER INCX, N CHARACTER*1 UPLO * .. Array Arguments .. DOUBLE PRECISION AP( * ), X( * ) * .. * * Purpose * ======= * * DSPR performs the symmetric rank 1 operation * * A := alpha*x*x' + A, * * where alpha is a real scalar, x is an n element vector and A is an * n by n symmetric matrix, supplied in packed form. * * Parameters * ========== * * UPLO - CHARACTER*1. * On entry, UPLO specifies whether the upper or lower * triangular part of the matrix A is supplied in the packed * array AP as follows: * * UPLO = 'U' or 'u' The upper triangular part of A is * supplied in AP. * * UPLO = 'L' or 'l' The lower triangular part of A is * supplied in AP. * * Unchanged on exit. * * N - INTEGER. * On entry, N specifies the order of the matrix A. * N must be at least zero. * Unchanged on exit. * * ALPHA - DOUBLE PRECISION. * On entry, ALPHA specifies the scalar alpha. * Unchanged on exit. * * X - DOUBLE PRECISION array of dimension at least * ( 1 + ( n - 1 )*abs( INCX ) ). * Before entry, the incremented array X must contain the n * element vector x. * Unchanged on exit. * * INCX - INTEGER. * On entry, INCX specifies the increment for the elements of * X. INCX must not be zero. * Unchanged on exit. * * AP - DOUBLE PRECISION array of DIMENSION at least * ( ( n*( n + 1 ) )/2 ). * Before entry with UPLO = 'U' or 'u', the array AP must * contain the upper triangular part of the symmetric matrix * packed sequentially, column by column, so that AP( 1 ) * contains a( 1, 1 ), AP( 2 ) and AP( 3 ) contain a( 1, 2 ) * and a( 2, 2 ) respectively, and so on. On exit, the array * AP is overwritten by the upper triangular part of the * updated matrix. * Before entry with UPLO = 'L' or 'l', the array AP must * contain the lower triangular part of the symmetric matrix * packed sequentially, column by column, so that AP( 1 ) * contains a( 1, 1 ), AP( 2 ) and AP( 3 ) contain a( 2, 1 ) * and a( 3, 1 ) respectively, and so on. On exit, the array * AP is overwritten by the lower triangular part of the * updated matrix. * * * Level 2 Blas routine. * * -- Written on 22-October-1986. * Jack Dongarra, Argonne National Lab. * Jeremy Du Croz, Nag Central Office. * Sven Hammarling, Nag Central Office. * Richard Hanson, Sandia National Labs. * * * .. Parameters .. DOUBLE PRECISION ZERO PARAMETER ( ZERO = 0.0D+0 ) * .. Local Scalars .. DOUBLE PRECISION TEMP INTEGER I, INFO, IX, J, JX, K, KK, KX * .. External Functions .. LOGICAL LSAME EXTERNAL LSAME * .. External Subroutines .. EXTERNAL XERBLA * .. * .. Executable Statements .. * * Test the input parameters. * INFO = 0 IF ( .NOT.LSAME( UPLO, 'U' ).AND. $ .NOT.LSAME( UPLO, 'L' ) )THEN INFO = 1 ELSE IF( N.LT.0 )THEN INFO = 2 ELSE IF( INCX.EQ.0 )THEN INFO = 5 END IF IF( INFO.NE.0 )THEN CALL XERBLA( 'DSPR ', INFO ) RETURN END IF * * Quick return if possible. * IF( ( N.EQ.0 ).OR.( ALPHA.EQ.ZERO ) ) $ RETURN * * Set the start point in X if the increment is not unity. * IF( INCX.LE.0 )THEN KX = 1 - ( N - 1 )*INCX ELSE IF( INCX.NE.1 )THEN KX = 1 END IF * * Start the operations. In this version the elements of the array AP * are accessed sequentially with one pass through AP. * KK = 1 IF( LSAME( UPLO, 'U' ) )THEN * * Form A when upper triangle is stored in AP. * IF( INCX.EQ.1 )THEN DO 20, J = 1, N IF( X( J ).NE.ZERO )THEN TEMP = ALPHA*X( J ) K = KK DO 10, I = 1, J AP( K ) = AP( K ) + X( I )*TEMP K = K + 1 10 CONTINUE END IF KK = KK + J 20 CONTINUE ELSE JX = KX DO 40, J = 1, N IF( X( JX ).NE.ZERO )THEN TEMP = ALPHA*X( JX ) IX = KX DO 30, K = KK, KK + J - 1 AP( K ) = AP( K ) + X( IX )*TEMP IX = IX + INCX 30 CONTINUE END IF JX = JX + INCX KK = KK + J 40 CONTINUE END IF ELSE * * Form A when lower triangle is stored in AP. * IF( INCX.EQ.1 )THEN DO 60, J = 1, N IF( X( J ).NE.ZERO )THEN TEMP = ALPHA*X( J ) K = KK DO 50, I = J, N AP( K ) = AP( K ) + X( I )*TEMP K = K + 1 50 CONTINUE END IF KK = KK + N - J + 1 60 CONTINUE ELSE JX = KX DO 80, J = 1, N IF( X( JX ).NE.ZERO )THEN TEMP = ALPHA*X( JX ) IX = JX DO 70, K = KK, KK + N - J AP( K ) = AP( K ) + X( IX )*TEMP IX = IX + INCX 70 CONTINUE END IF JX = JX + INCX KK = KK + N - J + 1 80 CONTINUE END IF END IF * RETURN * * End of DSPR . * END subroutine dswap (n,dx,incx,dy,incy) c c interchanges two vectors. c uses unrolled loops for increments equal one. c jack dongarra, linpack, 3/11/78. c modified 12/3/93, array(1) declarations changed to array(*) c double precision dx(*),dy(*),dtemp integer i,incx,incy,ix,iy,m,mp1,n c if(n.le.0)return if(incx.eq.1.and.incy.eq.1)go to 20 c c code for unequal increments or equal increments not equal c to 1 c ix = 1 iy = 1 if(incx.lt.0)ix = (-n+1)*incx + 1 if(incy.lt.0)iy = (-n+1)*incy + 1 do 10 i = 1,n dtemp = dx(ix) dx(ix) = dy(iy) dy(iy) = dtemp ix = ix + incx iy = iy + incy 10 continue return c c code for both increments equal to 1 c c c clean-up loop c 20 m = mod(n,3) if( m .eq. 0 ) go to 40 do 30 i = 1,m dtemp = dx(i) dx(i) = dy(i) dy(i) = dtemp 30 continue if( n .lt. 3 ) return 40 mp1 = m + 1 do 50 i = mp1,n,3 dtemp = dx(i) dx(i) = dy(i) dy(i) = dtemp dtemp = dx(i + 1) dx(i + 1) = dy(i + 1) dy(i + 1) = dtemp dtemp = dx(i + 2) dx(i + 2) = dy(i + 2) dy(i + 2) = dtemp 50 continue return end integer function idamax(n,dx,incx) c c finds the index of element having max. absolute value. c jack dongarra, linpack, 3/11/78. c modified 3/93 to return if incx .le. 0. c modified 12/3/93, array(1) declarations changed to array(*) c double precision dx(*),dmax integer i,incx,ix,n c idamax = 0 if( n.lt.1 .or. incx.le.0 ) return idamax = 1 if(n.eq.1)return if(incx.eq.1)go to 20 c c code for increment not equal to 1 c ix = 1 dmax = dabs(dx(1)) ix = ix + incx do 10 i = 2,n if(dabs(dx(ix)).le.dmax) go to 5 idamax = i dmax = dabs(dx(ix)) 5 ix = ix + incx 10 continue return c c code for increment equal to 1 c 20 dmax = dabs(dx(1)) do 30 i = 2,n if(dabs(dx(i)).le.dmax) go to 30 idamax = i dmax = dabs(dx(i)) 30 continue return end SUBROUTINE DGEMV ( TRANS, M, N, ALPHA, A, LDA, X, INCX, $ BETA, Y, INCY ) * .. Scalar Arguments .. DOUBLE PRECISION ALPHA, BETA INTEGER INCX, INCY, LDA, M, N CHARACTER*1 TRANS * .. Array Arguments .. DOUBLE PRECISION A( LDA, * ), X( * ), Y( * ) * .. * * Purpose * ======= * * DGEMV performs one of the matrix-vector operations * * y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y, * * where alpha and beta are scalars, x and y are vectors and A is an * m by n matrix. * * Parameters * ========== * * TRANS - CHARACTER*1. * On entry, TRANS specifies the operation to be performed as * follows: * * TRANS = 'N' or 'n' y := alpha*A*x + beta*y. * * TRANS = 'T' or 't' y := alpha*A'*x + beta*y. * * TRANS = 'C' or 'c' y := alpha*A'*x + beta*y. * * Unchanged on exit. * * M - INTEGER. * On entry, M specifies the number of rows of the matrix A. * M must be at least zero. * Unchanged on exit. * * N - INTEGER. * On entry, N specifies the number of columns of the matrix A. * N must be at least zero. * Unchanged on exit. * * ALPHA - DOUBLE PRECISION. * On entry, ALPHA specifies the scalar alpha. * Unchanged on exit. * * A - DOUBLE PRECISION array of DIMENSION ( LDA, n ). * Before entry, the leading m by n part of the array A must * contain the matrix of coefficients. * Unchanged on exit. * * LDA - INTEGER. * On entry, LDA specifies the first dimension of A as declared * in the calling (sub) program. LDA must be at least * max( 1, m ). * Unchanged on exit. * * X - DOUBLE PRECISION array of DIMENSION at least * ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n' * and at least * ( 1 + ( m - 1 )*abs( INCX ) ) otherwise. * Before entry, the incremented array X must contain the * vector x. * Unchanged on exit. * * INCX - INTEGER. * On entry, INCX specifies the increment for the elements of * X. INCX must not be zero. * Unchanged on exit. * * BETA - DOUBLE PRECISION. * On entry, BETA specifies the scalar beta. When BETA is * supplied as zero then Y need not be set on input. * Unchanged on exit. * * Y - DOUBLE PRECISION array of DIMENSION at least * ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n' * and at least * ( 1 + ( n - 1 )*abs( INCY ) ) otherwise. * Before entry with BETA non-zero, the incremented array Y * must contain the vector y. On exit, Y is overwritten by the * updated vector y. * * INCY - INTEGER. * On entry, INCY specifies the increment for the elements of * Y. INCY must not be zero. * Unchanged on exit. * * * Level 2 Blas routine. * * -- Written on 22-October-1986. * Jack Dongarra, Argonne National Lab. * Jeremy Du Croz, Nag Central Office. * Sven Hammarling, Nag Central Office. * Richard Hanson, Sandia National Labs. * * * .. Parameters .. DOUBLE PRECISION ONE , ZERO PARAMETER ( ONE = 1.0D+0, ZERO = 0.0D+0 ) * .. Local Scalars .. DOUBLE PRECISION TEMP INTEGER I, INFO, IX, IY, J, JX, JY, KX, KY, LENX, LENY * .. External Functions .. LOGICAL LSAME EXTERNAL LSAME * .. External Subroutines .. EXTERNAL XERBLA * .. Intrinsic Functions .. INTRINSIC MAX * .. * .. Executable Statements .. * * Test the input parameters. * INFO = 0 IF ( .NOT.LSAME( TRANS, 'N' ).AND. $ .NOT.LSAME( TRANS, 'T' ).AND. $ .NOT.LSAME( TRANS, 'C' ) )THEN INFO = 1 ELSE IF( M.LT.0 )THEN INFO = 2 ELSE IF( N.LT.0 )THEN INFO = 3 ELSE IF( LDA.LT.MAX( 1, M ) )THEN INFO = 6 ELSE IF( INCX.EQ.0 )THEN INFO = 8 ELSE IF( INCY.EQ.0 )THEN INFO = 11 END IF IF( INFO.NE.0 )THEN CALL XERBLA( 'DGEMV ', INFO ) RETURN END IF * * Quick return if possible. * IF( ( M.EQ.0 ).OR.( N.EQ.0 ).OR. $ ( ( ALPHA.EQ.ZERO ).AND.( BETA.EQ.ONE ) ) ) $ RETURN * * Set LENX and LENY, the lengths of the vectors x and y, and set * up the start points in X and Y. * IF( LSAME( TRANS, 'N' ) )THEN LENX = N LENY = M ELSE LENX = M LENY = N END IF IF( INCX.GT.0 )THEN KX = 1 ELSE KX = 1 - ( LENX - 1 )*INCX END IF IF( INCY.GT.0 )THEN KY = 1 ELSE KY = 1 - ( LENY - 1 )*INCY END IF * * Start the operations. In this version the elements of A are * accessed sequentially with one pass through A. * * First form y := beta*y. * IF( BETA.NE.ONE )THEN IF( INCY.EQ.1 )THEN IF( BETA.EQ.ZERO )THEN DO 10, I = 1, LENY Y( I ) = ZERO 10 CONTINUE ELSE DO 20, I = 1, LENY Y( I ) = BETA*Y( I ) 20 CONTINUE END IF ELSE IY = KY IF( BETA.EQ.ZERO )THEN DO 30, I = 1, LENY Y( IY ) = ZERO IY = IY + INCY 30 CONTINUE ELSE DO 40, I = 1, LENY Y( IY ) = BETA*Y( IY ) IY = IY + INCY 40 CONTINUE END IF END IF END IF IF( ALPHA.EQ.ZERO ) $ RETURN IF( LSAME( TRANS, 'N' ) )THEN * * Form y := alpha*A*x + y. * JX = KX IF( INCY.EQ.1 )THEN DO 60, J = 1, N IF( X( JX ).NE.ZERO )THEN TEMP = ALPHA*X( JX ) DO 50, I = 1, M Y( I ) = Y( I ) + TEMP*A( I, J ) 50 CONTINUE END IF JX = JX + INCX 60 CONTINUE ELSE DO 80, J = 1, N IF( X( JX ).NE.ZERO )THEN TEMP = ALPHA*X( JX ) IY = KY DO 70, I = 1, M Y( IY ) = Y( IY ) + TEMP*A( I, J ) IY = IY + INCY 70 CONTINUE END IF JX = JX + INCX 80 CONTINUE END IF ELSE * * Form y := alpha*A'*x + y. * JY = KY IF( INCX.EQ.1 )THEN DO 100, J = 1, N TEMP = ZERO DO 90, I = 1, M TEMP = TEMP + A( I, J )*X( I ) 90 CONTINUE Y( JY ) = Y( JY ) + ALPHA*TEMP JY = JY + INCY 100 CONTINUE ELSE DO 120, J = 1, N TEMP = ZERO IX = KX DO 110, I = 1, M TEMP = TEMP + A( I, J )*X( IX ) IX = IX + INCX 110 CONTINUE Y( JY ) = Y( JY ) + ALPHA*TEMP JY = JY + INCY 120 CONTINUE END IF END IF * RETURN * * End of DGEMV . * END SUBROUTINE DGER ( M, N, ALPHA, X, INCX, Y, INCY, A, LDA ) * .. Scalar Arguments .. DOUBLE PRECISION ALPHA INTEGER INCX, INCY, LDA, M, N * .. Array Arguments .. DOUBLE PRECISION A( LDA, * ), X( * ), Y( * ) * .. * * Purpose * ======= * * DGER performs the rank 1 operation * * A := alpha*x*y' + A, * * where alpha is a scalar, x is an m element vector, y is an n element * vector and A is an m by n matrix. * * Parameters * ========== * * M - INTEGER. * On entry, M specifies the number of rows of the matrix A. * M must be at least zero. * Unchanged on exit. * * N - INTEGER. * On entry, N specifies the number of columns of the matrix A. * N must be at least zero. * Unchanged on exit. * * ALPHA - DOUBLE PRECISION. * On entry, ALPHA specifies the scalar alpha. * Unchanged on exit. * * X - DOUBLE PRECISION array of dimension at least * ( 1 + ( m - 1 )*abs( INCX ) ). * Before entry, the incremented array X must contain the m * element vector x. * Unchanged on exit. * * INCX - INTEGER. * On entry, INCX specifies the increment for the elements of * X. INCX must not be zero. * Unchanged on exit. * * Y - DOUBLE PRECISION array of dimension at least * ( 1 + ( n - 1 )*abs( INCY ) ). * Before entry, the incremented array Y must contain the n * element vector y. * Unchanged on exit. * * INCY - INTEGER. * On entry, INCY specifies the increment for the elements of * Y. INCY must not be zero. * Unchanged on exit. * * A - DOUBLE PRECISION array of DIMENSION ( LDA, n ). * Before entry, the leading m by n part of the array A must * contain the matrix of coefficients. On exit, A is * overwritten by the updated matrix. * * LDA - INTEGER. * On entry, LDA specifies the first dimension of A as declared * in the calling (sub) program. LDA must be at least * max( 1, m ). * Unchanged on exit. * * * Level 2 Blas routine. * * -- Written on 22-October-1986. * Jack Dongarra, Argonne National Lab. * Jeremy Du Croz, Nag Central Office. * Sven Hammarling, Nag Central Office. * Richard Hanson, Sandia National Labs. * * * .. Parameters .. DOUBLE PRECISION ZERO PARAMETER ( ZERO = 0.0D+0 ) * .. Local Scalars .. DOUBLE PRECISION TEMP INTEGER I, INFO, IX, J, JY, KX * .. External Subroutines .. EXTERNAL XERBLA * .. Intrinsic Functions .. INTRINSIC MAX * .. * .. Executable Statements .. * * Test the input parameters. * INFO = 0 IF ( M.LT.0 )THEN INFO = 1 ELSE IF( N.LT.0 )THEN INFO = 2 ELSE IF( INCX.EQ.0 )THEN INFO = 5 ELSE IF( INCY.EQ.0 )THEN INFO = 7 ELSE IF( LDA.LT.MAX( 1, M ) )THEN INFO = 9 END IF IF( INFO.NE.0 )THEN CALL XERBLA( 'DGER ', INFO ) RETURN END IF * * Quick return if possible. * IF( ( M.EQ.0 ).OR.( N.EQ.0 ).OR.( ALPHA.EQ.ZERO ) ) $ RETURN * * Start the operations. In this version the elements of A are * accessed sequentially with one pass through A. * IF( INCY.GT.0 )THEN JY = 1 ELSE JY = 1 - ( N - 1 )*INCY END IF IF( INCX.EQ.1 )THEN DO 20, J = 1, N IF( Y( JY ).NE.ZERO )THEN TEMP = ALPHA*Y( JY ) DO 10, I = 1, M A( I, J ) = A( I, J ) + X( I )*TEMP 10 CONTINUE END IF JY = JY + INCY 20 CONTINUE ELSE IF( INCX.GT.0 )THEN KX = 1 ELSE KX = 1 - ( M - 1 )*INCX END IF DO 40, J = 1, N IF( Y( JY ).NE.ZERO )THEN TEMP = ALPHA*Y( JY ) IX = KX DO 30, I = 1, M A( I, J ) = A( I, J ) + X( IX )*TEMP IX = IX + INCX 30 CONTINUE END IF JY = JY + INCY 40 CONTINUE END IF * RETURN * * End of DGER . * END