! D03PJA Example Program Text
! Mark 25 Release. NAG Copyright 2014.
Module d03pjae_mod
! D03PJA 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 :: one = 1.0_nag_wp
Integer, Parameter, Public :: itrace = 0, ncode = 1, nin = 5, &
nout = 6, npde = 1, nxi = 1
Contains
Subroutine uvinit(npde,npts,x,u,ncode,v,iuser,ruser)
! Routine for PDE initial values (start time is 0.1D-6)
! .. Scalar Arguments ..
Integer, Intent (In) :: ncode, npde, npts
! .. Array Arguments ..
Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
Real (Kind=nag_wp), Intent (Out) :: u(npde,npts), v(ncode)
Real (Kind=nag_wp), Intent (In) :: x(npts)
Integer, Intent (Inout) :: iuser(*)
! .. Local Scalars ..
Real (Kind=nag_wp) :: ts
Integer :: i
! .. Intrinsic Procedures ..
Intrinsic :: exp
! .. Executable Statements ..
ts = ruser(1)
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,ncode,v,vdot,nxi,xi,ucp,ucpx,rcp,ucpt,ucptx,f, &
ires,iuser,ruser)
! .. 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) :: f(ncode)
Real (Kind=nag_wp), Intent (In) :: rcp(npde,*), ucp(npde,*), &
ucpt(npde,*), ucptx(npde,*), &
ucpx(npde,*), v(ncode), &
vdot(ncode), xi(nxi)
Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
Integer, Intent (Inout) :: iuser(*)
! .. 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,ncode,v,vdot,p,q,r,ires,iuser, &
ruser)
! .. Scalar Arguments ..
Real (Kind=nag_wp), Intent (In) :: t
Integer, Intent (Inout) :: ires
Integer, Intent (In) :: ncode, npde, nptl
! .. Array Arguments ..
Real (Kind=nag_wp), Intent (Out) :: p(npde,npde,nptl), &
q(npde,nptl), r(npde,nptl)
Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
Real (Kind=nag_wp), Intent (In) :: u(npde,nptl), ux(npde,nptl), &
v(ncode), vdot(ncode), x(nptl)
Integer, Intent (Inout) :: iuser(*)
! .. 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,ncode,v,vdot,ibnd,beta,gamma,ires,iuser, &
ruser)
! .. Scalar Arguments ..
Real (Kind=nag_wp), Intent (In) :: t
Integer, Intent (In) :: ibnd, ncode, npde
Integer, Intent (Inout) :: ires
! .. Array Arguments ..
Real (Kind=nag_wp), Intent (Out) :: beta(npde), gamma(npde)
Real (Kind=nag_wp), Intent (Inout) :: ruser(*)
Real (Kind=nag_wp), Intent (In) :: u(npde), ux(npde), v(ncode), &
vdot(ncode)
Integer, Intent (Inout) :: iuser(*)
! .. Intrinsic Procedures ..
Intrinsic :: exp
! .. Executable Statements ..
beta(1) = one
If (ibnd==0) Then
gamma(1) = -v(1)*exp(t)
Else
gamma(1) = -v(1)*vdot(1)
End If
Return
End Subroutine bndary
Subroutine exact(time,npts,x,u)
! Exact solution (for comparison purposes)
! .. Scalar Arguments ..
Real (Kind=nag_wp), Intent (In) :: time
Integer, Intent (In) :: npts
! .. Array Arguments ..
Real (Kind=nag_wp), Intent (Out) :: u(npts)
Real (Kind=nag_wp), Intent (In) :: x(npts)
! .. Local Scalars ..
Integer :: i
! .. Intrinsic Procedures ..
Intrinsic :: exp
! .. Executable Statements ..
Do i = 1, npts
u(i) = exp(time*(one-x(i))) - one
End Do
Return
End Subroutine exact
End Module d03pjae_mod
Program d03pjae
! D03PJA Example Main Program
! .. Use Statements ..
Use nag_library, Only: d03pja, nag_wp
Use d03pjae_mod, Only: bndary, exact, itrace, ncode, nin, nout, npde, &
nxi, odedef, pdedef, uvinit
! .. Implicit None Statement ..
Implicit None
! .. Local Scalars ..
Real (Kind=nag_wp) :: tout, ts
Integer :: i, ifail, ind, it, itask, itol, &
latol, lenode, lrtol, m, nbkpts, &
nel, neqn, nfuncs, niters, niw, &
njacs, npl1, npoly, npts, &
nsteps, nw, nwkres
Logical :: theta
Character (1) :: laopt, norm
! .. Local Arrays ..
Real (Kind=nag_wp) :: algopt(30), ruser(1), &
rwsav(1100), xi(nxi)
Real (Kind=nag_wp), Allocatable :: atol(:), exy(:), rtol(:), u(:), &
w(:), x(:), xbkpts(:)
Integer :: iuser(1), iwsav(505)
Integer, Allocatable :: iw(:)
Logical :: lwsav(100)
Character (80) :: cwsav(10)
! .. Intrinsic Procedures ..
Intrinsic :: mod, real
! .. Executable Statements ..
Write (nout,*) 'D03PJA 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 + ncode
npl1 = npoly + 1
nwkres = 3*npl1*npl1 + npl1*(npde*npde+6*npde+nbkpts+1) + 8*npde + &
nxi*(5*npde+1) + ncode + 3
lenode = 11*neqn + 50
nw = neqn*neqn + neqn + nwkres + lenode
niw = 25*neqn + 24
Allocate (exy(nbkpts),u(neqn),w(nw),x(npts),xbkpts(nbkpts),iw(niw))
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)
! 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.0E0_nag_wp
If (theta) Then
algopt(1) = 2.0E0_nag_wp
End If
ruser(1) = 1.0E-4_nag_wp
Write (nout,99998)
Write (nout,99997) atol(1)
Write (nout,99996) npoly
Write (nout,99995) nel
Write (nout,99994) npts
Write (nout,99999)(xbkpts(i),i=1,7,2), xbkpts(nbkpts)
! Output value solution at t = 0.1*(2**k) for k=1,2,...,5
tout = 0.1E0_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 d03pja(npde,m,ruser(1),tout,pdedef,bndary,u,nbkpts,xbkpts,npoly, &
npts,x,ncode,odedef,nxi,xi,neqn,uvinit,rtol,atol,itol,norm,laopt, &
algopt,w,nw,iw,niw,itask,itrace,ind,iuser,ruser,cwsav,lwsav,iwsav, &
rwsav,ifail)
! Check against the exact solution at alternate breakpoints.
ts = ruser(1)
Call exact(tout,nbkpts,xbkpts,exy)
Write (nout,99993) ts, (u(i),i=1,6*npoly+1,2*npoly), u(npts:neqn)
Write (nout,99992) ts, (exy(i),i=1,7,2), exy(nbkpts), ts
Write (nout,99991)
End Do
! Print integration statistics (reasonably rounded)
nsteps = 5*((iw(1)+2)/5)
nfuncs = 50*((iw(2)+25)/50)
njacs = 10*((iw(3)+5)/10)
niters = 50*((iw(5)+25)/50)
Write (nout,99990) nsteps
Write (nout,99989) nfuncs
Write (nout,99988) njacs
Write (nout,99987) niters
99999 Format (38X,'x'/12X,5F9.3,' | V(t) / v(t)'/)
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,'U(x,t=',F4.1,')',5F9.3,' | ',F8.3)
99992 Format (1X,'u(x,t=',F4.1,')',5F9.3,' | ',F8.3)
99991 Format (1X,'-----------+',44X,'-+-')
99990 Format (/' Number of time steps (nearest 5) = ',I6)
99989 Format (' Number of function evaluations (nearest 50) = ',I6)
99988 Format (' Number of Jacobian evaluations (neaerst 10) = ',I6)
99987 Format (' Number of iterations (nearest 50) = ',I6)
End Program d03pjae