C********************************************************************** C C Copyright (C) 1991-1992 Roland W. Freund and Noel M. Nachtigal C All rights reserved. C C This code is part of a copyrighted package. For details, see the C file "cpyrit.doc" in the top-level directory. C C ***************************************************************** C ANY USE OF THIS CODE CONSTITUTES ACCEPTANCE OF THE TERMS OF THE C COPYRIGHT NOTICE C ***************************************************************** C C********************************************************************** C C This file contains the routines for the QMR algorithm for C unsymmetric matrices, using the three-term recurrence variant of C the look-ahead Lanczos algorithm. C C********************************************************************** C SUBROUTINE SUQMR (NDIM,NLEN,NLIM,MAXVW,M,NORM,DWK,IDX,IWK, $ VECS,TOL,INFO) C C Purpose: C This subroutine uses the QMR algorithm based on the three-term C variant of the look-ahead Lanczos process to solve linear C systems. It runs the QMR algorithm to convergence or until a C user-specified iteration limit is reached. C C The code is set up to solve the system A x = b with initial C guess x_0 = 0. Here A x = b denotes the preconditioned system, C and it is connected with the original system as follows. Let C B y = c be the original unpreconditioned system to be solved, and C let y_0 be an arbitrary initial guess for its solution. Then: C A x = b, where A = M_1^{-1} B M_2^{-1}, C x = M_2 (y - y_0), b = M_1^{-1} (c - B y_0). C Here M = M_1 M_2 is the preconditioner. C C To recover the final iterate y_n for the original system B y = c C from the final iterate x_n for the preconditioned system A x = b, C set C y_n = y_0 + M_2^{-1} x_n. C C The algorithm was first described in the RIACS Technical Report C 90.46, `An Implementation of the Look-Ahead Lanczos Algorithm for C Non-Hermitian Matrices, Part II`, by R.W. Freund and N.M. C Nachtigal, November 1990. C C Parameters: C For a description of the parameters, see the file "suqmr.doc" in C the current directory. C C External routines used: C subroutine saxpby(n,z,a,x,b,y) C Library routine, computes z = a * x + b * y. C single precision sdot(n,dx,incx,dy,incy) C BLAS-1 routine, computes y^H * x. C single precision slamch(ch) C LAPACK routine, computes machine-related constants. C single precision snrm2(n,dx,incx) C BLAS-1 routine, computes the 2-norm of x. C subroutine srandn(n,x,seed) C Library routine, fills x with random numbers. C subroutine srotg(da,db,cos,sin) C BLAS-1 routine, computes the Givens rotation which rotates the C vector [a; b] into [ sqrt(a**2 + b**2); 0 ]. C subroutine suqmr1(ndim,nlen,m,maxvw,numchk,n,l,lstar,nl,nlstar,vf, C ierr,adjust,norm1,norma,tmax,tnrm,dwk,iwk,idx,vecs) C Low-level routine, builds the vectors v_{n+1} and w_{n+1}. C single precision suqmro(n) C User-supplied routine, specifies the QMR scaling factors. C C Noel M. Nachtigal C October 24, 1990 C C********************************************************************** C INTRINSIC ABS, FLOAT, AMAX1, SQRT, MAX0, MIN0, MOD EXTERNAL SAXPBY, SDOT, SLAMCH, SNRM2, SRANDN, SROTG, SUQMR1 EXTERNAL SUQMRO REAL SDOT, SLAMCH, SNRM2, SUQMRO C INTEGER INFO(4), M, MAXVW, NDIM, NLEN, NLIM INTEGER IDX(3,NLIM+2), IWK(M,4) REAL DWK(M,5*M+13), NORM, TOL REAL VECS(NDIM,4*MAXVW+4+4) C C Common block variables. C C C Common block SUQMRX. C REAL NORMA COMMON /SUQMRX/NORMA C C Miscellaneous parameters. C REAL DHUN, DONE, DTEN, SZERO PARAMETER (DHUN = 1.0E2,DONE = 1.0E0,DTEN = 1.0E1,SZERO = 0.0E0) C C Local variables, permanent. C INTEGER IERR, IEND, L, LSTAR, MVWBLT, N, NL, NLSTAR, NMAX, NQMR SAVE IERR, IEND, L, LSTAR, MVWBLT, N, NL, NLSTAR, NMAX, NQMR INTEGER NUMCHK, RETLBL, TF, TRES, VF SAVE NUMCHK, RETLBL, TF, TRES, VF REAL ADJUST, MAXOMG, NORM1, R0, RESN, TMAX, TMIN SAVE ADJUST, MAXOMG, NORM1, R0, RESN, TMAX, TMIN REAL TNRM, UCHK, UNRM SAVE TNRM, UCHK, UNRM C C Local variables, transient. C INTEGER I, IBASE, INIT, ISJ, ISN, J, REVCOM REAL STMP, STMP1, STMP2 C C Initialize some of the permanent variables. C DATA RETLBL /0/ C C Check the reverse communication flag to see where to branch. C REVCOM RETLBL Comment C 0 0 first call, go to label 10 C 1 40 returning from AXB, go to label 40 C 1 100 returning from AXB, go to label 100 C 2 50 returning from ATXB, go to label 50 C REVCOM = INFO(2) INFO(2) = 0 IF (REVCOM.EQ.0) THEN N = 0 NQMR = 0 MVWBLT = 0 IF (RETLBL.EQ.0) GO TO 10 ELSE IF (REVCOM.EQ.1) THEN IF (RETLBL.EQ.40) THEN GO TO 40 ELSE IF (RETLBL.EQ.100) THEN GO TO 100 END IF ELSE IF (REVCOM.EQ.2) THEN IF (RETLBL.EQ.50) GO TO 50 END IF IERR = 1 GO TO 130 C C Check whether the inputs are valid. C 10 IERR = 0 IF (NDIM.LT.1) IERR = 2 IF (NLEN.LT.1) IERR = 2 IF (NLIM.LT.1) IERR = 2 IF (MAXVW.LT.1) IERR = 2 IF (NLEN.GT.NDIM) IERR = 2 IF (M.LT.2*MAXVW+2) IERR = 2 IF (IERR.NE.0) GO TO 130 C C Extract from INFO the output units TF and VF, the true residual C flag TRES, and the left starting vector flag INIT. C VF = MAX0(INFO(1),0) INIT = VF / 100000 VF = VF - INIT * 100000 TRES = VF / 10000 VF = VF - TRES * 10000 TF = VF / 100 VF = VF - TF * 100 C C Extract the norms. C NORMA = SZERO NORM1 = ABS(NORM) C C Set the adjustment parameters. C NUMCHK = 25 ADJUST = DTEN C C Extract and check the various tolerances and norm estimates. C TNRM = SLAMCH('E') * DTEN TMIN = SQRT(SQRT(SLAMCH('S'))) TMAX = DONE / TMIN IF (TOL.LE.SZERO) TOL = SQRT(SLAMCH('E')) C C Start the trace messages and convergence history. C IF (VF.NE.0) WRITE (VF,'(I8,2E11.4)') 0, DONE, DONE IF (TF.NE.0) WRITE (TF,'(I8,2E11.4)') 0, DONE, DONE C C Set up wrapped indices. The following indices are used: C IDX(1,I) = indices used for the work arrays row-dimensioned M; C IDX(2,I) = indices used for v_i; C IDX(3,I) = indices used for w_i. C DO 20 I = 1, NLIM+2 IDX(1,I) = MOD(I-1,M) + 1 IDX(2,I) = 4 + 1 + 0*(MAXVW+2) + MOD(I-1,MAXVW+2) IDX(3,I) = 4 + 1 + 1*(MAXVW+2) + MOD(I-1,MAXVW+2) 20 CONTINUE C C Set x_0 = 0 and compute the norm of the initial residual. C CALL SAXPBY (NLEN,VECS(1,IDX(2,1)),DONE,VECS(1,2),SZERO,VECS(1, $IDX(2,1))) CALL SAXPBY (NLEN,VECS(1,1),SZERO,VECS(1,1),SZERO,VECS(1,1)) R0 = SNRM2(NLEN,VECS(1,IDX(2,1)),1) IF ((TOL.GE.DONE).OR.(R0.EQ.SZERO)) GO TO 130 C C Check whether the auxiliary vector must be supplied. C IF (INIT.EQ.0) CALL SRANDN (NLEN,VECS(1,3),1) CALL SAXPBY (NLEN,VECS(1,IDX(3,1)),DONE,VECS(1,3),SZERO,VECS(1, $IDX(3,1))) C C Scale the first pair of Lanczos vectors and check for invariant C subspaces. C STMP1 = R0 STMP2 = SNRM2(NLEN,VECS(1,IDX(3,1)),1) IF (STMP1.LT.TNRM) IERR = IERR + 16 IF (STMP2.LT.TNRM) IERR = IERR + 32 IF (IERR.NE.0) GO TO 130 DWK(IDX(1,1),5*M+8) = DONE DWK(IDX(1,1),IDX(1,1)) = SDOT(NLEN,VECS(1,IDX(2,1)),1,VECS(1,IDX(3 $,1)),1) / ( STMP1 * STMP2 ) IF ((STMP1.GE.TMAX).OR.(STMP1.LE.TMIN)) THEN STMP = DONE / STMP1 CALL SAXPBY (NLEN,VECS(1,IDX(2,1)),STMP,VECS(1,IDX(2,1)),SZERO, $VECS(1,IDX(2,1))) STMP1 = DONE END IF IF ((STMP2.GE.TMAX).OR.(STMP2.LE.TMIN)) THEN STMP = DONE / STMP2 CALL SAXPBY (NLEN,VECS(1,IDX(3,1)),STMP,VECS(1,IDX(3,1)),SZERO, $VECS(1,IDX(3,1))) STMP2 = DONE END IF DWK(IDX(1,1),5*M+9) = DONE / STMP1 DWK(IDX(1,1),5*M+10) = DONE / STMP2 C C Initialize the counters. C L = 1 N = 1 NMAX = 0 IWK(IDX(1,0+1),1) = 1 IWK(IDX(1,1+1),1) = 1 MVWBLT = 1 NL = IWK(IDX(1,L+1),1) C C Set up the QMR iteration. C RESN = DONE DWK(IDX(1,1),5*M+12) = SUQMRO(1) DWK(IDX(1,1),5*M+4) = DWK(IDX(1,1),5*M+12) * R0 MAXOMG = DONE / DWK(IDX(1,1),5*M+12) C C This is one step of the three-term Lanczos algorithm. C V_n, W_n, D_{n-1}, F_{n-1}, H_{n-1}, and D_{nn} are given. C C Note that the previous step was not necessarily step N-1, as it C could have been a restart. C Except at the first step, the following hold: C NL.LE.N, N-NL.LE.MAXVW-1. C 30 IWK(IDX(1,N),4) = L IWK(IDX(1,N),3) = N C C Build v_{n+1} and w_{n+1}. C Have the caller carry out AXB, then return here. C CALL AXB (VECS(1,IDX(2,N)),VECS(1,3)) C INFO(2) = 1 INFO(3) = IDX(2,N) INFO(4) = 3 RETLBL = 40 RETURN C C Have the caller carry out ATXB, then return here. C CALL ATXB (VECS(1,IDX(3,N)),VECS(1,4)) C 40 INFO(2) = 2 INFO(3) = IDX(3,N) INFO(4) = 4 RETLBL = 50 RETURN 50 IF ((VF.NE.0).AND.(N.LE.NMAX)) WRITE (VF,'(A11,I8)') $ 'Rebuilding:', N+1 CALL SUQMR1 (NDIM,NLEN,M,MAXVW,NUMCHK,N,L,LSTAR,NL,NLSTAR,VF,IERR, $ ADJUST,NORM1,NORMA,TMAX,TMIN,TNRM,DWK,IWK,IDX,VECS) IF ((IERR.NE.0).AND.(IERR.NE.16).AND.(IERR.NE.32) $ .AND.(IERR.NE.48)) GO TO 130 MVWBLT = MAX0(MVWBLT,NL-IWK(IDX(1,L-1+1),1)) IWK(IDX(1,N),2) = NLSTAR C C Update the counter for steps taken. C NMAX = MAX0(NMAX,N) IF (N.LT.NL-1) GO TO 120 C C The QMR code starts here. C At this point, steps up to NL-1 are guaranteed not to be rebuilt. C No errors are allowed in IERR, other than possibly having found C one or both invariant subspaces, in which case all remaining C iterates are computed. C IEND = NL-1 IF (IERR.NE.0) IEND = N 60 IF (NQMR.GT.IEND-1) GO TO 120 C C Update the QMR iteration counter. C NQMR = NQMR + 1 C C Get the next scaling factor omega(i) and update MAXOMG. C DWK(IDX(1,NQMR+1),5*M+12) = SUQMRO(NQMR+1) MAXOMG = AMAX1(MAXOMG,DONE/DWK(IDX(1,NQMR+1),5*M+12)) C C Compute the starting index IBASE for the column of \hat{R}. C IBASE = MAX0(1,IWK(IDX(1,NQMR),2)-1) DWK(IDX(1,IBASE),5*M+5) = SZERO C C Multiply the new column by the previous omegas. C DO 70 J = IWK(IDX(1,NQMR),2), NQMR+1 DWK(IDX(1,J),5*M+5) = DWK(IDX(1,J),5*M+12) * DWK(IDX(1,J),2*M+ $IDX(1,NQMR)) 70 CONTINUE C C Apply the previous rotations. C The loop below explicitly implements a call to SROT: C CALL SROT (1,DWK(IDX(1,J-1),5*M+5),1,DWK(IDX(1,J),5*M+5),1,DWK(IDX(1,J),5*M+13),DWK(IDX(1,J),5*M+6)) C DO 80 J = IBASE+1, NQMR STMP1 = DWK(IDX(1,J),5*M+5) STMP2 = DWK(IDX(1,J-1),5*M+5) DWK(IDX(1,J-1),5*M+5) = DWK(IDX(1,J),5*M+13) * STMP2 + $DWK(IDX(1,J),5*M+6) * STMP1 DWK(IDX(1,J),5*M+5) = DWK(IDX(1,J),5*M+13) * STMP1 - $DWK(IDX(1,J),5*M+6) * STMP2 80 CONTINUE C C Compute the rotation for the last element (this also applies it). C CALL SROTG (DWK(IDX(1,NQMR),5*M+5),DWK(IDX(1,NQMR+1),5*M+5), $DWK(IDX(1,NQMR+1),5*M+13),DWK(IDX(1,NQMR+1),5*M+6)) C C Apply the new rotation to the right-hand side vector. C Could be replaced with: C DWK(IDX(1,NQMR+1),5*M+4) = SZERO C CALL SROT (1,DWK(IDX(1,NQMR),5*M+4),1,DWK(IDX(1,NQMR+1),5*M+4),1,DWK(IDX(1,NQMR+1),5*M+13),DWK(IDX(1,NQMR+1),5*M+6)) C DWK(IDX(1,NQMR+1),5*M+4) = -DWK(IDX(1,NQMR+1),5*M+6) * DWK(IDX(1, $NQMR),5*M+4) DWK(IDX(1,NQMR),5*M+4) = DWK(IDX(1,NQMR+1),5*M+13) * DWK(IDX(1, $NQMR),5*M+4) C C Compute the next search direction s_i. C This is more complicated than it might have to be because storage C for the vectors VECS(1,ISN) is minimized. C STMP2 = SZERO STMP = DWK(IDX(1,NQMR),5*M+9) ISN = 4 + 1 + 2*(MAXVW+2) + MOD(NQMR-1,2*MAXVW) DO 90 J = IBASE, NQMR-1 STMP1 = DWK(IDX(1,J),5*M+5) * DWK(IDX(1,J),5*M+7) / STMP IF (STMP1.EQ.SZERO) GO TO 90 ISJ = 4 + 1 + 2*(MAXVW+2) + MOD(J-1,2*MAXVW) CALL SAXPBY (NLEN,VECS(1,ISN),STMP2,VECS(1,ISN),-STMP1,VECS(1, $ISJ)) STMP2 = DONE 90 CONTINUE CALL SAXPBY (NLEN,VECS(1,ISN),STMP2,VECS(1,ISN),DONE,VECS(1,IDX(2, $NQMR))) STMP = STMP / DWK(IDX(1,NQMR),5*M+5) C C Compute the new QMR iterate, then scale the search direction. C STMP1 = STMP * DWK(IDX(1,NQMR),5*M+4) CALL SAXPBY (NLEN,VECS(1,1),DONE,VECS(1,1),STMP1,VECS(1,ISN)) DWK(IDX(1,NQMR),5*M+7) = STMP STMP = ABS(STMP) IF ((STMP.GE.TMAX).OR.(STMP.LE.TMIN)) THEN DWK(IDX(1,NQMR),5*M+7) = DONE CALL SAXPBY (NLEN,VECS(1,ISN),STMP,VECS(1,ISN),SZERO,VECS(1,ISN $)) END IF C C Compute the residual norm upper bound. C If the scaled upper bound is within one order of magnitude of the C target convergence norm, compute the true residual norm. C UNRM = SQRT(FLOAT(NQMR+1)) * MAXOMG * ABS(DWK(IDX(1,NQMR+1),5*M+4) $) / R0 UCHK = UNRM IF ((TRES.EQ.0).AND.(UNRM/TOL.GT.DTEN).AND.(N.LT.NLIM)) GO TO 110 C C Have the caller carry out AXB, then return here. C CALL AXB (VECS(1,1),VECS(1,3)) C INFO(2) = 1 INFO(3) = 1 INFO(4) = 3 RETLBL = 100 RETURN 100 CALL SAXPBY (NLEN,VECS(1,3),DONE,VECS(1,2),-DONE,VECS(1,3)) RESN = SNRM2(NLEN,VECS(1,3),1) / R0 UCHK = RESN C C Output the convergence history. C 110 IF (VF.NE.0) WRITE (VF,'(I8,2E11.4)') NQMR, UNRM, RESN IF (TF.NE.0) WRITE (TF,'(I8,2E11.4)') NQMR, UNRM, RESN C C Check for convergence or termination. Stop if: C 1. algorithm converged; C 2. there is an error condition; C 3. the residual norm upper bound is smaller than the computed C residual norm by a factor of at least 100; C 4. algorithm exceeded the iterations limit. C IF (RESN.LE.TOL) THEN IERR = 0 GO TO 130 ELSE IF (IERR.NE.0) THEN GO TO 130 ELSE IF (UNRM.LT.UCHK/DHUN) THEN IERR = 4 GO TO 130 ELSE IF (NQMR.GE.NLIM) THEN IERR = 4 GO TO 130 END IF GO TO 60 C C Update the running counter. C 120 IF (IERR.NE.0) GO TO 130 IERR = 4 IF (N.GE.NLIM) GO TO 130 N = N + 1 GO TO 30 C C Done. C 130 RETLBL = 0 NLIM = NQMR INFO(1) = IERR NORM = NORM1 MAXVW = MVWBLT C RETURN END C C********************************************************************** C REAL FUNCTION SUQMRH(I,J) C C Purpose: C Returns the recurrence coefficients for inner vectors. C C Parameters: C I = row index of the coefficient (input). C J = column index of the coefficient (input). C C Noel M. Nachtigal C July 9, 1993 C C********************************************************************** C C C Common block SUQMRX. C REAL NORMA COMMON /SUQMRX/NORMA C C INTEGER I, J C IF ((I.LT.1).OR.(J.LT.1)) THEN SUQMRH = 0.0E0 ELSE IF (I.EQ.J) THEN SUQMRH = NORMA ELSE SUQMRH = 0.0E0 END IF C RETURN END C C********************************************************************** C REAL FUNCTION SUQMRO (I) C C Purpose: C Returns the scaling parameter OMEGA(I). C C Parameters: C I = the index of the parameter OMEGA (input). C C Noel M. Nachtigal C October 7, 1990 C C********************************************************************** C INTEGER I C SUQMRO = 1.0E0 C RETURN END C C********************************************************************** C SUBROUTINE SUQMR1 (NDIM,NLEN,M,MAXVW,NUMCHK,N,L,LSTAR,NL,NLSTAR, $ VF,IERR,ADJUST,NORM1,NORMA,TMAX,TMIN,TNRM,DWK, $ IWK,IDX,VECS) C C Purpose: C This subroutine builds a new pair of vectors v_{n+1} and w_{n+1}. C It is called internally by the Lanczos code, and is not meant to C to be called directly by the user code. C C Parameters: C See descriptions in "sulal.doc" or "suqmr.doc". C C External routines used: C subroutine saxpby(n,z,a,x,b,y) C Library routine, computes z = a * x + b * y. C single precision sdot(n,x,incx,y,incy) C BLAS-1 routine, computes y^H * x. C single precision slamch(ch) C LAPACK routine, computes machine-related constants. C single precision snrm2(n,x,incx) C BLAS-1 routine, computes the 2-norm of x. C subroutine sqrdc(x,ldx,n,p,qraux,jpvt,work,job) C LINPACK routine, computes the QR factorization of x. C subroutine sqrsl(x,ldx,n,k,qraux,y,qy,qty,b,rsd,xb,job,info) C LINPACK routine, applies the QR factorization of x. C subroutine ssvdc(x,ldx,n,p,s,e,u,ldu,v,ldv,work,job,info) C LINPACK routine, computes the SVD of x. C subroutine suqmr2 (ndim,nlen,m,dwk,iwk,idx,vecs) C Forces closure of an inner block. C single precision suqmrh(i,n) C User-supplied routine, computes inner recurrence coefficients. C C Noel M. Nachtigal C July 9, 1993 C C********************************************************************** C INTRINSIC ABS, AMAX1, MAX0 EXTERNAL SAXPBY, SDOT, SLAMCH, SNRM2, SQRDC, SQRSL, SSVDC EXTERNAL SUQMR2, SUQMRH REAL SDOT, SLAMCH, SNRM2, SUQMRH C INTEGER IERR, L, LSTAR, M, MAXVW, N, NDIM, NL, NLEN, NLSTAR INTEGER IDX(3,3), IWK(M,4), NUMCHK, VF REAL ADJUST, NORM1, NORMA, TMAX, TMIN, TNRM REAL DWK(M,5*M+13), VECS(NDIM,4*MAXVW+4+4) C C Miscellaneous parameters. C REAL DONE, SZERO PARAMETER (DONE = 1.0E0,SZERO = 0.0E0) C C Local variables. C INTEGER I, IJ, J, LBLKSZ, NP1 REAL STMP, STMP1, STMP2, STMP3, STMP4, STMP5 REAL ISTMP1, ISTMP2, ISTMP3, ISTMP4 LOGICAL IBUILT, INNER, RERUN C C Initialize local variables. C NP1 = N + 1 DWK(IDX(1,N),5*M+11) = SZERO RERUN = .FALSE. IBUILT = .FALSE. IERR = 0 LBLKSZ = N - NL + 1 C C Update the norm estimates. C IF (N.LE.NUMCHK) THEN STMP1 = SNRM2(NLEN,VECS(1,3),1) * DWK(IDX(1,N),5*M+9) STMP2 = SNRM2(NLEN,VECS(1,4),1) * DWK(IDX(1,N),5*M+10) NORMA = AMAX1(STMP1,STMP2,NORMA) END IF NORM1 = AMAX1(ADJUST*NORMA,NORM1) C C Clear current column of H. C DO 100 I = 1, M DWK(IDX(1,I),2*M+IDX(1,N)) = SZERO 100 CONTINUE C C Update D^{(n-1)} to D^{(n)}. C DO 120 I = NL, N-1 STMP = SZERO DO 110 J = NL, N-1 STMP = STMP + DWK(IDX(1,I),IDX(1,J)) * DWK(IDX(1,J),2*M+ $IDX(1,N-1)) 110 CONTINUE STMP = ( DWK(IDX(1,I),M+IDX(1,N-1)) - STMP ) / DWK(IDX(1,N) $,2*M+IDX(1,N-1)) DWK(IDX(1,I),IDX(1,N)) = STMP DWK(IDX(1,N),IDX(1,I)) = STMP * DWK(IDX(1,N),5*M+8) / DWK(IDX(1 $,I),5*M+8) 120 CONTINUE C C Compute l^\star. C LSTAR = MAX0(1,L-1) NLSTAR = IWK(IDX(1,LSTAR+1),1) C C Compute F_{n,1:n-1} and F_{1:n-1,n}. Delay computing F_{nn}. C DO 140 I = MAX0(1,NL-1), N-1 STMP = SZERO DO 130 J = NL, I+1 STMP = STMP + DWK(IDX(1,N),IDX(1,J)) * DWK(IDX(1,J),2*M+ $IDX(1,I)) 130 CONTINUE DWK(IDX(1,N),M+IDX(1,I)) = STMP DWK(IDX(1,I),M+IDX(1,N)) = STMP * DWK(IDX(1,I),5*M+8) / $DWK(IDX(1,N),5*M+8) NORMA = AMAX1(NORMA,ABS(DWK(IDX(1,N),M+IDX(1,I))), $ABS(DWK(IDX(1,I),M+IDX(1,N)))) NORM1 = AMAX1(ADJUST*NORMA,NORM1) 140 CONTINUE C C Check whether D_l is nonsingular. C IF (LBLKSZ.EQ.1) THEN INNER = ABS(DWK(IDX(1,NL),IDX(1,NL))).EQ.SZERO ELSE DO 160 I = NL, N DO 150 J = NL, N DWK(I-NL+1,4*M+J-NL+1) = DWK(IDX(1,I),IDX(1,J)) 150 CONTINUE 160 CONTINUE CALL SSVDC (DWK(1,4*M+1),M,LBLKSZ,LBLKSZ,DWK(1,5*M+1),DWK(1,5*M $+3),SZERO,0,SZERO,0,DWK(1,5*M+2),0,I) IF (I.NE.0) THEN IERR = -I GOTO 370 END IF STMP = SZERO IF (DWK(1,5*M+1).NE.SZERO) STMP = DWK(LBLKSZ,5*M+1) / DWK(1,5*M $+1) INNER = STMP.LT.(FLOAT(NLEN) * SLAMCH('E')) END IF IF ((VF.NE.0).AND.INNER) $ WRITE (VF,'(A31)') '... moment matrix D is singular' C C Compute H_{n_{l-1}:n_l-1,n}, for n_l > 1. C Also, build the part common to both inner and regular vectors. C Assume that VECS(1,3) = A v_n and VECS(1,4) = A^T w_n. C DWK(IDX(1,NP1),2*M+IDX(1,N)) = SZERO IF (NL.EQ.1) GO TO 180 IF (NL.EQ.NLSTAR+1) THEN DWK(IDX(1,NLSTAR),2*M+IDX(1,N)) = DWK(IDX(1,NLSTAR),M+IDX(1,N)) $ / DWK(IDX(1,NLSTAR),IDX(1,NLSTAR)) STMP = DWK(IDX(1,NLSTAR),2*M+IDX(1,N)) * DWK(IDX(1,NLSTAR),5*M+ $9) / DWK(IDX(1,N),5*M+9) CALL SAXPBY (NLEN,VECS(1,3),DONE,VECS(1,3),-STMP,VECS(1,IDX(2, $NLSTAR))) STMP = DWK(IDX(1,NLSTAR),2*M+IDX(1,N)) * DWK(IDX(1,NLSTAR),5*M+ $10) / DWK(IDX(1,N),5*M+10) * DWK(IDX(1,N),5*M+8) / DWK(IDX(1, $NLSTAR),5*M+8) CALL SAXPBY (NLEN,VECS(1,4),DONE,VECS(1,4),-STMP,VECS(1,IDX(3, $NLSTAR))) ELSE IF (NL.GT.NLSTAR+1) THEN C C If a block of size larger than 1 was just closed, combine all its C vectors into one. C STMP2 = SZERO DO 170 I = NL-1, NLSTAR, -1 DWK(IDX(1,I),2*M+IDX(1,N)) = DWK(IDX(1,NL-1),M+IDX(1,N)) * $DWK(IDX(1,I),3*M+IDX(1,L-1)) IF (N.NE.NL) GO TO 170 STMP1 = DWK(IDX(1,I),3*M+IDX(1,L-1)) * DWK(IDX(1,I),5*M+9) CALL SAXPBY (NLEN,VECS(1,IDX(2,NL-1)),STMP2,VECS(1,IDX(2,NL- $1)),STMP1,VECS(1,IDX(2,I))) STMP1 = DWK(IDX(1,I),3*M+IDX(1,L-1)) * DWK(IDX(1,I),5*M+10) $/ DWK(IDX(1,I),5*M+8) CALL SAXPBY (NLEN,VECS(1,IDX(3,NL-1)),STMP2,VECS(1,IDX(3,NL- $1)),STMP1,VECS(1,IDX(3,I))) STMP2 = DONE 170 CONTINUE STMP1 = DWK(IDX(1,NL-1),M+IDX(1,N)) / DWK(IDX(1,N),5*M+9) CALL SAXPBY (NLEN,VECS(1,3),DONE,VECS(1,3),-STMP1,VECS(1,IDX(2, $NL-1))) STMP1 = DWK(IDX(1,NL-1),M+IDX(1,N)) * DWK(IDX(1,N),5*M+8) / $DWK(IDX(1,N),5*M+10) CALL SAXPBY (NLEN,VECS(1,4),DONE,VECS(1,4),-STMP1,VECS(1,IDX(3, $NL-1))) END IF C C Compute F_{nn}. C 180 DWK(IDX(1,N),M+IDX(1,N)) = SDOT(NLEN,VECS(1,IDX(3,N)),1,VECS(1,3), $1) * DWK(IDX(1,N),5*M+9) * DWK(IDX(1,N),5*M+10) NORMA = AMAX1(NORMA,ABS(DWK(IDX(1,N),M+IDX(1,N)))) NORM1 = AMAX1(ADJUST*NORMA,NORM1) IF (INNER) GO TO 320 C C Compute H_{n_l:n,n}. C IWK(IDX(1,L+1+1),1) = NP1 IF (LBLKSZ.EQ.1) THEN DWK(IDX(1,NL),2*M+IDX(1,N)) = DWK(IDX(1,NL),M+IDX(1,N)) / $DWK(IDX(1,NL),IDX(1,NL)) ELSE DO 200 J = NL, N DWK(J-NL+1,5*M+1) = DWK(IDX(1,J),M+IDX(1,N)) DO 190 IJ = NL, N DWK(J-NL+1,4*M+IJ-NL+1) = DWK(IDX(1,J),IDX(1,IJ)) 190 CONTINUE 200 CONTINUE CALL SQRDC (DWK(1,4*M+1),M,LBLKSZ,LBLKSZ,DWK(1,5*M+3),0,SZERO,0 $) CALL SQRSL (DWK(1,4*M+1),M,LBLKSZ,LBLKSZ,DWK(1,5*M+3),DWK(1,5*M $+1), $ SZERO,DWK(1,5*M+1),DWK(1,5*M+2),DWK(1,5*M+1),SZERO, $100,J) DO 210 J = NL, N DWK(IDX(1,J),2*M+IDX(1,N)) = DWK(J-NL+1,5*M+2) 210 CONTINUE END IF C C Either D_l is nonsingular, or VW is being rerun. For the latter, C just finish building the vectors. For the former, the check for C H_{n_{l-1}:n,n} could be done. However, the look-ahead strategy C requires the smaller of the norm checks for H_{n_{l-1}:n,n} and C H_{n_l:n,n+1}, so both norms are computed first and then checked. C C Build regular vectors. C STMP = DWK(IDX(1,N),2*M+IDX(1,N)) CALL SAXPBY (NLEN,VECS(1,IDX(2,NP1)),DONE,VECS(1,3),-STMP,VECS(1, $IDX(2,N))) CALL SAXPBY (NLEN,VECS(1,IDX(3,NP1)),DONE,VECS(1,4),-STMP,VECS(1, $IDX(3,N))) DO 220 I = NL, N-1 STMP = DWK(IDX(1,I),2*M+IDX(1,N)) * DWK(IDX(1,I),5*M+9) / $DWK(IDX(1,N),5*M+9) CALL SAXPBY (NLEN,VECS(1,IDX(2,NP1)),DONE,VECS(1,IDX(2,NP1)),- $STMP,VECS(1,IDX(2,I))) STMP = DWK(IDX(1,I),2*M+IDX(1,N)) * DWK(IDX(1,I),5*M+10) / $DWK(IDX(1,N),5*M+10) * DWK(IDX(1,N),5*M+8) / DWK(IDX(1,I),5*M+8) CALL SAXPBY (NLEN,VECS(1,IDX(3,NP1)),DONE,VECS(1,IDX(3,NP1)),- $STMP,VECS(1,IDX(3,I))) 220 CONTINUE C C Compute scale factors for the new vectors. C 230 STMP3 = SNRM2(NLEN,VECS(1,IDX(2,NP1)),1) STMP4 = SNRM2(NLEN,VECS(1,IDX(3,NP1)),1) STMP1 = DWK(IDX(1,N),5*M+9) * STMP3 STMP2 = DWK(IDX(1,N),5*M+10) * STMP4 DWK(IDX(1,NP1),2*M+IDX(1,N)) = STMP1 IF (STMP1.LT.TNRM) IERR = IERR + 16 IF (STMP2.LT.TNRM) IERR = IERR + 32 IF (IERR.NE.0) GOTO 360 DWK(IDX(1,NP1),5*M+8) = DWK(IDX(1,N),5*M+8) * STMP1 / STMP2 DWK(IDX(1,NP1),IDX(1,NP1)) = SDOT(NLEN,VECS(1,IDX(3,NP1)),1,VECS(1 $,IDX(2,NP1)),1) / ( STMP3 * STMP4 ) IF ((STMP3.GE.TMAX).OR.(STMP3.LE.TMIN)) THEN STMP = DONE / STMP3 CALL SAXPBY (NLEN,VECS(1,IDX(2,NP1)),STMP,VECS(1,IDX(2,NP1)), $SZERO,VECS(1,IDX(2,NP1))) STMP3 = DONE END IF IF ((STMP4.GE.TMAX).OR.(STMP4.LE.TMIN)) THEN STMP = DONE / STMP4 CALL SAXPBY (NLEN,VECS(1,IDX(3,NP1)),STMP,VECS(1,IDX(3,NP1)), $SZERO,VECS(1,IDX(3,NP1))) STMP4 = DONE END IF DWK(IDX(1,NP1),5*M+9) = DONE / STMP3 DWK(IDX(1,NP1),5*M+10) = DONE / STMP4 C C Compute the last column of D_l^{-1}. C IF (LBLKSZ.EQ.1) THEN DWK(IDX(1,NL),3*M+IDX(1,L)) = DONE / DWK(IDX(1,NL),IDX(1,NL)) ELSE DO 240 I = 1, LBLKSZ-1 DWK(I,5*M+1) = SZERO 240 CONTINUE DWK(LBLKSZ,5*M+1) = DONE CALL SQRSL (DWK(1,4*M+1),M,LBLKSZ,LBLKSZ,DWK(1,5*M+3),DWK(1,5*M $+1), $ SZERO,DWK(1,5*M+1),DWK(1,5*M+2),DWK(1,5*M+1),SZERO, $100,J) DO 250 I = NL, N DWK(IDX(1,I),3*M+IDX(1,L)) = DWK(I-NL+1,5*M+2) 250 CONTINUE END IF C C If VW is being rerun, then skip to the end. C IF (RERUN) GO TO 360 C C Compute the norm of H_{n_{l-1}:n,n}. C STMP1 = SZERO STMP2 = SZERO DO 260 I = IWK(IDX(1,L-1+1),1), N STMP = ABS(DWK(IDX(1,I),2*M+IDX(1,N))) STMP1 = STMP1 + STMP STMP2 = STMP2 + STMP / DWK(IDX(1,I),5*M+8) 260 CONTINUE STMP2 = STMP2 * DWK(IDX(1,N),5*M+8) C C Build the 2nd term for the next step, H_{n_l:n,n+1}, regular. C Compute the norm of H_{n_l:n,n+1}. C STMP3 = SZERO STMP4 = SZERO STMP5 = DWK(IDX(1,NP1),IDX(1,NP1)) * DWK(IDX(1,NP1),2*M+IDX(1,N)) $* DWK(IDX(1,N),5*M+8) / DWK(IDX(1,NP1),5*M+8) DO 270 I = NL, N STMP = ABS(STMP5 * DWK(IDX(1,I),3*M+IDX(1,L))) STMP3 = STMP3 + STMP STMP4 = STMP4 + STMP / DWK(IDX(1,I),5*M+8) 270 CONTINUE STMP4 = STMP4 * DWK(IDX(1,NP1),5*M+8) C C Check H_{n_{l-1}:n,n} and H_{n_l:n-1,n+1}. C STMP = AMAX1(STMP1,STMP2,STMP3,STMP4) INNER = STMP.GT.NORM1 IF (.NOT.INNER) GO TO 360 DWK(IDX(1,N),5*M+11) = STMP C C If H_{n_{l-1}:n,n} is bad, build inner vectors. C IF (AMAX1(STMP1,STMP2).GT.NORM1) GO TO 320 C C If H_{n_l:n-1,n+1} is bad, check the inner vectors. C This only applies if n_l > 1. C IF (NL.LE.1) GO TO 320 C C Build the inner vectors to compute the 2nd term at the next step. C Get the coefficients H_{n_l:n+1,n} in a temporary location. C Build inner vectors in VECS(1,3) and VECS(1,4). C IBUILT = .TRUE. DO 280 I = NL, N DWK(IDX(1,I),2*M+IDX(1,NP1)) = SUQMRH(I,N) STMP = DWK(IDX(1,I),2*M+IDX(1,NP1)) * DWK(IDX(1,I),5*M+9) $ / DWK(IDX(1,N),5*M+9) CALL SAXPBY (NLEN,VECS(1,3),DONE,VECS(1,3),-STMP,VECS(1,IDX(2,I $))) STMP = DWK(IDX(1,I),2*M+IDX(1,NP1)) * DWK(IDX(1,I),5*M+10 $) / DWK(IDX(1,N),5*M+10) * DWK(IDX(1,N),5*M+8) / DWK(IDX(1,I),5*M+ $8) CALL SAXPBY (NLEN,VECS(1,4),DONE,VECS(1,4),-STMP,VECS(1,IDX(3,I $))) 280 CONTINUE ISTMP3 = SNRM2(NLEN,VECS(1,3),1) ISTMP4 = SNRM2(NLEN,VECS(1,4),1) ISTMP1 = DWK(IDX(1,N),5*M+9) * ISTMP3 ISTMP2 = DWK(IDX(1,N),5*M+10) * ISTMP4 DWK(IDX(1,NP1),2*M+IDX(1,NP1)) = ISTMP1 IF (ISTMP1.LT.TNRM) IERR = IERR + 16 IF (ISTMP2.LT.TNRM) IERR = IERR + 32 IF (IERR.NE.0) GOTO 330 C C Build the 2nd term for the next step, H_{n_{l-1}:n_l-1,n+1}, C inner. Compute the norm of H_{n_{l-1}:n_l-1,n+1}. C STMP = SZERO DO 290 J = NLSTAR, NL-1 STMP = STMP + DWK(IDX(1,NL),IDX(1,J)) * DWK(IDX(1,I),2*M+IDX(1, $N)) 290 CONTINUE DO 300 J = NL, N STMP = STMP + DWK(IDX(1,NL),IDX(1,J)) * DWK(IDX(1,I),2*M+IDX(1, $NP1)) 300 CONTINUE STMP1 = SZERO STMP2 = SZERO STMP3 = AMAX1(STMP3,STMP4) STMP = ( DWK(IDX(1,NL),M+IDX(1,N)) - STMP ) / DWK(IDX(1,NP1),2*M+ $IDX(1,NP1)) STMP5 = STMP * DWK(IDX(1,NL),2*M+IDX(1,NL-1)) * DWK(IDX(1,NL-1), $5*M+8) / DWK(IDX(1,NL),5*M+8) DO 310 I = IWK(IDX(1,L-1+1),1), NL-1 STMP = ABS(DWK(IDX(1,I),3*M+IDX(1,L-1)) * STMP5) STMP1 = STMP1 + STMP STMP2 = STMP2 + STMP / DWK(IDX(1,I),5*M+8) 310 CONTINUE STMP2 = STMP2 * DWK(IDX(1,N),5*M+8) C C Compare the inner and regular versions of the 2nd term at the next C step. Build the vector corresponding to the smaller term. C INNER = STMP3.GT.AMAX1(STMP1,STMP2) IF (.NOT.INNER) GO TO 360 C C Build inner vectors. C Check whether the block has to be forced to close. C 320 IF (VF.NE.0) WRITE (VF,'(A7,I5,A9)') 'Vector ',NP1,' is inner' IF (NP1-NL.EQ.MAXVW) THEN CALL SUQMR2 (NDIM,NLEN,M,N,L,LSTAR,NL,NLSTAR,VF,IERR,ADJUST, $ NORM1,DWK,IWK,IDX,VECS) IF (IERR.NE.0) GO TO 370 LBLKSZ = N - NL + 1 INNER = .FALSE. RERUN = .TRUE. NP1 = N + 1 GO TO 230 END IF C C The temporary vectors contain either just partial inner vectors, C or the completed ones, depending on whether IBUILT is TRUE or C FALSE. In either case, replace the regular vectors. C 330 CALL SAXPBY (NLEN,VECS(1,IDX(2,NP1)),DONE,VECS(1,3),SZERO,VECS(1, $IDX(2,NP1))) CALL SAXPBY (NLEN,VECS(1,IDX(3,NP1)),DONE,VECS(1,4),SZERO,VECS(1, $IDX(3,NP1))) C C Get the coefficients H_{n_l:n+1,n} and build inner vectors. C IF (IBUILT) THEN DO 340 I = NL, NP1 DWK(IDX(1,I),2*M+IDX(1,N)) = DWK(IDX(1,I),2*M+IDX(1,NP1)) 340 CONTINUE STMP1 = ISTMP1 STMP2 = ISTMP2 STMP3 = ISTMP3 STMP4 = ISTMP4 ELSE DO 350 I = NL, N DWK(IDX(1,I),2*M+IDX(1,N)) = SUQMRH(I,N) STMP = DWK(IDX(1,I),2*M+IDX(1,N)) * DWK(IDX(1,I),5*M+9) $/ DWK(IDX(1,N),5*M+9) CALL SAXPBY (NLEN,VECS(1,IDX(2,NP1)),DONE,VECS(1,IDX(2,NP1)) $,-STMP,VECS(1,IDX(2,I))) STMP = DWK(IDX(1,I),2*M+IDX(1,N)) * DWK(IDX(1,I),5*M+10) $ / DWK(IDX(1,N),5*M+10) * DWK(IDX(1,N),5*M+8) / DWK(IDX(1,I),5*M+8 $) CALL SAXPBY (NLEN,VECS(1,IDX(3,NP1)),DONE,VECS(1,IDX(3,NP1)) $,-STMP,VECS(1,IDX(3,I))) 350 CONTINUE STMP3 = SNRM2(NLEN,VECS(1,IDX(2,NP1)),1) STMP4 = SNRM2(NLEN,VECS(1,IDX(3,NP1)),1) STMP1 = DWK(IDX(1,N),5*M+9) * STMP3 STMP2 = DWK(IDX(1,N),5*M+10) * STMP4 END IF DWK(IDX(1,NP1),2*M+IDX(1,N)) = STMP1 IF (STMP1.LT.TNRM) IERR = IERR + 16 IF (STMP2.LT.TNRM) IERR = IERR + 32 IF (IERR.NE.0) GOTO 360 DWK(IDX(1,NP1),5*M+8) = DWK(IDX(1,N),5*M+8) * STMP1 / STMP2 DWK(IDX(1,NP1),IDX(1,NP1)) = SDOT(NLEN,VECS(1,IDX(3,NP1)),1,VECS(1 $,IDX(2,NP1)),1) / ( STMP3 * STMP4 ) IF ((STMP3.GE.TMAX).OR.(STMP3.LE.TMIN)) THEN STMP = DONE / STMP3 CALL SAXPBY (NLEN,VECS(1,IDX(2,NP1)),STMP,VECS(1,IDX(2,NP1)), $SZERO,VECS(1,IDX(2,NP1))) STMP3 = DONE END IF IF ((STMP4.GE.TMAX).OR.(STMP4.LE.TMIN)) THEN STMP = DONE / STMP4 CALL SAXPBY (NLEN,VECS(1,IDX(3,NP1)),STMP,VECS(1,IDX(3,NP1)), $SZERO,VECS(1,IDX(3,NP1))) STMP4 = DONE END IF DWK(IDX(1,NP1),5*M+9) = DONE / STMP3 DWK(IDX(1,NP1),5*M+10) = DONE / STMP4 C C If regular vectors were built, update the counters. C 360 IF (.NOT.INNER) THEN L = L + 1 NL = IWK(IDX(1,L+1),1) END IF C 370 RETURN END C C********************************************************************** C SUBROUTINE SUQMR2 (NDIM,NLEN,M,N,L,LSTAR,NL,NLSTAR,VF,IERR, $ ADJUST,NORM1,DWK,IWK,IDX,VECS) C C Purpose: C This subroutine rebuilds the data for vectors v_{n+1} and w_{n+1} C at the point where a block was forced to close. The routine is C called internally by the look-ahead Lanczos code and is not meant C to be called directly by the user code. C C Parameters: C See descriptions in "sulal.doc" or "suqmr.doc". C C External routines used: C subroutine saxpby(n,z,a,x,b,y) C Library routine, computes z = a * x + b * y. C subroutine sqrdc(x,ldx,n,p,qraux,jpvt,work,job) C LINPACK routine, computes the QR factorization of x. C subroutine sqrsl(x,ldx,n,k,qraux,y,qy,qty,b,rsd,xb,job,info) C LINPACK routine, applies the QR factorization of x. C C Noel M. Nachtigal C May 31, 1993 C C********************************************************************** C EXTERNAL SAXPBY, SQRDC, SQRSL C INTEGER IERR, L, LSTAR, M, N, NDIM, NL, NLEN, NLSTAR INTEGER IDX(3,3), IWK(M,4), VF REAL ADJUST, DWK(M,5*M+13), NORM1, VECS(NDIM,*) C C Miscellaneous parameters. C REAL DONE, SZERO PARAMETER (DONE = 1.0E0,SZERO = 0.0E0) C C Local variables. C INTEGER I, IJ, J, LBLKSZ, NP1 REAL STMP1, STMP2, SCALV, SCALW C C Find the index of the vector pair with the smallest pass value. C IERR = 0 IF (VF.NE.0) WRITE (VF,'(A20)') 'block did not close:' J = NL STMP1 = DWK(IDX(1,J),5*M+11) DO 100 I = NL+1, N STMP2 = DWK(IDX(1,I),5*M+11) IF (STMP2.GT.SZERO) THEN IF ((STMP1.EQ.SZERO).OR.(STMP2.LT.STMP1)) THEN J = I STMP1 = STMP2 END IF END IF 100 CONTINUE IF (STMP1.EQ.SZERO) THEN IF (VF.NE.0) WRITE (VF,'(A47)') $'... no new norm estimates available (aborting).' IERR = 8 RETURN END IF NORM1 = ADJUST * STMP1 IF (VF.NE.0) WRITE (VF,'(A40,I5,E11.4)') $ '... updated norms, restarting from step:', IWK(IDX(1,J),3), $ NORM1 IF (IWK(IDX(1,J),3).EQ.N) RETURN N = IWK(IDX(1,J),3) L = IWK(IDX(1,N),4) NL = IWK(IDX(1,L+1),1) C C Initialize local variables. C NP1 = N + 1 DWK(IDX(1,N),5*M+11) = SZERO LBLKSZ = N - NL + 1 C C Compute l^\star. C LSTAR = L - 1 NLSTAR = IWK(IDX(1,LSTAR+1),1) C C Compute H_{n_l:n,n}. Save the old coefficients. C IWK(IDX(1,L+1+1),1) = NP1 IF (LBLKSZ.EQ.1) THEN DWK(IDX(1,NL),2*M+IDX(1,NP1)) = DWK(IDX(1,NL),2*M+IDX(1,N)) DWK(IDX(1,NL),2*M+IDX(1,N)) = DWK(IDX(1,NL),M+IDX(1,N)) / $DWK(IDX(1,NL),IDX(1,NL)) ELSE DO 120 J = NL, N DWK(J-NL+1,5*M+1) = DWK(IDX(1,J),M+IDX(1,N)) DO 110 IJ = NL, N DWK(J-NL+1,4*M+IJ-NL+1) = DWK(IDX(1,J),IDX(1,IJ)) 110 CONTINUE 120 CONTINUE CALL SQRDC (DWK(1,4*M+1),M,LBLKSZ,LBLKSZ,DWK(1,5*M+3),0,SZERO,0 $) CALL SQRSL (DWK(1,4*M+1),M,LBLKSZ,LBLKSZ,DWK(1,5*M+3),DWK(1,5*M $+1), $ SZERO,DWK(1,5*M+1),DWK(1,5*M+2),DWK(1,5*M+1),SZERO, $100,J) DO 130 J = NL, N DWK(IDX(1,J),2*M+IDX(1,NP1)) = DWK(IDX(1,J),2*M+IDX(1,N)) DWK(IDX(1,J),2*M+IDX(1,N)) = DWK(J-NL+1,5*M+2) 130 CONTINUE END IF C C Convert inner vectors to regular vectors. C SCALV = DWK(IDX(1,NP1),2*M+IDX(1,N)) * DWK(IDX(1,NP1),5*M+9) SCALW = DWK(IDX(1,NP1),2*M+IDX(1,N)) * DWK(IDX(1,NP1),5*M+10) * $DWK(IDX(1,N),5*M+8) / DWK(IDX(1,NP1),5*M+8) DWK(IDX(1,NP1),2*M+IDX(1,N)) = SZERO DO 140 I = NL, N STMP1 = DWK(IDX(1,I),2*M+IDX(1,N)) - DWK(IDX(1,I),2*M+IDX(1,NP1 $)) STMP2 = STMP1 * DWK(IDX(1,I),5*M+9) / SCALV CALL SAXPBY (NLEN,VECS(1,IDX(2,NP1)),DONE,VECS(1,IDX(2,NP1)),- $STMP2,VECS(1,IDX(2,I))) STMP2 = STMP1 * DWK(IDX(1,I),5*M+10) / SCALW * DWK(IDX(1,N),5*M $+8) / DWK(IDX(1,I),5*M+8) CALL SAXPBY (NLEN,VECS(1,IDX(3,NP1)),DONE,VECS(1,IDX(3,NP1)),- $STMP2,VECS(1,IDX(3,I))) 140 CONTINUE C RETURN END C C**********************************************************************