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

    Module f12aefe_mod

!     F12AEF 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 :: four = 4.0_nag_wp
      Real (Kind=nag_wp), Parameter, Public :: three = 3.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       :: nin = 5, nout = 6
    Contains
      Subroutine mv(n,v)
!       Compute the in-place matrix vector multiplication X<---M*X,
!       where M is mass matrix formed by using piecewise linear elements
!       on [0,1].

!       .. Scalar Arguments ..
        Integer, Intent (In)           :: n
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Inout) :: v(n)
!       .. Local Scalars ..
        Real (Kind=nag_wp)             :: vm1, vv
        Integer                        :: j
!       .. Executable Statements ..
        vm1 = v(1)
        v(1) = four*v(1) + v(2)
        Do j = 2, n - 1
          vv = v(j)
          v(j) = vm1 + four*vv + v(j+1)
          vm1 = vv
        End Do
        v(n) = vm1 + four*v(n)
        Return
      End Subroutine mv

      Subroutine av(n,v,w)

!       .. 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 ..
        Integer                        :: j
!       .. Executable Statements ..
        w(1) = two*v(1) + three*v(2)
        Do j = 2, n - 1
          w(j) = -two*v(j-1) + two*v(j) + three*v(j+1)
        End Do
        w(n) = -two*v(n-1) + two*v(n)
        Return
      End Subroutine av
    End Module f12aefe_mod

    Program f12aefe

!     F12AEF Example Main Program

!     .. Use Statements ..
      Use f12aefe_mod, Only: av, four, mv, nin, nout, three, two, zero
      Use nag_library, Only: ddot, dnrm2, f06bnf, f12aaf, f12abf, f12acf,      &
                             f12adf, f12aef, nag_wp, zgttrf, zgttrs
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Complex (Kind=nag_wp)            :: c1, c2, c3, csig
      Real (Kind=nag_wp)               :: deni, denr, nev_nrm, numi, numr,     &
                                          sigmai, sigmar
      Integer                          :: ifail, ifail1, info, irevcm, j,      &
                                          lcomm, ldv, licomm, n, nconv, ncv,   &
                                          nev, niter, nshift
      Logical                          :: first
!     .. Local Arrays ..
      Complex (Kind=nag_wp), Allocatable :: cdd(:), cdl(:), cdu(:), cdu2(:),   &
                                          ctemp(:)
      Real (Kind=nag_wp), Allocatable  :: ax(:), comm(:), d(:,:), mx(:),       &
                                          resid(:), v(:,:), x(:)
      Integer, Allocatable             :: icomm(:), ipiv(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: cmplx, real
!     .. Executable Statements ..
      Write (nout,*) 'F12AEF Example Program Results'
      Write (nout,*)
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) n, nev, ncv, sigmar, sigmai
      ldv = n
      licomm = 140
      lcomm = 3*n + 3*ncv*ncv + 6*ncv + 60

      Allocate (cdd(n),cdl(n),cdu(n),cdu2(n),ctemp(n),ax(n),comm(lcomm),       &
        d(ncv,3),mx(n),resid(n),v(ldv,ncv),x(n),icomm(licomm),ipiv(n))

!     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('SHIFTED REAL',icomm,comm,ifail)

!     Set problem type
      Call f12adf('GENERALIZED',icomm,comm,ifail)

