Example description
!   F12ANF Example Program Text
!   Mark 27.0 Release. NAG Copyright 2019.

    Module f12anfe_mod

!     F12ANF Example Program Module:
!            Parameters and User-defined Routines

!     .. Use Statements ..
      Use nag_library, Only: nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Accessibility Statements ..
      Private
      Public                           :: av
!     .. Parameters ..
      Integer, Parameter, Public       :: imon = 0, ipoint = 0, licomm = 140,  &
                                          nin = 5, nout = 6
    Contains
      Subroutine tv(nx,x,y)
!       Compute the matrix vector multiplication y<---T*x where T is a nx
!       by nx tridiagonal matrix.

!       .. Parameters ..
        Complex (Kind=nag_wp), Parameter :: four = (4.0_nag_wp,0.0_nag_wp)
        Complex (Kind=nag_wp), Parameter :: half = (0.5_nag_wp,0.0_nag_wp)
        Complex (Kind=nag_wp), Parameter :: rho = (100.0_nag_wp,0.0_nag_wp)
!       .. Scalar Arguments ..
        Integer, Intent (In)           :: nx
!       .. Array Arguments ..
        Complex (Kind=nag_wp), Intent (In) :: x(nx)
        Complex (Kind=nag_wp), Intent (Out) :: y(nx)
!       .. Local Scalars ..
        Complex (Kind=nag_wp)          :: dd, dl, du, h, h2
        Integer                        :: j
!       .. Intrinsic Procedures ..
        Intrinsic                      :: cmplx
!       .. Executable Statements ..
        h = cmplx(nx+1,kind=nag_wp)
        h2 = h*h
        dd = four*h2
        dl = -h2 - half*rho*h
        du = -h2 + half*rho*h

        y(1) = dd*x(1) + du*x(2)
        Do j = 2, nx - 1
          y(j) = dl*x(j-1) + dd*x(j) + du*x(j+1)
        End Do
        y(nx) = dl*x(nx-1) + dd*x(nx)
        Return
      End Subroutine tv
      Subroutine av(nx,v,w)

!       .. Use Statements ..
        Use nag_library, Only: zaxpy
!       .. Scalar Arguments ..
        Integer, Intent (In)           :: nx
!       .. Array Arguments ..
        Complex (Kind=nag_wp), Intent (In) :: v(nx*nx)
        Complex (Kind=nag_wp), Intent (Out) :: w(nx*nx)
!       .. Local Scalars ..
        Complex (Kind=nag_wp)          :: h2
        Integer                        :: j, lo
!       .. Intrinsic Procedures ..
        Intrinsic                      :: cmplx
!       .. Executable Statements ..
        h2 = cmplx(-(nx+1)*(nx+1),kind=nag_wp)

        Call tv(nx,v(1),w(1))
!       The NAG name equivalent of zaxpy is f06gcf
        Call zaxpy(nx,h2,v(nx+1),1,w(1),1)

        Do j = 2, nx - 1
          lo = (j-1)*nx
          Call tv(nx,v(lo+1),w(lo+1))
          Call zaxpy(nx,h2,v(lo-nx+1),1,w(lo+1),1)
          Call zaxpy(nx,h2,v(lo+nx+1),1,w(lo+1),1)
        End Do

        lo = (nx-1)*nx
        Call tv(nx,v(lo+1),w(lo+1))
        Call zaxpy(nx,h2,v(lo-nx+1),1,w(lo+1),1)

        Return
      End Subroutine av
    End Module f12anfe_mod
    Program f12anfe

!     F12ANF Example Main Program

!     .. Use Statements ..
      Use f12anfe_mod, Only: av, imon, ipoint, licomm, nin, nout
      Use nag_library, Only: dznrm2, f12anf, f12apf, f12aqf, f12arf, f12asf,   &
                             nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Complex (Kind=nag_wp)            :: sigma
      Integer                          :: i, ifail, ifail1, irevcm, lcomm,     &
                                          ldv, n, nconv, ncv, nev, niter,      &
                                          nshift, nx
!     .. Local Arrays ..
      Complex (Kind=nag_wp), Allocatable :: ax(:), comm(:), d(:,:), mx(:),     &
                                          resid(:), v(:,:), x(:)
      Integer                          :: icomm(licomm)
!     .. Executable Statements ..
      Write (nout,*) 'F12ANF Example Program Results'
      Write (nout,*)
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) nx, nev, ncv

      n = nx*nx
      ldv = n
      lcomm = 3*n + 3*ncv*ncv + 5*ncv + 60
      Allocate (ax(n),comm(lcomm),d(ncv,2),mx(n),resid(n),v(n,ncv),x(n))

      ifail = 0
      Call f12anf(n,nev,ncv,icomm,licomm,comm,lcomm,ifail)

      If (ipoint/=0) Then
!       Use pointers to Workspace in calculating matrix vector
!       products rather than interfacing through the array X.
        ifail = 0
        Call f12arf('POINTERS=YES',icomm,comm,ifail)
      End If

      irevcm = 0
      ifail = -1

revcm: Do
        Call f12apf(irevcm,resid,v,ldv,x,mx,nshift,comm,icomm,ifail)
        If (irevcm==5) Then
          Exit revcm
        Else If (irevcm==-1 .Or. irevcm==1) Then
!         Perform matrix vector multiplication y <--- Op*x.
          If (ipoint==0) Then
            Call av(nx,x,ax)
            x(1:n) = ax(1:n)
          Else
            Call av(nx,comm(icomm(1)),comm(icomm(2)))
          End If
        Else If (irevcm==4 .And. imon/=0) Then
!         Output monitoring information.
          Call f12asf(niter,nconv,d,d(1,2),icomm,comm)
!         The NAG name equivalent of dznrm2 is f06jjf
          Write (6,99999) niter, nconv, dznrm2(nev,d(1,2),1)
        End If
      End Do revcm

      If (ifail==0) Then
!       Post-Process using F12AQF to compute eigenvalues/vectors.
        ifail1 = 0
        Call f12aqf(nconv,d,v,ldv,sigma,resid,v,ldv,comm,icomm,ifail1)

        Write (nout,99998) nconv
        Write (nout,99997)(i,d(i,1),i=1,nconv)
      End If

99999 Format (1X,'Iteration',1X,I3,', No. converged =',1X,I3,', norm o',       &
        'f estimates =',E16.8)
99998 Format (1X,/,' The ',I4,' Ritz values of largest magnitude are:',/)
99997 Format (1X,I8,5X,'( ',F12.4,' , ',F12.4,' )')
    End Program f12anfe