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

    Module d03prfe_mod

!     D03PRF 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                               :: bndary, exact, monitf, pdedef,   &
                                              uvinit
!     .. Parameters ..
      Real (Kind=nag_wp), Parameter        :: half = 0.5_nag_wp
      Real (Kind=nag_wp), Parameter        :: two = 2.0_nag_wp
      Integer, Parameter, Public           :: itrace = 0, ncode = 0, nin = 5,  &
                                              nleft = 1, nout = 6, npde = 2,   &
                                              nxfix = 0, nxi = 0
    Contains
      Subroutine uvinit(npde,npts,nxi,x,xi,u,ncode,v)

!       .. Use Statements ..
        Use nag_library, Only: x01aaf
!       .. Scalar Arguments ..
        Integer, Intent (In)                 :: ncode, npde, npts, nxi
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: u(npde,npts), v(ncode)
        Real (Kind=nag_wp), Intent (In)      :: x(npts), xi(nxi)
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: pi
        Integer                              :: i
!       .. Intrinsic Procedures ..
        Intrinsic                            :: exp, sin
!       .. Executable Statements ..
        pi = x01aaf(pi)
        Do i = 1, npts
          u(1,i) = exp(x(i))
          u(2,i) = x(i)**2 + sin(two*pi*x(i)**2)
        End Do
        Return
      End Subroutine uvinit
      Subroutine pdedef(npde,t,x,u,udot,ux,ncode,v,vdot,res,ires)

!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In)      :: t, x
        Integer, Intent (Inout)              :: ires
        Integer, Intent (In)                 :: ncode, npde
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: res(npde)
        Real (Kind=nag_wp), Intent (In)      :: u(npde), udot(npde), ux(npde), &
                                                v(ncode), vdot(ncode)
!       .. Executable Statements ..
        If (ires==-1) Then
          res(1) = udot(1)
          res(2) = udot(2)
        Else
          res(1) = udot(1) + ux(1) + ux(2)
          res(2) = udot(2) + 4.0_nag_wp*ux(1) + ux(2)
        End If
        Return
      End Subroutine pdedef
      Subroutine bndary(npde,t,ibnd,nobc,u,udot,ncode,v,vdot,res,ires)

!       .. Use Statements ..
        Use nag_library, Only: x01aaf
!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In)      :: t
        Integer, Intent (In)                 :: ibnd, ncode, nobc, npde
        Integer, Intent (Inout)              :: ires
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: res(nobc)
        Real (Kind=nag_wp), Intent (In)      :: u(npde), udot(npde), v(ncode), &
                                                vdot(ncode)
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: pp, ppt1, ppt3, t1, t3
!       .. Intrinsic Procedures ..
        Intrinsic                            :: exp, sin
!       .. Executable Statements ..
        If (ires==-1) Then
          res(1) = 0.0_nag_wp
        Else
          pp = two*x01aaf(pp)
          t1 = t
          t3 = -3.0_nag_wp*t
          If (ibnd==0) Then
            ppt3 = sin(pp*t3**2)
            ppt1 = sin(pp*t1**2)
            res(1) = u(1) - half*(exp(t3)+exp(t1)+half*(ppt3-ppt1))
            res(1) = res(1) - 2.0_nag_wp*t**2
          Else
            t3 = t3 + 1.0_nag_wp
            t1 = t1 + 1.0_nag_wp
            ppt3 = sin(pp*t3**2)
            ppt1 = sin(pp*t1**2)
            res(1) = u(2) - (exp(t3)-exp(t1)+half*(ppt3+ppt1))
            res(1) = res(1) - (1.0_nag_wp+5.0_nag_wp*t**2-2.0_nag_wp*t)
          End If
        End If
        Return
      End Subroutine bndary
      Subroutine monitf(t,npts,npde,x,u,fmon)

!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In)      :: t
        Integer, Intent (In)                 :: npde, npts
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: fmon(npts)
        Real (Kind=nag_wp), Intent (In)      :: u(npde,npts), x(npts)
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: d2x1, d2x2, h1, h2, h3
        Integer                              :: i
!       .. Intrinsic Procedures ..
        Intrinsic                            :: abs, max
!       .. Executable Statements ..
        Do i = 2, npts - 1
          h1 = x(i) - x(i-1)
          h2 = x(i+1) - x(i)
          h3 = half*(x(i+1)-x(i-1))
!         Second derivatives ..
          d2x1 = abs(((u(1,i+1)-u(1,i))/h2-(u(1,i)-u(1,i-1))/h1)/h3)
          d2x2 = abs(((u(2,i+1)-u(2,i))/h2-(u(2,i)-u(2,i-1))/h1)/h3)
          fmon(i) = max(d2x1,d2x2)
        End Do
        fmon(1) = fmon(2)
        fmon(npts) = fmon(npts-1)
        Return
      End Subroutine monitf
      Subroutine exact(t,npde,npts,x,u)

!       Exact solution (for comparison purposes)

!       .. Use Statements ..
        Use nag_library, Only: x01aaf
!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In)      :: t
        Integer, Intent (In)                 :: npde, npts
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: u(npde,npts)
        Real (Kind=nag_wp), Intent (In)      :: x(npts)
!       .. Local Scalars ..
        Real (Kind=nag_wp)                   :: pp, ppt1, ppt3, x1, x3
        Integer                              :: i
!       .. Intrinsic Procedures ..
        Intrinsic                            :: exp, sin
