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

    Module f12fdfe_mod

!     F12FDF 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 ..
      Real (Kind=nag_wp), Parameter, Public :: four = 4.0_nag_wp
      Real (Kind=nag_wp), Parameter, Public :: one = 1.0_nag_wp
      Real (Kind=nag_wp), Parameter, Public :: six = 6.0_nag_wp
      Real (Kind=nag_wp), Parameter, Public :: two = 2.0_nag_wp
      Real (Kind=nag_wp), Parameter, Public :: zero = 0.0_nag_wp
      Integer, Parameter, Public       :: imon = 0, ipoint = 1, licomm = 140,  &
                                          nin = 5, nout = 6
    Contains
      Subroutine mv(n,v,w)

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

!     F12FDF Example Main Program

!     .. Use Statements ..
      Use f12fdfe_mod, Only: four, imon, ipoint, licomm, mv, nin, nout, one,   &
                             six, two, zero
      Use nag_library, Only: dcopy, dgttrf, dgttrs, dnrm2, f12faf, f12fbf,     &
                             f12fcf, f12fdf, f12fef, nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Real (Kind=nag_wp)               :: h, r1, r2, sigma
      Integer                          :: ifail, info, irevcm, j, lcomm, ldv,  &
                                          n, nconv, ncv, nev, niter, nshift
!     .. Local Arrays ..
      Real (Kind=nag_wp), Allocatable  :: ad(:), adl(:), adu(:), adu2(:),      &
                                          comm(:), d(:,:), mx(:), resid(:),    &
                                          v(:,:), x(:)
      Integer                          :: icomm(licomm)
      Integer, Allocatable             :: ipiv(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: real
!     .. Executable Statements ..
      Write (nout,*) 'F12FDF Example Program Results'
      Write (nout,*)
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) n, nev, ncv

      lcomm = 3*n + ncv*ncv + 8*ncv + 60
      ldv = n
      Allocate (ad(n),adl(n),adu(n),adu2(n),comm(lcomm),d(ncv,2),mx(n),        &
        resid(n),v(ldv,ncv),x(n),ipiv(n))

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

!     We are solving a generalized problem
      ifail = 0
      Call f12fdf('GENERALIZED',icomm,comm,ifail)
!     Indicate that we are using the shift and invert mode.
      Call f12fdf('SHIFTED INVERSE',icomm,comm,ifail)
      If (ipoint==1) Then
!       Use pointers to Workspace in calculating matrix vector
!       products rather than interfacing through the array X
        Call f12fdf('POINTERS=YES',icomm,comm,ifail)
      End If

      h = one/real(n+1,kind=nag_wp)
      r1 = (four/six)*h
      r2 = (one/six)*h
      sigma = zero
      ad(1:n) = two/h - sigma*r1
      adl(1:n) = -one/h - sigma*r2
      adu(1:n) = adl(1:n)
!     The NAG name equivalent of dgttrf is f07cdf
      Call dgttrf(n,adl,ad,adu,adu2,ipiv,info)

      irevcm = 0
      ifail = -1
revcm: Do
        Call f12fbf(irevcm,resid,v,ldv,x,mx,nshift,comm,icomm,ifail)
        If (irevcm==5) Then
          Exit revcm
        Else If (irevcm==-1) Then
!         Perform  y <--- OP*x = inv[A-SIGMA*M]*M*x
!         The NAG name equivalent of dgttrs is f07cef
          If (ipoint==0) Then
            Call mv(n,x,mx)
            x(1:n) = mx(1:n)
            Call dgttrs('N',n,1,adl,ad,adu,adu2,ipiv,x,n,info)
          Else
            Call mv(n,comm(icomm(1)),comm(icomm(2)))
            Call dgttrs('N',n,1,adl,ad,adu,adu2,ipiv,comm(icomm(2)),n,info)
          End If
        Else If (irevcm==1) Then
!         Perform y <-- OP*x = inv[A-sigma*M]*M*x
!         The NAG name equivalent of dgttrs is f07cef.
!         M*x has been saved in COMM(ICOMM(3)) or MX.
          If (ipoint==0) Then
            x(1:n) = mx(1:n)
            Call dgttrs('N',n,1,adl,ad,adu,adu2,ipiv,x,n,info)
          Else
!           The NAG name equivalent of dcopy is f06eff
            Call dcopy(n,comm(icomm(3)),1,comm(icomm(2)),1)
            Call dgttrs('N',n,1,adl,ad,adu,adu2,ipiv,comm(icomm(2)),n,info)
          End If
        Else If (irevcm==2) Then
!         Perform  y <--- M*x.
          If (ipoint==0) Then
            Call mv(n,x,mx)
          Else
            Call mv(n,comm(icomm(1)),comm(icomm(2)))
          End If
        Else If (irevcm==4 .And. imon/=0) Then
!         Output monitoring information
          Call f12fef(niter,nconv,d,d(1,2),icomm,comm)
!         The NAG name equivalent of dnrm2 is f06ejf
          Write (6,99999) niter, nconv, dnrm2(nev,d(1,2),1)
        End If
      End Do revcm

      If (ifail==0) Then
!       Post-Process using F12FCF to compute eigenvalues/values.
        ifail = 0
        Call f12fcf(nconv,d,v,ldv,sigma,resid,v,ldv,comm,icomm,ifail)
        Write (nout,99998) nconv, sigma
        Write (nout,99997)(j,d(j,1),j=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 closest to ',F8.4,' are:',/)
99997 Format (1X,I8,5X,F12.4)
    End Program f12fdfe