!   F12FBF Example Program Text
!   Mark 25 Release. NAG Copyright 2014.

    Module f12fbfe_mod

!     F12FBF 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                               :: atv, av
!     .. Parameters ..
      Real (Kind=nag_wp), Parameter, Public :: one = 1.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, licomm = 140, nin = 5, &
                                              nout = 6
    Contains
      Subroutine av(m,n,x,w)

!       Computes  w <- A*x.

!       .. Scalar Arguments ..
        Integer, Intent (In)                 :: m, n
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: w(m)
        Real (Kind=nag_wp), Intent (In)      :: x(n)
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: h, k, s, t
        Integer                              :: i, j
!       .. Intrinsic Procedures ..
        Intrinsic                            :: real
!       .. Executable Statements ..
        h = one/real(m+1,kind=nag_wp)
        k = one/real(n+1,kind=nag_wp)
        w(1:m) = zero
        t = zero

        Do j = 1, n
          t = t + k
          s = zero
          Do i = 1, j
            s = s + h
            w(i) = w(i) + k*s*(t-one)*x(j)
          End Do
          Do i = j + 1, m
            s = s + h
            w(i) = w(i) + k*t*(s-one)*x(j)
          End Do
        End Do

        Return
      End Subroutine av

      Subroutine atv(m,n,w,y)

!       Computes  y <- A'*w.

!       .. Scalar Arguments ..
        Integer, Intent (In)                 :: m, n
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (In)      :: w(m)
        Real (Kind=nag_wp), Intent (Out)     :: y(n)
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: h, k, s, t
        Integer                              :: i, j
!       .. Intrinsic Procedures ..
        Intrinsic                            :: real
!       .. Executable Statements ..
        h = one/real(m+1,kind=nag_wp)
        k = one/real(n+1,kind=nag_wp)
        y(1:n) = zero
        t = zero

        Do j = 1, n
          t = t + k
          s = zero
          Do i = 1, j
            s = s + h
            y(j) = y(j) + k*s*(t-one)*w(i)
          End Do
          Do i = j + 1, m
            s = s + h
            y(j) = y(j) + k*t*(s-one)*w(i)
          End Do
        End Do

        Return
      End Subroutine atv
    End Module f12fbfe_mod
    Program f12fbfe

!     F12FBF Example Main Program

!     .. Use Statements ..
      Use f12fbfe_mod, Only: nout
!     .. Implicit None Statement ..
      Implicit None
!     .. Executable Statements ..
      Write (nout,*) 'F12FBF Example Program Results'

      Call ex1

      Call ex2

    Contains
      Subroutine ex1

!       .. Use Statements ..
        Use nag_library, Only: dgttrf, dgttrs, dnrm2, f12faf, f12fbf, f12fcf,  &
                               f12fdf, f12fef, nag_wp
        Use f12fbfe_mod, Only: imon, licomm, nin, one, two, zero
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: h2, 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,*)
        Write (nout,*) 'F12FBF Example 1'
        Write (nout,*)
        Flush (nout)
!       Skip heading in data file
        Read (nin,*)
        Read (nin,*)
        Read (nin,*) n, nev, ncv

        ldv = n
        lcomm = 3*n + ncv*ncv + 8*ncv + 60
        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)

!       Set the region of the spectrum that is required.
        Call f12fdf('LARGEST MAGNITUDE',icomm,comm,ifail)
!       Use the Shifted Inverse mode.
        Call f12fdf('SHIFTED INVERSE',icomm,comm,ifail)

        h2 = one/real((n+1)*(n+1),kind=nag_wp)
        sigma = zero
        ad(1:n) = two/h2 - sigma
        adl(1:n) = -one/h2
        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 .Or. irevcm==1) Then
!           Perform matrix vector multiplication
!           y <--- inv[A-sigma*I]*x
!           The NAG name equivalent of dgttrs is f07cef
            Call dgttrs('N',n,1,adl,ad,adu,adu2,ipiv,x,n,info)
          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)
            Flush (nout)
          End If
        End Do revcm

        If (ifail==0) Then
!         Post-Process using F12FCF to compute eigenvalues/vectors.
          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

        Deallocate (ad,adl,adu,adu2,comm,d,mx,resid,v,x,ipiv)

        Return

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 Subroutine ex1
      Subroutine ex2

!       .. Use Statements ..
        Use nag_library, Only: daxpy, dnrm2, dscal, f12faf, f12fbf, f12fcf,    &
                               nag_wp, x04caf
        Use f12fbfe_mod, Only: atv, av, licomm, nin, one
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: sigma, temp
        Integer                              :: ifail, ifail1, irevcm, j,      &
                                                lcomm, ldu, ldv, m, n, nconv,  &
                                                ncv, nev, nshift
!       .. Local Arrays ..
        Real (Kind=nag_wp), Allocatable      :: ax(:), comm(:), d(:,:), mx(:), &
                                                resid(:), u(:,:), v(:,:), x(:)
        Integer                              :: icomm(licomm)
!       .. Intrinsic Procedures ..
        Intrinsic                            :: sqrt
!       .. Executable Statements ..
        Write (nout,*)
        Write (nout,*) 'F12FBF Example 2'
        Write (nout,*)
        Flush (nout)
!       Skip heading in data file
        Read (nin,*)
        Read (nin,*) m, n, nev, ncv

        ldu = m
        ldv = n
        lcomm = 3*n + ncv*ncv + 8*ncv + 60
        Allocate (ax(m),comm(lcomm),d(ncv,2),mx(m),resid(n),u(ldu,nev), &
          v(ldv,ncv),x(m))

!       Initialize for problem.
        ifail = 0
        Call f12faf(n,nev,ncv,icomm,licomm,comm,lcomm,ifail)

!       Main reverse communication loop.
        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 .Or. irevcm==1) Then
!           Perform the operation X <- A'AX
            Call av(m,n,x,ax)
            Call atv(m,n,ax,x)
          End If
        End Do revcm

        If (ifail==0) Then
!         Post-Process using F12FCF.
!         Computed singular values may be extracted.
!         Singular vectors may also be computed now if desired.
          ifail1 = 0
          Call f12fcf(nconv,d,v,ldv,sigma,resid,v,ldv,comm,icomm,ifail1)

!         Singular values (squared) are returned in the first column
!         of D and the corresponding right singular vectors are
!         returned in the first NEV columns of V.

          Do j = 1, nconv
            d(j,1) = sqrt(d(j,1))

!           Compute the left singular vectors from the formula
!           u = Av/sigma/norm(Av).

            Call av(m,n,v(1,j),ax)
            u(1:m,j) = ax(1:m)
!           The NAG name equivalent of dnrm2 is f06ejf
            temp = one/dnrm2(m,u(1,j),1)
!           The NAG name equivalent of dscal is f06edf
            Call dscal(m,temp,u(1,j),1)

!           Compute the residual norm ||A*v - sigma*u|| for the nconv
!           accurately computed singular values and vectors.
!           Store the result in 2nd column of array D.

!           The NAG name equivalent of daxpy is f06ecf
            Call daxpy(m,-d(j,1),u(1,j),1,ax,1)
            d(j,2) = dnrm2(m,ax,1)

          End Do

!         Print computed residuals

          Call x04caf('G','N',nconv,2,d,ncv, &
            ' Singular values and direct residuals',ifail1)
        End If

        Deallocate (ax,comm,d,mx,resid,u,v,x)

        Return

      End Subroutine ex2
    End Program f12fbfe