!       .. Executable Statements ..
        pp = 2.0_nag_wp*x01aaf(pp)
        Do i = 1, npts
          x1 = x(i) + t
          x3 = x(i) - 3.0_nag_wp*t
          ppt3 = sin(pp*x3**2)
          ppt1 = sin(pp*x1**2)
          u(1,i) = half*(exp(x3)+exp(x1)+half*(ppt3-ppt1)) - two*x(i)*t + &
            two*t**2
          u(2,i) = (exp(x3)-exp(x1)+half*(ppt3+ppt1)) - two*x(i)*t + x(i)**2 + &
            5.0_nag_wp*t**2
        End Do
        Return
      End Subroutine exact
    End Module d03prfe_mod
    Program d03prfe

!     D03PRF Example Main Program

!     .. Use Statements ..
      Use nag_library, Only: d03pek, d03prf, d03pzf, nag_wp
      Use d03prfe_mod, Only: bndary, exact, itrace, monitf, ncode, nin, nleft, &
                             nout, npde, nxfix, nxi, pdedef, uvinit
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Real (Kind=nag_wp)                   :: con, dxmesh, tout, trmesh, ts,   &
                                              xratio
      Integer                              :: i, ifail, ind, intpts, ipminf,   &
                                              it, itask, itol, itype, lenode,  &
                                              lisave, lrsave, neqn, npts,      &
                                              nrmesh, nwkres
      Logical                              :: remesh, theta
      Character (1)                        :: laopt, norm
!     .. Local Arrays ..
      Real (Kind=nag_wp)                   :: algopt(30), atol(1), rtol(1),    &
                                              xfix(1), xi(1)
      Real (Kind=nag_wp), Allocatable      :: rsave(:), u(:,:), ue(:,:),       &
                                              uout(:,:,:), x(:), xout(:)
      Integer, Allocatable                 :: isave(:)
!     .. Intrinsic Procedures ..
      Intrinsic                            :: real
!     .. Executable Statements ..
      Write (nout,*) 'D03PRF Example Program Results'
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) npts, intpts, itype
      lisave = 25 + nxfix
      neqn = npde*npts + ncode
      nwkres = npde*(npts+21+3*npde) + 7*npts + nxfix + 3
      lenode = 11*neqn + 50
      lrsave = neqn*neqn + neqn + nwkres + lenode

      Allocate (rsave(lrsave),u(npde,npts),ue(npde,npts), &
        uout(npde,intpts,itype),x(npts),xout(intpts),isave(lisave))
      Read (nin,*) itol
      Read (nin,*) atol(1), rtol(1)

!     Set remesh parameters
      remesh = .True.
      nrmesh = 3
      dxmesh = 0.0_nag_wp
      con = 5.0_nag_wp/real(npts-1,kind=nag_wp)
      xratio = 1.2_nag_wp
      ipminf = 0

!     Initialise mesh
      Do i = 1, npts
        x(i) = real(i-1,kind=nag_wp)/real(npts-1,kind=nag_wp)
      End Do
      xi(1) = 0.0_nag_wp

      Read (nin,*) xout(1:intpts)
      Read (nin,*) norm, laopt
      ind = 0
      itask = 1

!     Set theta to .TRUE. if the Theta integrator is required
      theta = .False.
      algopt(1:30) = 0.0_nag_wp
      If (theta) Then
        algopt(1) = 2.0_nag_wp
        algopt(6) = 2.0_nag_wp
        algopt(7) = 1.0_nag_wp
      End If

!     Loop over output value of t
      ts = 0.0_nag_wp
      tout = 0.0_nag_wp

      Do it = 1, 5
        tout = 0.05_nag_wp*real(it,kind=nag_wp)

!       ifail: behaviour on error exit   
!              =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft  
        ifail = 0
        Call d03prf(npde,ts,tout,pdedef,bndary,uvinit,u,npts,x,nleft,ncode, &
          d03pek,nxi,xi,neqn,rtol,atol,itol,norm,laopt,algopt,remesh,nxfix, &
          xfix,nrmesh,dxmesh,trmesh,ipminf,xratio,con,monitf,rsave,lrsave, &
          isave,lisave,itask,itrace,ind,ifail)

        If (it==1) Then
          Write (nout,99996) atol, npts
          Write (nout,99999) nrmesh
          Write (nout,99998) xout(1:intpts)
        End If

!       Interpolate at output points ..
        ifail = 0
        Call d03pzf(npde,0,u,npts,x,xout,intpts,itype,uout,ifail)

!       Check against exact solution ..
        Call exact(ts,npde,intpts,xout,ue)

        Write (nout,99997) ts
        Write (nout,99994) uout(1,1:intpts,1)
        Write (nout,99993) ue(1,1:intpts)
        Write (nout,99992) uout(2,1:intpts,1)
        Write (nout,99991) ue(2,1:intpts)

      End Do
      Write (nout,99995) isave(1), isave(2), isave(3), isave(5)

99999 Format (' Remeshing every ',I3,' time steps'/)
99998 Format (' X        ',5F10.4/)
99997 Format (' T = ',F6.3)
99996 Format (//'  Accuracy requirement =',E10.3,' Number of points = ',I3/)
99995 Format (' Number of integration steps in time = ',I6/' Number o', &
        'f function evaluations = ',I6/' Number of Jacobian eval','uations =', &
        I6/' Number of iterations = ',I6)
99994 Format (' Approx U1',5F10.4)
99993 Format (' Exact  U1',5F10.4)
99992 Format (' Approx U2',5F10.4)
99991 Format (' Exact  U2',5F10.4/)
    End Program d03prfe