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

    Module d03pkfe_mod

!     D03PKF 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, odedef, pdedef,   &
                                              uvinit
!     .. Parameters ..
      Real (Kind=nag_wp), Parameter, Public :: one = 1.0_nag_wp
      Integer, Parameter, Public           :: itrace = 0, ncode = 1, nin = 5,  &
                                              nleft = 1, nout = 6, npde = 2,   &
                                              nxi = 1
!     .. Local Scalars ..
      Real (Kind=nag_wp), Public, Save     :: ts
    Contains
      Subroutine odedef(npde,t,ncode,v,vdot,nxi,xi,ucp,ucpx,ucpt,r,ires)

!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In)      :: t
        Integer, Intent (Inout)              :: ires
        Integer, Intent (In)                 :: ncode, npde, nxi
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: r(ncode)
        Real (Kind=nag_wp), Intent (In)      :: ucp(npde,*), ucpt(npde,*),     &
                                                ucpx(npde,*), v(ncode),        &
                                                vdot(ncode), xi(nxi)
!       .. Executable Statements ..
        If (ires==-1) Then
          r(1) = vdot(1)
        Else
          r(1) = vdot(1) - v(1)*ucp(1,1) - ucp(2,1) - one - t
        End If
        Return
      End Subroutine odedef
      Subroutine pdedef(npde,t,x,u,ut,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), ut(npde), ux(npde),   &
                                                v(ncode), vdot(ncode)
!       .. Executable Statements ..
        If (ires==-1) Then
          res(1) = v(1)*v(1)*ut(1) - x*u(2)*v(1)*vdot(1)
          res(2) = 0.0_nag_wp
        Else
          res(1) = v(1)*v(1)*ut(1) - x*u(2)*v(1)*vdot(1) - ux(2)
          res(2) = u(2) - ux(1)
        End If
        Return
      End Subroutine pdedef
      Subroutine bndary(npde,t,ibnd,nobc,u,ut,ncode,v,vdot,res,ires)

!       .. 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), ut(npde), v(ncode),   &
                                                vdot(ncode)
!       .. Intrinsic Procedures ..
        Intrinsic                            :: exp
!       .. Executable Statements ..
        If (ibnd==0) Then
          If (ires==-1) Then
            res(1) = 0.0_nag_wp
          Else
            res(1) = u(2) + v(1)*exp(t)
          End If
        Else
          If (ires==-1) Then
            res(1) = v(1)*vdot(1)
          Else
            res(1) = u(2) + v(1)*vdot(1)
          End If
        End If
        Return
      End Subroutine bndary
      Subroutine uvinit(npde,npts,x,u,ncode,neqn)

!       Routine for PDE initial values

!       .. Scalar Arguments ..
        Integer, Intent (In)                 :: ncode, neqn, npde, npts
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: u(neqn)
        Real (Kind=nag_wp), Intent (In)      :: x(npts)
!       .. Local Scalars ..
        Integer                              :: i, k
!       .. Intrinsic Procedures ..
        Intrinsic                            :: exp
!       .. Executable Statements ..
        k = 1
        Do i = 1, npts
          u(k) = exp(ts*(one-x(i))) - one
          u(k+1) = -ts*exp(ts*(one-x(i)))
          k = k + 2
        End Do
        u(neqn) = ts
        Return
      End Subroutine uvinit
      Subroutine exact(time,neqn,npts,x,u)
!       Exact solution (for comparison purposes)

!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In)      :: time
        Integer, Intent (In)                 :: neqn, npts
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out)     :: u(neqn)
        Real (Kind=nag_wp), Intent (In)      :: x(npts)
!       .. Local Scalars ..
        Integer                              :: i, k
!       .. Intrinsic Procedures ..
        Intrinsic                            :: exp
!       .. Executable Statements ..
        k = 1
        Do i = 1, npts
          u(k) = exp(time*(one-x(i))) - one
          k = k + 2
        End Do
        Return
      End Subroutine exact
    End Module d03pkfe_mod
    Program d03pkfe

