Example description
!   F12ACF Example Program Text
!   Mark 26.2 Release. NAG Copyright 2017.

    Module f12acfe_mod

!     F12ACF 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, mv
!     .. Parameters ..
      Real (Kind=nag_wp), Parameter, Public :: one = 1.0_nag_wp
      Integer, Parameter, Public       :: imon = 0, nin = 5, nout = 6
    Contains
      Subroutine av(nx,rho,v,w)

!       .. Parameters ..
        Real (Kind=nag_wp), Parameter  :: two = 2.0_nag_wp
!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In) :: rho
        Integer, Intent (In)           :: nx
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (In) :: v(nx*nx)
        Real (Kind=nag_wp), Intent (Out) :: w(nx*nx)
!       .. Local Scalars ..
        Real (Kind=nag_wp)             :: dd, dl, du, h, s
        Integer                        :: j, n
!       .. Intrinsic Procedures ..
        Intrinsic                      :: real
!       .. Executable Statements ..
        n = nx*nx
        h = one/real(n+1,kind=nag_wp)
        s = rho/two
        dd = two/h
        dl = -one/h - s
        du = -one/h + s
        w(1) = dd*v(1) + du*v(2)
        Do j = 2, n - 1
          w(j) = dl*v(j-1) + dd*v(j) + du*v(j+1)
        End Do
        w(n) = dl*v(n-1) + dd*v(n)
        Return
      End Subroutine av

      Subroutine mv(nx,v,w)

!       .. Use Statements ..
        Use nag_library, Only: dscal
!       .. Parameters ..
        Real (Kind=nag_wp), Parameter  :: four = 4.0_nag_wp
!       .. Scalar Arguments ..
        Integer, Intent (In)           :: nx
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (In) :: v(nx*nx)
        Real (Kind=nag_wp), Intent (Out) :: w(nx*nx)
!       .. Local Scalars ..
        Real (Kind=nag_wp)             :: h
        Integer                        :: j, n
!       .. Intrinsic Procedures ..
        Intrinsic                      :: real
!       .. Executable Statements ..
        n = nx*nx
        w(1) = four*v(1) + one*v(2)
        Do j = 2, n - 1
          w(j) = one*v(j-1) + four*v(j) + one*v(j+1)
        End Do
        w(n) = one*v(n-1) + four*v(n)
        h = one/real(n+1,kind=nag_wp)
!       The NAG name equivalent of dscal is f06edf
        Call dscal(n,h,w,1)
        Return
      End Subroutine mv
    End Module f12acfe_mod

    Program f12acfe

!     F12ACF Example Main Program

!     .. Use Statements ..
      Use f12acfe_mod, Only: av, imon, mv, nin, nout, one
      Use nag_library, Only: dnrm2, dpttrf, dpttrs, f12aaf, f12abf, f12acf,    &
                             f12adf, f12aef, nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Real (Kind=nag_wp)               :: h, rho, sigmai, sigmar
      Integer                          :: ifail, ifail1, info, irevcm, j,      &
                                          lcomm, ldv, licomm, n, nconv, ncv,   &
                                          nev, niter, nshift, nx
!     .. Local Arrays ..
      Real (Kind=nag_wp), Allocatable  :: comm(:), d(:,:), md(:), me(:),       &
                                          mx(:), resid(:), v(:,:), x(:)
      Integer, Allocatable             :: icomm(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: real
!     .. Executable Statements ..
      Write (nout,*) 'F12ACF Example Program Results'
      Write (nout,*)
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) nx, nev, ncv, rho
      n = nx*nx
      ldv = n
      licomm = 140
      lcomm = 3*n + 3*ncv*ncv + 6*ncv + 60
      Allocate (comm(lcomm),d(ncv,3),md(n),me(n-1),mx(n),resid(n),v(ldv,ncv),  &
        x(n),icomm(licomm))

!     ifail: behaviour on error exit
!             =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft
      ifail = 0
      Call f12aaf(n,nev,ncv,icomm,licomm,comm,lcomm,ifail)

!     Set the mode.
      ifail = 0
      Call f12adf('REGULAR INVERSE',icomm,comm,ifail)
!     Set problem type.
      Call f12adf('GENERALIZED',icomm,comm,ifail)
!     Use pointers to Workspace in calculating matrix vector
!     products rather than interfacing through the array X
      Call f12adf('POINTERS=YES',icomm,comm,ifail)

!     Construct M, and factorize using DPTTRF/F07JDF.
      h = one/real(n+1,kind=nag_wp)
      md(1:n-1) = 4.0_nag_wp*h
      me(1:n-1) = h
      md(n) = 4.0_nag_wp*h

!     The NAG name equivalent of dpttrf is f07jdf
      Call dpttrf(n,md,me,info)

      irevcm = 0

      ifail = -1
loop: Do
        Call f12abf(irevcm,resid,v,ldv,x,mx,nshift,comm,icomm,ifail)

        If (irevcm/=5) Then
          Select Case (irevcm)
          Case (-1,1)
!           Perform  y <--- OP*x = inv[M]*A*x using DPTTRS/F07JEF.
            Call av(nx,rho,comm(icomm(1)),comm(icomm(2)))
!           The NAG name equivalent of dpttrs is f07jef
            Call dpttrs(n,1,md,me,comm(icomm(2)),n,info)
          Case (2)
!           Perform  y <--- M*x.
            Call mv(nx,comm(icomm(1)),comm(icomm(2)))
          Case (4)
            If (imon/=0) Then
!             Output monitoring information if required.
              Call f12aef(niter,nconv,d,d(1,2),d(1,3),icomm,comm)
!             The NAG name equivalent of dnrm2 is f06ejf
              Write (6,99999) niter, nconv, dnrm2(nev,d(1,3),1)
            End If
          End Select
        Else
          Exit loop
        End If
      End Do loop
      If (ifail==0) Then
!       Post-Process using F12ACF to compute eigenvalues/vectors.
        ifail1 = 0
        Call f12acf(nconv,d,d(1,2),v,ldv,sigmar,sigmai,resid,v,ldv,comm,icomm, &
          ifail1)
!       Print computed eigenvalues.
        Write (nout,99998) nconv
        Do j = 1, nconv
          Write (nout,99997) j, d(j,1), d(j,2)
        End Do
      End If

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