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

    Module f12asfe_mod

!     F12ASF 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                           :: mv
!     .. Parameters ..
      Complex (Kind=nag_wp), Parameter, Public :: four = (4.0_nag_wp,          &
                                          0.0_nag_wp)
      Complex (Kind=nag_wp), Parameter, Public :: one = (1.0_nag_wp,0.0_nag_wp &
                                          )
      Complex (Kind=nag_wp), Parameter, Public :: six = (6.0_nag_wp,0.0_nag_wp &
                                          )
      Complex (Kind=nag_wp), Parameter, Public :: two = (2.0_nag_wp,0.0_nag_wp &
                                          )
      Integer, Parameter, Public       :: imon = 1, licomm = 140, nerr = 6,    &
                                          nin = 5, nout = 6
    Contains
      Subroutine mv(nx,v,w)
!       Compute the out-of--place matrix vector multiplication Y<---M*X,
!       where M is mass matrix formed by using piecewise linear elements
!       on [0,1].

!       .. Use Statements ..
        Use nag_library, Only: zscal
!       .. 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)          :: h
        Integer                        :: j, n
!       .. Intrinsic Procedures ..
        Intrinsic                      :: cmplx
!       .. Executable Statements ..
        n = nx*nx
        w(1) = (four*v(1)+v(2))/six
        Do j = 2, n - 1
          w(j) = (v(j-1)+four*v(j)+v(j+1))/six
        End Do
        w(n) = (v(n-1)+four*v(n))/six

        h = one/cmplx(n+1,kind=nag_wp)
!       The NAG name equivalent of zscal is f06gdf
        Call zscal(n,h,w,1)
        Return
      End Subroutine mv
    End Module f12asfe_mod

    Program f12asfe

!     F12ASF Example Main Program

!     .. Use Statements ..
      Use f12asfe_mod, Only: four, imon, licomm, mv, nerr, nin, nout, one,     &
                             six, two
      Use nag_library, Only: dznrm2, f12anf, f12apf, f12aqf, f12arf, f12asf,   &
                             nag_wp, zgttrf, zgttrs
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Complex (Kind=nag_wp)            :: h, rho, s, s1, s2, s3, sigma
      Real (Kind=nag_wp)               :: nev_nrm
      Integer                          :: ifail, ifail1, info, irevcm, j,      &
                                          lcomm, ldv, n, nconv, ncv, nev,      &
                                          niter, nshift, nx
!     .. Local Arrays ..
      Complex (Kind=nag_wp), Allocatable :: ax(:), comm(:), d(:,:), dd(:),     &
                                          dl(:), du(:), du2(:), mx(:),         &
                                          resid(:), v(:,:), x(:)
      Integer                          :: icomm(licomm)
      Integer, Allocatable             :: ipiv(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: cmplx
!     .. Executable Statements ..
      Write (nout,*) 'F12ASF Example Program Results'
      Write (nout,*)
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) nx, nev, ncv

      n = nx*nx
      lcomm = 3*n + 3*ncv*ncv + 5*ncv + 60
      ldv = n
      Allocate (ax(n),comm(lcomm),d(ncv,2),dd(n),dl(n),du(n),du2(n),mx(n),     &
        resid(n),v(ldv,ncv),x(n),ipiv(n))

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

!     Set the mode.
      ifail = 0
      Call f12arf('SHIFTED INVERSE',icomm,comm,ifail)
!     Set problem type.
      Call f12arf('GENERALIZED',icomm,comm,ifail)
      sigma = (5000.0_nag_wp,0.0_nag_wp)
      rho = (10.0_nag_wp,0.0_nag_wp)
      h = one/cmplx(n+1,kind=nag_wp)
      s = rho/two
      s1 = -one/h - s - sigma*h/six
      s2 = two/h - four*sigma*h/six
      s3 = -one/h + s - sigma*h/six

      dl(1:n-1) = s1
      dd(1:n-1) = s2
      du(1:n-1) = s3
      dd(n) = s2

!     The NAG name equivalent of zgttrf is f07crf
      Call zgttrf(n,dl,dd,du,du2,ipiv,info)
      If (info/=0) Then
        Write (nerr,99999) info
        Go To 100
      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) Then
!         Perform  x <--- OP*x = inv[A-SIGMA*M]*M*x
          Call mv(nx,x,ax)
          x(1:n) = ax(1:n)
!         The NAG name equivalent of zgttrs is f07csf
          Call zgttrs('N',n,1,dl,dd,du,du2,ipiv,x,n,info)
          If (info/=0) Then
            Write (nerr,99998) info
            Exit revcm
          End If
        Else If (irevcm==1) Then
!         Perform  x <--- OP*x = inv[A-SIGMA*M]*M*x,
!         MX stored in COMM from location IPNTR(3)
!         The NAG name equivalent of zgttrs is f07csf
          Call zgttrs('N',n,1,dl,dd,du,du2,ipiv,mx,n,info)
          x(1:n) = mx(1:n)
          If (info/=0) Then
            Write (nerr,99998) info
            Exit revcm
          End If
        Else If (irevcm==2) Then
!         Perform  y <--- M*x
          Call mv(nx,x,ax)
          x(1:n) = ax(1:n)
        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
          nev_nrm = dznrm2(nev,d(1,2),1)
          Write (6,99997) niter, nconv, nev_nrm
        End If
      End Do revcm

      If (ifail==0 .And. info==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,99996) nconv, sigma
        Write (nout,99995)(j,d(j,1),j=1,nconv)
      End If
100   Continue

99999 Format (1X,'** Error status returned by ZGTTRF, INFO =',I12)
99998 Format (1X,'** Error status returned by ZGTTRS, INFO =',I12)
99997 Format (1X,'Iteration',1X,I3,', No. converged =',1X,I3,', norm o',       &
        'f estimates =',E12.4)
99996 Format (1X,/,' The ',I4,' generalized Ritz values closest to (',F8.3,    &
        ',',F8.3,') are:',/)
99995 Format (1X,I8,5X,'( ',F10.4,' , ',F10.4,' )')
    End Program f12asfe