Example description
!   D03PJF Example Program Text
!   Mark 26.1 Release. NAG Copyright 2017.

    Module d03pjfe_mod

!     D03PJF 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, odedef, pdedef, uvinit
!     .. Parameters ..
      Real (Kind=nag_wp), Parameter    :: one = 1.0_nag_wp
      Integer, Parameter, Public       :: itrace = 0, nin = 5, nout = 6,       &
                                          npde = 1, nv = 1, nxi = 1
      Logical, Parameter, Public       :: print_stat = .False.
!     .. Local Scalars ..
      Real (Kind=nag_wp), Public, Save :: ts
    Contains
      Subroutine uvinit(npde,npts,x,u,nv,v)

!       Routine for PDE initial values (start time is 0.1D-6)

!       .. Scalar Arguments ..
        Integer, Intent (In)           :: npde, npts, nv
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out) :: u(npde,npts), v(nv)
        Real (Kind=nag_wp), Intent (In) :: x(npts)
!       .. Local Scalars ..
        Integer                        :: i
!       .. Intrinsic Procedures ..
        Intrinsic                      :: exp
!       .. Executable Statements ..
        v(1) = ts
        Do i = 1, npts
          u(1,i) = exp(ts*(one-x(i))) - one
        End Do
        Return
      End Subroutine uvinit
      Subroutine odedef(npde,t,nv,v,vdot,nxi,xi,ucp,ucpx,rcp,ucpt,ucptx,f,     &
        ires)

!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In) :: t
        Integer, Intent (Inout)        :: ires
        Integer, Intent (In)           :: npde, nv, nxi
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out) :: f(nv)
        Real (Kind=nag_wp), Intent (In) :: rcp(npde,nxi), ucp(npde,nxi),       &
                                          ucpt(npde,nxi), ucptx(npde,nxi),     &
                                          ucpx(npde,nxi), v(nv), vdot(nv),     &
                                          xi(nxi)
!       .. Executable Statements ..
        If (ires==1) Then
          f(1) = vdot(1) - v(1)*ucp(1,1) - ucpx(1,1) - one - t
        Else If (ires==-1) Then
          f(1) = vdot(1)
        End If
        Return
      End Subroutine odedef
      Subroutine pdedef(npde,t,x,nptl,u,ux,nv,v,vdot,p,q,r,ires)

!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In) :: t
        Integer, Intent (Inout)        :: ires
        Integer, Intent (In)           :: npde, nptl, nv
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out) :: p(npde,npde,nptl), q(npde,nptl),   &
                                          r(npde,nptl)
        Real (Kind=nag_wp), Intent (In) :: u(npde,nptl), ux(npde,nptl), v(nv), &
                                          vdot(nv), x(nptl)
!       .. Local Scalars ..
        Integer                        :: i
!       .. Executable Statements ..
        Do i = 1, nptl
          p(1,1,i) = v(1)*v(1)
          r(1,i) = ux(1,i)
          q(1,i) = -x(i)*ux(1,i)*v(1)*vdot(1)
        End Do
        Return
      End Subroutine pdedef
      Subroutine bndary(npde,t,u,ux,nv,v,vdot,ibnd,beta,gamma,ires)

!       .. Scalar Arguments ..
        Real (Kind=nag_wp), Intent (In) :: t
        Integer, Intent (In)           :: ibnd, npde, nv
        Integer, Intent (Inout)        :: ires
!       .. Array Arguments ..
        Real (Kind=nag_wp), Intent (Out) :: beta(npde), gamma(npde)
        Real (Kind=nag_wp), Intent (In) :: u(npde), ux(npde), v(nv), vdot(nv)
!       .. Intrinsic Procedures ..
        Intrinsic                      :: exp
!       .. Executable Statements ..
        beta(1) = 1.0E0_nag_wp
        If (ibnd==0) Then
          gamma(1) = -v(1)*exp(t)
        Else
          gamma(1) = -v(1)*vdot(1)
        End If
        Return
      End Subroutine bndary
    End Module d03pjfe_mod
    Program d03pjfe