!     D03PKF Example Program Text

!     .. Use Statements ..
      Use nag_library, Only: d03pkf, nag_wp
      Use d03pkfe_mod, Only: bndary, exact, itrace, ncode, nin, nleft, nout,   &
                             npde, nxi, odedef, one, pdedef, ts, uvinit
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Real (Kind=nag_wp)                   :: tout
      Integer                              :: i, ifail, ind, it, itask, itol,  &
                                              latol, lenode, lisave, lrsave,   &
                                              lrtol, neqn, npts, nwkres
      Character (1)                        :: laopt, norm
!     .. Local Arrays ..
      Real (Kind=nag_wp)                   :: algopt(30), xi(nxi)
      Real (Kind=nag_wp), Allocatable      :: atol(:), exy(:), rsave(:),       &
                                              rtol(:), u(:), x(:)
      Integer, Allocatable                 :: isave(:)
!     .. Intrinsic Procedures ..
      Intrinsic                            :: mod, real
!     .. Executable Statements ..
      Write (nout,*) 'D03PKF Example Program Results'
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) npts
      neqn = npde*npts + ncode
      nwkres = npde*(npts+6*nxi+3*npde+15) + ncode + nxi + 7*npts + 2
      lenode = 11*neqn + 50
      lisave = 25*neqn + 24
      lrsave = neqn*neqn + neqn + nwkres + lenode
      Allocate (exy(neqn),u(neqn),rsave(lrsave),x(npts),isave(lisave))

      Read (nin,*) itol
      latol = 1
      lrtol = 1
      If (itol>2) latol = neqn
      If (mod(itol,2)==0) lrtol = neqn
      Allocate (atol(latol),rtol(lrtol))
      Read (nin,*) atol(1:latol), rtol(1:lrtol)
      Read (nin,*) ts

!     Set spatial-mesh points
      Do i = 1, npts
        x(i) = real(i-1,kind=nag_wp)/real(npts-1,kind=nag_wp)
      End Do

      Read (nin,*) xi(1:nxi)
      Read (nin,*) norm, laopt
      ind = 0
      itask = 1

      algopt(1:30) = 0.0_nag_wp
      algopt(1) = one
      algopt(13) = 0.005_nag_wp

!     Loop over output value of t

      Call uvinit(npde,npts,x,u,ncode,neqn)

      tout = 0.2_nag_wp
      Do it = 1, 5
!       ifail: behaviour on error exit   
!              =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft  
        ifail = 0
        Call d03pkf(npde,ts,tout,pdedef,bndary,u,npts,x,nleft,ncode,odedef, &
          nxi,xi,neqn,rtol,atol,itol,norm,laopt,algopt,rsave,lrsave,isave, &
          lisave,itask,itrace,ind,ifail)

        If (it==1) Then
          Write (nout,99997) atol, npts
          Write (nout,99999)(x(i),i=1,13,4), x(npts)
        End If

!       Check against the exact solution

        Call exact(tout,neqn,npts,x,exy)

        Write (nout,99998) ts
        Write (nout,99995)(u(i),i=1,25,4*npde), u(neqn-2), u(neqn)
        Write (nout,99994)(exy(i),i=1,25,4*npde), exy(neqn-2), ts
        tout = 2.0_nag_wp*tout
      End Do
      Write (nout,99996) isave(1), isave(2), isave(3), isave(5)

99999 Format ('  X        ',5F9.3/)
99998 Format (' T = ',F6.3)
99997 Format (//'  Accuracy requirement =',E10.3,' Number of points = ',I3/)
99996 Format (' Number of integration steps in time = ',I6/' Number o', &
        'f function evaluations = ',I6/' Number of Jacobian eval','uations =', &
        I6/' Number of iterations = ',I6)
99995 Format (1X,'App.  sol.  ',F7.3,4F9.3,'  ODE sol. =',F8.3)
99994 Format (1X,'Exact sol.  ',F7.3,4F9.3,'  ODE sol. =',F8.3/)
    End Program d03pkfe