End of file
Contents
Index
SUBROUTINE CYTNPD (N,DM,DU,CR,RS,X,MARK)
C
C*****************************************************************
C *
C Solving a system of linear equations *
C A * X = RS *
C for a cyclically tridiagonal, symmetric, strongly *
C nonsingular matrix A. *
C The matrix A is given by two *
C N-vectors DM and DU. The system of equations has the form: *
C *
C DM(1) * X(1) + DU(1) * X(2) + DU(N) * X(N) = RS(1) *
C *
C DU(I-1) * X(I-1) + DM(I) * X(I) + DU(I) * X(I+1) = RS(I) *
C for I = 2, ..., N - 1, and *
C *
C DU(N) * X(1) + DU(N-1) * X(N-1) + DM(N) * X(N) = RS(N) *
C *
C *
C *
C INPUT PARAMETERS: *
C ================= *
C N : number of equations, N > 2 *
C DM : N-vector DM(1:N); main diagonal of A *
C DM(1), DM(2), ... , DM(N) *
C DU : N-vector DU(1:N); upper co-diagonal of A *
C DU(1), DU(2), ... , DU(N-1); the off-diagonal *
C element A(1,N) is stored in DU(N). *
C RS : N-vector RS(1:N); the right hand side *
C *
C *
C OUTPUT PARAMETERS: *
C ================== *
C DM :) *
C DU :) overwritten with intermediate vectors *
C CR :) *
C RS :) *
C X : N-vector X(1:N), containing the solution *
C MARK : error parameter *
C MARK=-2 : condition N > 2 is not satisfied *
C MARK=-1 : A is strongly nonsingular, but not *
C positive definite *
C MARK= 0 : numerically the matrix A is not strongly *
C nonsingular *
C MARK= 1 : A is positive definite *
C *
C NOTE : If MARK = +/- 1, the determinant of A is given as: *
C DET A = DM(1) * DM(2) * ... * DM(N) *
C *
C----------------------------------------------------------------*
C *
C subroutines required: CYTSYP, CYTSYS, MACHPD *
C *
C*****************************************************************
C *
C authors : Gisela Engeln-Muellges *
C date : 01.07.1992 *
C source : FORTRAN 77 *
C *
C*****************************************************************
C
IMPLICIT DOUBLE PRECISION (A-H,O-Z)
DOUBLE PRECISION DM(1:N),DU(1:N),CR(1:N),RS(1:N),X(1:N)
MARK = -2
IF (N .LT. 3) RETURN
C
C factorization of the matrix A
C
CALL CNPSYP (N,DM,DU,CR,MARK)
C
C if MARK = +/- 1, update and backsubstitute
C
IF ((MARK .EQ. 1) .OR. (MARK .EQ. -1)) THEN
CALL CYTSYS (N,DM,DU,CR,RS,X)
ENDIF
RETURN
END
C
C
SUBROUTINE CNPSYP (N,DM,DU,CR,MARK)
C
C*****************************************************************
C *
C Factoring a cyclically tridiagonal, symmetric and strongly *
C nonsingular matrix A, that is given by the two N-vectors *
C DM and DU, into its Cholesky factors *
C A = R(TRANSP) * D * R *
C by applying the root-free Cholesky-method for tridiagonal *
C cyclic matrices. The form of the system of equations is *
C identical to the one described in SUBROUTINE CYTSY. *
C *
C *
C INPUT PARAMETERS: *
C ================= *
C N : number of equations, N > 2 *
C DM : N-vector DM(1:N); main diagonal of A *
C DM(1), DM(2), ... , DM(N) *
C DU : N-vector DU(1:N); upper co-diagonal of A *
C DU(1), DU(2), ... , DU(N-1); the off-diagonal *
C element A(1,N) is stored in DU(N). *
C Due to symmetry the lower co-diagonal does not need *
C to be stored separately. *
C *
C *
C OUTPUT PARAMETER: *
C ================= *
C DM :) overwritten with auxiliary vectors from the *
C DU :) factorization of A. The co-diagonal of the unit *
C CR :) upper tridiagonal matrix R is stored in DU, the *
C diagonal matrix D appears in DM and the right hand *
C side in CR. *
C MARK : error parameter *
C MARK=-2 : condition N > 2 is not satisfied *
C MARK=-1 : A is strongly nonsingular, but not *
C positive definite *
C MARK= 0 : numerically A is not strongly *
C nonsingular *
C MARK= 1 : A is positive definite *
C *
C NOTE : If MARK = +/- 1, then the inertia of A, i. e., the *
C number of positive and negative eigenvalues of A, *
C is the same as the number of positive and negative *
C numbers among the components of DM. *
C *
C----------------------------------------------------------------*
C *
C subroutine required: MACHPD *
C *
C*****************************************************************
C *
C authors : Gisela Engeln-Muellges *
C date : 01.07.1992 *
C source : FORTRAN 77 *
C *
C*****************************************************************
C
IMPLICIT DOUBLE PRECISION (A-H,O-Z)
DOUBLE PRECISION DM(1:N),DU(1:N),CR(1:N)
C
C calculating the machine constant
C
FMACHP = 1.0D0
10 FMACHP = 0.5D0 * FMACHP
IF (MACHPD(1.0D0+FMACHP) .EQ. 1) GOTO 10
FMACHP = FMACHP * 2.0D0
C
C determinaing the relative error bound
C
EPS = 4.0D0 * FMACHP
C
C testing of condition N > 2
C
MARK = -2
IF (N .LT. 3) RETURN
MARK = 1
C
C checking for strong nonsingularity of A for N=1
C
ROW = DABS(DM(1)) + DABS(DU(1)) + DABS(DU(N))
IF (ROW .EQ. 0.0D0) THEN
MARK = 0
RETURN
END IF
D = 1.0D0/ROW
IF (DM(1) .LT. 0.0D0) MARK = -1
IF (DABS(DM(1))*D .LE. EPS) THEN
MARK = 0
RETURN
END IF
C
C factoring A while checking for strong nonsingularity
C
DUMMY = DU(1)
DU(1) = DU(1)/DM(1)
CR(1) = DU(N)/DM(1)
DO 20 I=2,N-1,1
ROW = DABS(DM(I)) + DABS(DU(I)) + DABS(DUMMY)
IF (ROW .EQ. 0.0D0) THEN
MARK = 0
RETURN
END IF
D = 1.0D0/ROW
DM(I) = DM(I) - DUMMY * DU(I-1)
IF (DM(I) .LT. 0.0D0) MARK = -1
IF (DABS(DM(I))*D .LE. EPS) THEN
MARK = 0
RETURN
ENDIF
IF (I .LT. (N-1)) THEN
CR(I) = -DUMMY * CR(I-1)/DM(I)
DUMMY = DU(I)
DU(I) = DU(I)/DM(I)
ELSE
DUMMY2 = DU(I)
DU(I) = (DU(I) - DUMMY * CR(I-1))/DM(I)
ENDIF
20 CONTINUE
ROW = DABS(DU(N)) + DABS(DM(N)) + DABS(DUMMY2)
IF (ROW .EQ. 0.0D0) THEN
MARK = 0
RETURN
END IF
D = 1.0D0/ROW
DM(N) = DM(N) - DM(N-1) * DU(N-1) * DU(N-1)
DUMMY = 0.0D0
DO 30 I=1,N-2,1
DUMMY = DUMMY + DM(I) * CR(I) * CR(I)
30 CONTINUE
DM(N) = DM(N) - DUMMY
IF (DM(N) .LT. 0) MARK = -1
IF (DABS(DM(N))*D .LE. EPS) THEN
MARK = 0
RETURN
ENDIF
RETURN
END
Begin of file
Contents
Index