!     D03PJF Example Main Program

!     .. Use Statements ..
      Use d03pjfe_mod, Only: bndary, itrace, nin, nout, npde, nv, nxi, odedef, &
                             pdedef, print_stat, ts, uvinit
      Use nag_library, Only: d03pjf, nag_wp
!     .. Implicit None Statement ..
      Implicit None
!     .. Local Scalars ..
      Real (Kind=nag_wp)               :: tout
      Integer                          :: i, ifail, ind, it, itask, itol,      &
                                          latol, lenode, lisave, lrsave,       &
                                          lrtol, m, nbkpts, nel, neqn, np1,    &
                                          npoly, npts, nwkres
      Logical                          :: theta
      Character (1)                    :: laopt, norm
!     .. Local Arrays ..
      Real (Kind=nag_wp)               :: algopt(30), xi(nxi)
      Real (Kind=nag_wp), Allocatable  :: atol(:), rsave(:), rtol(:), u(:),    &
                                          x(:), xbkpts(:)
      Integer, Allocatable             :: isave(:)
!     .. Intrinsic Procedures ..
      Intrinsic                        :: abs, mod, real
!     .. Executable Statements ..
      Write (nout,*) 'D03PJF Example Program Results'
!     Skip heading in data file
      Read (nin,*)
      Read (nin,*) m, nbkpts, npoly

      nel = nbkpts - 1
      npts = nel*npoly + 1
      neqn = npde*npts + nv
      np1 = npoly + 1
      nwkres = np1*(3*np1+npde*npde+6*npde+nbkpts+1)
      nwkres = nwkres + 8*npde + nxi*(5*npde+1) + nv + 3
      lenode = 11*neqn + 50
      lisave = 25*neqn + 24
      lrsave = neqn*neqn + neqn + nwkres + lenode
      Allocate (u(neqn),rsave(lrsave),x(npts),xbkpts(nbkpts),isave(lisave))

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

      Read (nin,*) ts

!     Set break-points
      Do i = 1, nbkpts
        xbkpts(i) = real(i-1,kind=nag_wp)/real(nbkpts-1,kind=nag_wp)
      End Do

      Read (nin,*) xi(1:nxi)
      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
      End If

      Write (nout,99998)
      Write (nout,99997) atol
      Write (nout,99996) npoly
      Write (nout,99995) nel
      Write (nout,99994) npts
      Write (nout,99999)

!     Output value solution at t = 0.1*(2**k) for k=1,2,...,5

      tout = 0.1_nag_wp
      Do it = 1, 5
        tout = tout + tout

!       ifail: behaviour on error exit
!              =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft
        ifail = 0
        Call d03pjf(npde,m,ts,tout,pdedef,bndary,u,nbkpts,xbkpts,npoly,npts,x, &
          nv,odedef,nxi,xi,neqn,uvinit,rtol,atol,itol,norm,laopt,algopt,rsave, &
          lrsave,isave,lisave,itask,itrace,ind,ifail)
        If (abs(u(1))<10.0_nag_wp) Then
          Write (nout,99993) ts, u(1)
        Else
          Write (nout,99992) ts, u(1)
        End If
      End Do

      If (print_stat) Then
        Write (nout,*)
        Write (nout,99991) 'time steps', isave(1)
        Write (nout,99991) 'function evaluations', isave(2)
        Write (nout,99991) 'Jacobian evaluations', isave(3)
        Write (nout,99991) 'iterations', isave(5)
      End If

99999 Format (3X,'time',8X,'solution at x=0')
99998 Format (/,/,' Simple coupled PDE using BDF')
99997 Format (' Accuracy requirement  = ',1P,E12.3)
99996 Format (' Degree of Polynomial  = ',I4)
99995 Format (' Number of elements    = ',I4)
99994 Format (' Number of mesh points = ',I4,/)
99993 Format (1X,F6.1,10X,F7.2)
99992 Format (1X,F6.1,10X,F6.1)
99991 Format (' Number of ',A20,' = ',I6)
    End Program d03pjfe