!   C02AGF Example Program Text
!   Mark 23 Release. NAG Copyright 2011.

    MODULE c02agfe_mod

!      C02AGF Example Program Module: Parameters

!      .. Implicit None Statement ..
       IMPLICIT NONE
!      .. Parameters ..
       INTEGER, PARAMETER              :: nin = 5, nout = 6
       LOGICAL, PARAMETER              :: scal = .TRUE.
    END MODULE c02agfe_mod
    PROGRAM c02agfe

!      C02AGF Example Main Program

!      .. Use Statements ..
       USE c02agfe_mod, ONLY : nout
!      .. Implicit None Statement ..
       IMPLICIT NONE
!      .. External Subroutines ..
       EXTERNAL                           ex1, ex2
!      .. Executable Statements ..
       WRITE (nout,*) 'C02AGF Example Program Results'

       CALL ex1
       CALL ex2

    END PROGRAM c02agfe
    SUBROUTINE ex1

!      .. Use Statements ..
       USE nag_library, ONLY : c02agf, nag_wp
       USE c02agfe_mod, ONLY : nin, nout, scal
!      .. Implicit None Statement ..
       IMPLICIT NONE
!      .. Local Scalars ..
       REAL (KIND=nag_wp)              :: zi, zr
       INTEGER                         :: i, ifail, n, nroot
!      .. Local Arrays ..
       REAL (KIND=nag_wp), ALLOCATABLE :: a(:), w(:), z(:,:)
!      .. Intrinsic Functions ..
       INTRINSIC                          abs
!      .. Executable Statements ..
       WRITE (nout,*)
       WRITE (nout,*)
       WRITE (nout,*) 'Example 1'

!      Skip heading in data file
       READ (nin,*)
       READ (nin,*)
       READ (nin,*)

       READ (nin,*) n
       ALLOCATE (a(0:n),w(2*(n+1)),z(2,n))

       READ (nin,*) (a(i),i=0,n)

       WRITE (nout,*)
       WRITE (nout,99999) 'Degree of polynomial = ', n

       ifail = 0
       CALL c02agf(a,n,scal,z,w,ifail)

       WRITE (nout,99998) 'Computed roots of polynomial'

       nroot = 1

       DO WHILE (nroot<=n)

          zr = z(1,nroot)
          zi = z(2,nroot)
          IF (zi==0.0E0_nag_wp) THEN
             WRITE (nout,99997) 'z = ', zr
             nroot = nroot + 1
          ELSE
             WRITE (nout,99997) 'z = ', zr, ' +/- ', abs(zi), '*i'
             nroot = nroot + 2
          END IF

       END DO

99999  FORMAT (/1X,A,I4)
99998  FORMAT (/1X,A/)
99997  FORMAT (1X,A,1P,E12.4,A,E12.4,A)
    END SUBROUTINE ex1
    SUBROUTINE ex2

!      .. Use Statements ..
       USE nag_library, ONLY : a02abf, c02agf, nag_wp, x02ajf, x02alf
       USE c02agfe_mod, ONLY : nin, nout, scal
!      .. Implicit None Statement ..
       IMPLICIT NONE
!      .. Local Scalars ..
       REAL (KIND=nag_wp)              :: deltac, deltai, di, eps, epsbar, f,  &
                                          r1, r2, r3, rmax
       INTEGER                         :: i, ifail, j, jmin, n
!      .. Local Arrays ..
       REAL (KIND=nag_wp), ALLOCATABLE :: a(:), abar(:), r(:), w(:), z(:,:),   &
                                          zbar(:,:)
       INTEGER, ALLOCATABLE            :: m(:)
!      .. Intrinsic Functions ..
       INTRINSIC                          abs, max, min
!      .. Executable Statements ..
       WRITE (nout,*)
       WRITE (nout,*)
       WRITE (nout,*) 'Example 2'

!      Skip heading in data file
       READ (nin,*)
       READ (nin,*)

       READ (nin,*) n
       ALLOCATE (a(0:n),abar(0:n),r(n),w(2*(n+1)),z(2,n),zbar(2,n),m(n))

!      Read in the coefficients of the original polynomial.

       READ (nin,*) (a(i),i=0,n)

!      Compute the roots of the original polynomial.

       ifail = 0
       CALL c02agf(a,n,scal,z,w,ifail)

!      Form the coefficients of the perturbed polynomial.

       eps = x02ajf()
       epsbar = 3.0_nag_wp*eps

       DO i = 0, n

          IF (a(i)/=0.0_nag_wp) THEN
             f = 1.0_nag_wp + epsbar
             epsbar = -epsbar
             abar(i) = f*a(i)
          ELSE
             abar(i) = 0.0E0_nag_wp
          END IF

       END DO

!      Compute the roots of the perturbed polynomial.

       ifail = 0
       CALL c02agf(abar,n,scal,zbar,w,ifail)

!      Perform error analysis.

!      Initialize markers to 0 (unmarked).

       m(1:n) = 0

       rmax = x02alf()

!      Loop over all unperturbed roots (stored in Z).

       DO i = 1, n
          deltai = rmax
          r1 = a02abf(z(1,i),z(2,i))

!         Loop over all perturbed roots (stored in ZBAR).

          DO j = 1, n

!            Compare the current unperturbed root to all unmarked
!            perturbed roots.

             IF (m(j)==0) THEN
                r2 = a02abf(zbar(1,j),zbar(2,j))
                deltac = abs(r1-r2)

                IF (deltac<deltai) THEN
                   deltai = deltac
                   jmin = j
                END IF

             END IF

          END DO

!         Mark the selected perturbed root.

          m(jmin) = 1

!         Compute the relative error.

          IF (r1/=0.0E0_nag_wp) THEN
             r3 = a02abf(zbar(1,jmin),zbar(2,jmin))
             di = min(r1,r3)
             r(i) = max(deltai/max(di,deltai/rmax),eps)
          ELSE
             r(i) = 0.0_nag_wp
          END IF

       END DO

       WRITE (nout,*)
       WRITE (nout,99999) 'Degree of polynomial = ', n
       WRITE (nout,*)
       WRITE (nout,*) 'Computed roots of polynomial   ', '  Error estimates'
       WRITE (nout,*) '                            ', '   (machine-dependent)'
       WRITE (nout,*)

       DO i = 1, n
          WRITE (nout,99998) 'z = ', z(1,i), z(2,i), '*i', r(i)
       END DO

99999  FORMAT (1X,A,I4)
99998  FORMAT (1X,A,1P,E12.4,SP,E12.4,A,5X,SS,E9.1)
    END SUBROUTINE ex2