Example description
    Program f08tcfe

!     F08TCF Example Program Text

!     Mark 27.1 Release. NAG Copyright 2020.

!     .. Use Statements ..
      Use nag_library, Only: dspgvd, dtpcon, f06rdf, nag_wp, x02ajf
!     .. Implicit None Statement ..
      Implicit None
!     .. Parameters ..
      Integer, Parameter               :: nin = 5, nout = 6
      Character (1), Parameter         :: uplo = 'U'
!     .. Local Scalars ..
      Real (Kind=nag_wp)               :: anorm, bnorm, eps, rcond, rcondb, t1
      Integer                          :: aplen, i, info, j, liwork, lwork, n
!     .. Local Arrays ..
      Real (Kind=nag_wp), Allocatable  :: ap(:), bp(:), eerbnd(:), w(:),       &
                                          work(:)
      Real (Kind=nag_wp)               :: dummy(1,1)
      Integer                          :: idum(1)
      Integer, Allocatable             :: iwork(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: abs, max, nint
!     .. Executable Statements ..
      Write (nout,*) 'F08TCF Example Program Results'
      Write (nout,*)
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) n
      aplen = (n*(n+1))/2
      Allocate (ap(aplen),bp(aplen),eerbnd(n),w(n))

!     Use routine workspace query to get optimal workspace.
      lwork = -1
      liwork = -1
!     The NAG name equivalent of dspgvd is f08tcf
      Call dspgvd(2,'No vectors',uplo,n,ap,bp,w,dummy,n,dummy,lwork,idum,      &
        liwork,info)

!     Make sure that there is at least minimum workspace.
      lwork = max(3*n,nint(dummy(1,1)))
      liwork = max(n,idum(1))
      Allocate (work(lwork),iwork(liwork))

!     Read the upper or lower triangular parts of the matrices A and
!     B from data file

      If (uplo=='U') Then
        Read (nin,*)((ap(i+(j*(j-1))/2),j=i,n),i=1,n)
        Read (nin,*)((bp(i+(j*(j-1))/2),j=i,n),i=1,n)
      Else If (uplo=='L') Then
        Read (nin,*)((ap(i+((2*n-j)*(j-1))/2),j=1,i),i=1,n)
        Read (nin,*)((bp(i+((2*n-j)*(j-1))/2),j=1,i),i=1,n)
      End If

!     Compute the one-norms of the symmetric matrices A and B

      anorm = f06rdf('One norm',uplo,n,ap,work)
      bnorm = f06rdf('One norm',uplo,n,bp,work)

!     Solve the generalized symmetric eigenvalue problem
!     A*B*x = lambda*x (itype = 2)

!     In the following call the 9th argument is set to n rather
!     than 1 to avoid an incorrect error message in some vendor
!     versions of LAPACK.
!     The NAG name equivalent of dspgvd is f08tcf
      Call dspgvd(2,'No vectors',uplo,n,ap,bp,w,dummy,n,work,lwork,iwork,      &
        liwork,info)

      If (info==0) Then

!       Print solution

        Write (nout,*) 'Eigenvalues'
        Write (nout,99999) w(1:n)

!       Call DTPCON (F07UGF) to estimate the reciprocal condition
!       number of the Cholesky factor of B.  Note that:
!       cond(B) = 1/rcond**2.  DTPCON requires WORK and IWORK to be
!       of length at least 3*n and n respectively

        Call dtpcon('One norm',uplo,'Non-unit',n,bp,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

!       Get the machine precision, eps, and if rcondb is not less
!       than eps**2, compute error estimates for the eigenvalues

        eps = x02ajf()
        If (rcond>=eps) Then
          t1 = anorm*bnorm
          Do i = 1, n
            eerbnd(i) = eps*(t1+abs(w(i))/rcondb)
          End Do

!         Print the approximate error bounds for the eigenvalues

          Write (nout,*)
          Write (nout,*) 'Error estimates for the eigenvalues'
          Write (nout,99998) eerbnd(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 DSPGVD. 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 f08tcfe