PROGRAM f08safe ! F08SAF Example Program Text ! Mark 23 Release. NAG Copyright 2011. ! .. Use Statements .. USE nag_library, ONLY : ddisna, dsygv, dtrcon, f06rcf, nag_wp, x02ajf, & x04caf ! .. Implicit None Statement .. IMPLICIT NONE ! .. Parameters .. INTEGER, PARAMETER :: nb = 64, nin = 5, nout = 6 ! .. Local Scalars .. REAL (KIND=nag_wp) :: anorm, bnorm, eps, rcond, rcondb, & t1, t2, t3 INTEGER :: i, ifail, info, lda, ldb, lwork, n ! .. Local Arrays .. REAL (KIND=nag_wp), ALLOCATABLE :: a(:,:), b(:,:), eerbnd(:), & rcondz(:), w(:), work(:), zerbnd(:) REAL (KIND=nag_wp) :: dummy(1) INTEGER, ALLOCATABLE :: iwork(:) ! .. Intrinsic Functions .. INTRINSIC abs, max, nint ! .. Executable Statements .. WRITE (nout,*) 'F08SAF Example Program Results' WRITE (nout,*) ! Skip heading in data file READ (nin,*) READ (nin,*) n lda = n ldb = n ALLOCATE (a(lda,n),b(ldb,n),eerbnd(n),rcondz(n),w(n),zerbnd(n), & iwork(n)) ! Use routine workspace query to get optimal workspace. lwork = -1 ! The NAG name equivalent of dsygv is f08saf CALL dsygv(1,'Vectors','Upper',n,a,lda,b,ldb,w,dummy,lwork,info) ! Make sure that there is enough workspace for blocksize nb. lwork = max((nb+2)*n,nint(dummy(1))) ALLOCATE (work(lwork)) ! Read the upper triangular parts of the matrices A and B READ (nin,*) (a(i,i:n),i=1,n) READ (nin,*) (b(i,i:n),i=1,n) ! Compute the one-norms of the symmetric matrices A and B anorm = f06rcf('One norm','Upper',n,a,lda,work) bnorm = f06rcf('One norm','Upper',n,b,ldb,work) ! Solve the generalized symmetric eigenvalue problem ! A*x = lambda*B*x (ITYPE = 1) ! The NAG name equivalent of dsygv is f08saf CALL dsygv(1,'Vectors','Upper',n,a,lda,b,ldb,w,work,lwork,info) IF (info==0) THEN ! Print solution WRITE (nout,*) 'Eigenvalues' WRITE (nout,99999) w(1:n) WRITE (nout,*) FLUSH (nout) ! Normalize the eigenvectors DO i = 1, n a(1:n,i) = a(1:n,i)/a(1,i) END DO ! ifail: behaviour on error exit ! =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft ifail = 0 CALL x04caf('General',' ',n,n,a,lda,'Eigenvectors',ifail) ! Call DTRCON (F07TGF) to estimate the reciprocal condition ! number of the Cholesky factor of B. Note that: ! cond(B) = 1/rcond**2 ! The NAG name equivalent of dtrcon is f07tgf CALL dtrcon('One norm','Upper','Non-unit',n,b,ldb,rcond,work,iwork, & info) ! Print the reciprocal condition number of B rcondb = rcond**2 WRITE (nout,*) WRITE (nout,*) 'Estimate of reciprocal condition number for B' WRITE (nout,99998) rcondb FLUSH (nout) ! Get the machine precision, eps, and if rcondb is not less ! than eps**2, compute error estimates for the eigenvalues and ! eigenvectors eps = x02ajf() IF (rcond>=eps) THEN ! Call DDISNA (F08FLF) to estimate reciprocal condition ! numbers for the eigenvectors of (A - lambda*B) CALL ddisna('Eigenvectors',n,n,w,rcondz,info) ! Compute the error estimates for the eigenvalues and ! eigenvectors t1 = eps/rcondb t2 = anorm/bnorm t3 = t2/rcond DO i = 1, n eerbnd(i) = t1*(t2+abs(w(i))) zerbnd(i) = t1*(t3+abs(w(i)))/rcondz(i) END DO ! Print the approximate error bounds for the eigenvalues ! and vectors WRITE (nout,*) WRITE (nout,*) 'Error estimates for the eigenvalues' WRITE (nout,99998) eerbnd(1:n) WRITE (nout,*) WRITE (nout,*) 'Error estimates for the eigenvectors' WRITE (nout,99998) zerbnd(1:n) ELSE WRITE (nout,*) WRITE (nout,*) 'B is very ill-conditioned, error ', & 'estimates have not been computed' END IF ELSE IF (info>n .AND. info<=2*n) THEN i = info - n WRITE (nout,99997) 'The leading minor of order ', i, & ' of B is not positive definite' ELSE WRITE (nout,99996) 'Failure in DSYGV. INFO =', info END IF 99999 FORMAT (3X,(6F11.4)) 99998 FORMAT (4X,1P,6E11.1) 99997 FORMAT (1X,A,I4,A) 99996 FORMAT (1X,A,I4) END PROGRAM f08safe