!     Solve A*x = lambda*B*x in shift-invert mode.
!     The shift, sigma, is a complex number (sigmar, sigmai).
!     OP = Real_Part{inv[A-(SIGMAR,SIGMAI)*M]*M and  B = M.
      csig = cmplx(sigmar,sigmai,kind=nag_wp)
      c1 = cmplx(-two,kind=nag_wp) - csig
      c2 = cmplx(two,kind=nag_wp) - cmplx(four,kind=nag_wp)*csig
      c3 = cmplx(three,kind=nag_wp) - csig

      cdl(1:n-1) = c1
      cdd(1:n-1) = c2
      cdu(1:n-1) = c3
      cdd(n) = c2

!     The NAG name equivalent of zgttrf is f07crf
      Call zgttrf(n,cdl,cdd,cdu,cdu2,ipiv,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)
!           Perform  x <--- OP*x = inv[A-SIGMA*M]*M*x
            Call mv(n,x)
            ctemp(1:n) = cmplx(x(1:n),kind=nag_wp)
!           The NAG name equivalent of zgttrs is f07csf
            Call zgttrs('N',n,1,cdl,cdd,cdu,cdu2,ipiv,ctemp,n,info)
            x(1:n) = real(ctemp(1:n))
          Case (1)
!           Perform  x <--- OP*x = inv[A-SIGMA*M]*M*x,
!           M*X stored in MX.
            ctemp(1:n) = cmplx(mx(1:n),kind=nag_wp)
!           The NAG name equivalent of zgttrs is f07csf
            Call zgttrs('N',n,1,cdl,cdd,cdu,cdu2,ipiv,ctemp,n,info)
            x(1:n) = real(ctemp(1:n))
          Case (2)
!           Perform  y <--- M*x
            Call mv(n,x)
          Case (4)
!           Output monitoring information
            Call f12aef(niter,nconv,d,d(1,2),d(1,3),icomm,comm)
!           The NAG name equivalent of dnrm2 is f06ejf
            nev_nrm = dnrm2(nev,d(1,3),1)
            Write (6,99999) niter, nconv, nev_nrm
          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)

        first = .True.
        Do j = 1, nconv
!         Use Rayleigh Quotient to recover eigenvalues of the original
!         problem.
!         The NAG name equivalent of ddot is f06eaf
          If (d(j,2)==zero) Then
!           Ritz value is real. x = v(:,j); eig = x'Ax/x'Mx.
            Call av(n,v(1,j),ax)
            numr = ddot(n,v(1,j),1,ax,1)
            mx(1:n) = v(1:n,j)
            Call mv(n,mx)
            denr = ddot(n,v(1,j),1,mx,1)
            d(j,1) = numr/denr
          Else If (first) Then
!           Ritz value is complex: x = v(:,j) - i v(:,j+1).

!           Compute x'(Ax):
!              first (xr,xi)'*(A xr)
            Call av(n,v(1,j),ax)
            numr = ddot(n,v(1,j),1,ax,1)
            numi = ddot(n,v(1,j+1),1,ax,1)
!           then add (xi,-xr)'*(A xi)
            Call av(n,v(1,j+1),ax)
            numr = numr + ddot(n,v(1,j+1),1,ax,1)
            numi = -numi + ddot(n,v(1,j),1,ax,1)

!           Compute x'(Mx) as above using mv in, place of av.
            mx(1:n) = v(1:n,j)
            Call mv(n,mx)
            denr = ddot(n,v(1,j),1,mx,1)
            deni = ddot(n,v(1,j+1),1,mx,1)
            mx(1:n) = v(1:n,j+1)
            Call mv(n,mx)
            denr = denr + ddot(n,v(1,j+1),1,mx,1)
            deni = -deni + ddot(n,v(1,j),1,mx,1)

!           Rayleigh quotient,  d=x'(Ax)/x'(Mx), (complex division).
            d(j,1) = (numr*denr+numi*deni)/f06bnf(denr,deni)
            d(j,2) = (numi*denr-numr*deni)/f06bnf(denr,deni)
            first = .False.
          Else
!           Second of complex conjugate pair.
            d(j,1) = d(j-1,1)
            d(j,2) = -d(j-1,2)
            first = .True.
          End If
        End Do
!       Print computed eigenvalues.
        Write (nout,99998) nconv, sigmar, sigmai
        Write (nout,99997)(j,d(j,1:2),j=1,nconv)
      End If

99999 Format (1X,'Iteration',1X,I3,', No. converged =',1X,I3,', norm o',       &
        'f estimates =',E12.4)
99998 Format (1X,/,' The ',I4,' generalized Ritz values closest to (',F8.4,    &
        ', ',F8.4,') are:',/)
99997 Format (1X,I8,5X,'(',F7.4,',',F7.4,')')
    End Program f12aefe