Keyword: 対流, 拡散, 偏微分方程式
概要
本サンプルは保存型のソース項を用いた対流・拡散偏微分方程式を解くFortranによるサンプルプログラムです。 本サンプルは以下に示される2つの例の対流・拡散偏微分方程式を解いて結果を出力します。
※本サンプルはnAG Fortranライブラリに含まれるルーチン d03plf() のExampleコードです。本サンプル及びルーチンの詳細情報は d03plf のマニュアルページをご参照ください。
ご相談やお問い合わせはこちらまで
入力データ
(本ルーチンの詳細はd03plf のマニュアルページを参照)このデータをダウンロード |
D03PLF Example Program Data 141 : npts 1 : itol '1' : norm 0.1E-4 0.25E-3 : atol(1), rtol(1) 141 14 : npts, nskip 2.5 0.25 1.4 1.0 0.125 : el0, er0, gamma, rl0, rr0 1 : itol '2' : norm 0.5E-2 0.5E-3 : atol(1), rtol(1) D, V, P at selected output pts. For t = 0.1: 1.0000 0.0000 1.0000 1.0000 0.0000 1.0000 0.8775 0.1527 0.8327 0.4263 0.9275 0.3031 0.2656 0.9275 0.3031 0.1250 0.0000 0.1000 0.1250 0.0000 0.1000 0.1250 0.0000 0.1000 For t = 0.2: 1.0000 0.0000 1.0000 0.8775 0.1527 0.8327 0.6029 0.5693 0.4925 0.4263 0.9275 0.3031 0.4263 0.9275 0.3031 0.2656 0.9275 0.3031 0.2656 0.9275 0.3031 0.1250 0.0000 0.1000
- 1行目はタイトル行で読み飛ばされます。
- 2〜5行目は例1の入力データです。
- 2行目はメッシュ点の数(npts)を指定しています。
- 3行目は局所誤差チェックの形式(itol=1:相対局所誤差の許容値がスカラーで絶対局所誤差の許容値もスカラー)を指定しています。
- 4行目は使用されるノルムの種類(norm='1':平均L1ノルム)を指定しています。
- 5行目は絶対局所誤差の許容値(atol(1))と相対局所誤差の許容値(rtol(1))を指定しています。
- 7〜29行目は例2の入力データです。
- 7行目はメッシュ点の数(npts)と刻み幅(nskip)を指定しています。
- 8行目はx<0.5 のときの比エネルギー(el0)、x>0.5 のときの比エネルギー(er0)、比熱比(gamma)、x<0.5 のときの密度(rl0)、x>0.5 のときの密度(r0)を指定しています。
- 9行目は局所誤差チェックの形式(itol=1:相対局所誤差の許容値がスカラーで絶的局所誤差の許容値もスカラー)を指定しています。
- 10行目は使用されるノルムの種類(norm='2':平均L2ノルム)を指定しています。
- 11行目は絶的局所誤差の許容値(atol(1))と相対局所誤差の許容値(rtol(1))を指定しています。
- 12行目は入力データの説明文で読み飛ばされます。t=0.1の場合のデータであることを示しています。
- 13〜20行目にはt=0.1の場合の密度(d)、速度(v)、圧力(p)を指定しています。
- 21行目には入力データの説明文で読み飛ばされます。t=0.2の場合のデータであることを示しています。
- 22〜29行目にはt=0.2の場合の密度(d)、速度(v)、圧力(p)を指定しています。
出力結果
(本ルーチンの詳細はd03plf のマニュアルページを参照)この出力例をダウンロード |
D03PLF Example Program Results Example 1 NPTS = 141 ATOL = 0.100E-04 RTOL = 0.250E-03 T = 0.500 X Approx U1 Exact U1 Approx U2 Exact U2 0.0000 -0.0432 -0.0432 0.0432 0.0432 0.1429 -0.0221 -0.0220 0.0001 0.0000 0.2857 -0.0199 -0.0199 -0.0231 -0.0231 0.4286 -0.0123 -0.0123 -0.0176 -0.0176 0.5714 0.0247 0.0245 0.0226 0.0224 0.7143 0.0834 0.0827 0.0830 0.0825 0.8571 0.1042 0.1036 0.1044 0.1039 1.0000 -0.0008 -0.0001 -0.0007 0.0001 Number of integration steps in time = 159 Number of function evaluations = 1199 Number of Jacobian evaluations = 19 Number of iterations = 426 Example 2 GAMMA = 1.400 EL0 = 2.500 ER0 = 0.250 RL0 = 1.000 RR0 = 0.125 NPTS = 141 ATOL = 0.500E-02 RTOL = 0.500E-03 X APPROX D EXACT D APPROX V EXACT V APPROX P EXACT P T = 0.100 0.2000 1.0000 1.0000 0.0000 0.0000 1.0000 1.0000 0.3000 1.0000 1.0000 0.0000 0.0000 1.0000 1.0000 0.4000 0.8668 0.8775 0.1665 0.1527 0.8188 0.8327 0.5000 0.4299 0.4263 0.9182 0.9275 0.3071 0.3031 0.6000 0.2969 0.2656 0.9274 0.9275 0.3028 0.3031 0.7000 0.1250 0.1250 0.0000 0.0000 0.1000 0.1000 0.8000 0.1250 0.1250 0.0000 0.0000 0.1000 0.1000 0.9000 0.1250 0.1250 0.0000 0.0000 0.1000 0.1000 T = 0.200 0.2000 1.0000 1.0000 0.0000 0.0000 1.0000 1.0000 0.3000 0.8718 0.8775 0.1601 0.1527 0.8253 0.8327 0.4000 0.6113 0.6029 0.5543 0.5693 0.5022 0.4925 0.5000 0.4245 0.4263 0.9314 0.9275 0.3014 0.3031 0.6000 0.4259 0.4263 0.9277 0.9275 0.3030 0.3031 0.7000 0.2772 0.2656 0.9272 0.9275 0.3031 0.3031 0.8000 0.2657 0.2656 0.9276 0.9275 0.3032 0.3031 0.9000 0.1250 0.1250 0.0000 0.0000 0.1000 0.1000 Number of integration steps in time = 170 Number of function evaluations = 411 Number of Jacobian evaluations = 1 Number of iterations = 2
- 4〜24行目には例1の計算結果が出力されています。
- 7行目にはメッシュ点の数、絶対誤差許容値、相対誤差許容値が出力されています。
- 9行目にはtの値が出力されています。
- 10〜19行目には空間変数xの値、U1(x,t)の近似解と厳密解、U2(x,t)の近似解と厳密解が出力されています。
- 21行目には積分のステップの数が出力されています。
- 22行目には関数評価の数が出力されています。
- 23行目にはヤコビアン評価の数が出力されています。
- 24行目には反復数が出力されています。
- 28〜61行目には例2の計算結果が出力されています。
- 30行目にはパラメータの値が出力されています。
- 32行目にはメッシュ点の数、絶対誤差許容値、相対誤差許容値が出力されています。
- 36〜45行目にはt=0.100の場合の空間変数xの値、密度の近似解と厳密解、速度の近似解と厳密解、圧力の近似解と厳密解が出力されています。
- 47〜56行目にはt=0.200の場合の空間変数xの値、密度の近似解と厳密解、速度の近似解と厳密解、圧力の近似解と厳密解が出力されています
- 58行目には積分のステップの数が出力されています。
- 59行目には関数評価の数が出力されています。
- 60行目にはヤコビアン評価の数が出力されています。
- 61行目には反復数が出力されています。
ソースコード
(本ルーチンの詳細はd03plf のマニュアルページを参照)
※本サンプルソースコードは科学技術・統計計算ライブラリである「nAG Fortranライブラリ」のルーチンを呼び出します。
サンプルのコンパイル及び実行方法
このソースコードをダウンロード |
! D03PLF Example Program Text ! Mark 23 Release. nAG Copyright 2011. MODULE d03plfe_mod ! D03PLF Example Program Module: ! Parameters and User-defined Routines ! .. Use Statements .. USE nag_library, ONLY : nag_wp ! .. Implicit None Statement .. IMPLICIT NONE ! .. Parameters .. INTEGER, PARAMETER :: itrace = 0, ncode1 = 2, & ncode2 = 0, nin = 5, nout = 6, & npde1 = 2, npde2 = 3, nxi1 = 2, & nxi2 = 0 ! .. Local Scalars .. REAL (KIND=nag_wp), SAVE :: el0, er0, gamma, rl0, rr0 CONTAINS SUBROUTINE exact(t,u,npde,x,npts) ! Exact solution (for comparison and b.c. purposes) ! .. Use Statements .. USE nag_library, ONLY : x01aaf ! .. Implicit None Statement .. IMPLICIT NONE ! .. 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(*) ! .. Local Scalars .. REAL (KIND=nag_wp) :: f, g, pi, pi2, x1, x3 INTEGER :: i ! .. Intrinsic Functions .. INTRINSIC cos, exp, sin ! .. Executable Statements .. f = 0.0_nag_wp pi = x01aaf(f) pi2 = 2.0_nag_wp*pi DO i = 1, npts x1 = x(i) + t x3 = x(i) - 3.0_nag_wp*t f = exp(pi*x3)*sin(pi2*x3) g = exp(-pi2*x1)*cos(pi2*x1) u(1,i) = f + g u(2,i) = f - g END DO RETURN END SUBROUTINE exact SUBROUTINE pdedef(npde,t,x,u,ux,ncode,v,vdot,p,c,d,s,ires) ! .. Implicit None Statement .. IMPLICIT NONE ! .. 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) :: c(npde), d(npde), & p(npde,npde), s(npde) REAL (KIND=nag_wp), INTENT (IN) :: u(npde), ux(npde), v(ncode), & vdot(ncode) ! .. Local Scalars .. INTEGER :: i ! .. Executable Statements .. c(1:npde) = 1.0_nag_wp d(1:npde) = 0.0_nag_wp s(1:npde) = 0.0_nag_wp p(1:npde,1:npde) = 0.0_nag_wp DO i = 1, npde p(i,i) = 1.0_nag_wp END DO RETURN END SUBROUTINE pdedef SUBROUTINE bndry1(npde,npts,t,x,u,ncode,v,vdot,ibnd,g,ires) ! .. Implicit None Statement .. IMPLICIT NONE ! .. Scalar Arguments .. REAL (KIND=nag_wp), INTENT (IN) :: t INTEGER, INTENT (IN) :: ibnd, ncode, npde, npts INTEGER, INTENT (INOUT) :: ires ! .. Array Arguments .. REAL (KIND=nag_wp), INTENT (OUT) :: g(npde) REAL (KIND=nag_wp), INTENT (IN) :: u(npde,npts), v(ncode), & vdot(ncode), x(npts) ! .. Local Scalars .. REAL (KIND=nag_wp) :: dudx INTEGER :: i ! .. Local Arrays .. REAL (KIND=nag_wp) :: ue(2,1) ! .. Executable Statements .. IF (ibnd==0) THEN i = 1 CALL exact(t,ue,npde,x(i),1) g(1) = u(1,i) + u(2,i) - ue(1,1) - ue(2,1) dudx = (u(1,i+1)-u(2,i+1)-u(1,i)+u(2,i))/(x(i+1)-x(i)) g(2) = vdot(1) - dudx ELSE i = npts CALL exact(t,ue,npde,x(i),1) g(1) = u(1,i) - u(2,i) - ue(1,1) + ue(2,1) dudx = (u(1,i)+u(2,i)-u(1,i-1)-u(2,i-1))/(x(i)-x(i-1)) g(2) = vdot(2) + 3.0_nag_wp*dudx END IF RETURN END SUBROUTINE bndry1 SUBROUTINE nmflx1(npde,t,x,ncode,v,uleft,uright,flux,ires) ! .. Implicit None Statement .. IMPLICIT NONE ! .. 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) :: flux(npde) REAL (KIND=nag_wp), INTENT (IN) :: uleft(npde), uright(npde), & v(ncode) ! .. Local Scalars .. REAL (KIND=nag_wp) :: tmpl, tmpr ! .. Executable Statements .. tmpl = 3.0_nag_wp*(uleft(1)+uleft(2)) tmpr = uright(1) - uright(2) flux(1) = 0.5_nag_wp*(tmpl-tmpr) flux(2) = 0.5_nag_wp*(tmpl+tmpr) RETURN END SUBROUTINE nmflx1 SUBROUTINE odedef(npde,t,ncode,v,vdot,nxi,xi,ucp,ucpx,ucpt,r,ires) ! .. Implicit None Statement .. IMPLICIT NONE ! .. 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) = 0.0_nag_wp r(2) = 0.0_nag_wp ELSE r(1) = v(1) - ucp(1,1) + ucp(2,1) r(2) = v(2) - ucp(1,2) - ucp(2,2) END IF RETURN END SUBROUTINE odedef SUBROUTINE bndry2(npde,npts,t,x,u,ncode,v,vdot,ibnd,g,ires) ! .. Implicit None Statement .. IMPLICIT NONE ! .. Scalar Arguments .. REAL (KIND=nag_wp), INTENT (IN) :: t INTEGER, INTENT (IN) :: ibnd, ncode, npde, npts INTEGER, INTENT (INOUT) :: ires ! .. Array Arguments .. REAL (KIND=nag_wp), INTENT (OUT) :: g(npde) REAL (KIND=nag_wp), INTENT (IN) :: u(npde,npts), v(ncode), & vdot(ncode), x(npts) ! .. Executable Statements .. IF (ibnd==0) THEN g(1) = u(1,1) - rl0 g(2) = u(2,1) g(3) = u(3,1) - el0 ELSE g(1) = u(1,npts) - rr0 g(2) = u(2,npts) g(3) = u(3,npts) - er0 END IF RETURN END SUBROUTINE bndry2 SUBROUTINE nmflx2(npde,t,x,ncode,v,uleft,uright,flux,ires) ! .. Use Statements .. USE nag_library, ONLY : d03puf, d03pvf ! .. Implicit None Statement .. IMPLICIT NONE ! .. 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) :: flux(npde) REAL (KIND=nag_wp), INTENT (IN) :: uleft(npde), uright(npde), & v(ncode) ! .. Local Scalars .. INTEGER :: ifail CHARACTER (1) :: path, solver ! .. Executable Statements .. ifail = 0 solver = 'R' IF (solver=='R') THEN ! ROE scheme .. CALL d03puf(uleft,uright,gamma,flux,ifail) ELSE ! OSHER scheme .. path = 'P' CALL d03pvf(uleft,uright,gamma,path,flux,ifail) END IF RETURN END SUBROUTINE nmflx2 SUBROUTINE uvinit(npde,npts,x,u) ! .. Implicit None Statement .. IMPLICIT NONE ! .. Scalar Arguments .. 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 .. INTEGER :: i ! .. Executable Statements .. DO i = 1, npts IF (x(i)<0.5_nag_wp) THEN u(1,i) = rl0 u(2,i) = 0.0_nag_wp u(3,i) = el0 ELSE IF (x(i)==0.5_nag_wp) THEN u(1,i) = 0.5_nag_wp*(rl0+rr0) u(2,i) = 0.0_nag_wp u(3,i) = 0.5_nag_wp*(el0+er0) ELSE u(1,i) = rr0 u(2,i) = 0.0_nag_wp u(3,i) = er0 END IF END DO RETURN END SUBROUTINE uvinit END MODULE d03plfe_mod PROGRAM d03plfe ! D03PLF Example Main Program ! .. Use Statements .. USE d03plfe_mod, ONLY : nout ! .. Implicit None Statement .. IMPLICIT NONE ! .. External Subroutines .. EXTERNAL ex1, ex2 ! .. Executable Statements .. WRITE (nout,*) 'D03PLF Example Program Results' CALL ex1 CALL ex2 END PROGRAM d03plfe SUBROUTINE ex1 ! .. Use Statements .. USE nag_library, ONLY : d03plf USE d03plfe_mod, ONLY : bndry1, exact, itrace, nag_wp, ncode1, nin, & nmflx1, nout, npde1, nxi1, odedef, pdedef ! .. Implicit None Statement .. IMPLICIT NONE ! .. Local Scalars .. REAL (KIND=nag_wp) :: tout, ts INTEGER :: i, ifail, ind, itask, itol, & lenode, lisave, lrsave, ncode, & neqn, nop, npde, npts, nwkres, & nxi, outpts CHARACTER (1) :: laopt, norm ! .. Local Arrays .. REAL (KIND=nag_wp) :: algopt(30), atol(1), rtol(1) REAL (KIND=nag_wp), ALLOCATABLE :: rsave(:), u(:), ue(:,:), & uout(:,:), x(:), xi(:), xout(:) INTEGER, ALLOCATABLE :: isave(:) ! .. Intrinsic Functions .. INTRINSIC real ! .. Executable Statements .. WRITE (nout,*) WRITE (nout,*) WRITE (nout,*) 'Example 1' WRITE (nout,*) ! Skip heading in data file READ (nin,*) READ (nin,*) npts npde = npde1 ncode = ncode1 nxi = nxi1 neqn = npde*npts + ncode lisave = 25*neqn + 24 nwkres = npde*(2*npts+6*nxi+3*npde+26) + nxi + ncode + 7*npts + 2 lenode = 11*neqn + 50 lrsave = 4*neqn + 11*neqn/2 + 1 + nwkres + lenode lisave = lisave*4 lrsave = lrsave*4 outpts = npts/20 + 1 ALLOCATE (rsave(lrsave),u(neqn),ue(npde,outpts),uout(npde,outpts), & x(npts),xi(nxi),xout(outpts),isave(lisave)) READ (nin,*) itol READ (nin,*) norm READ (nin,*) atol(1), rtol(1) ! 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 xi(2) = 1.0_nag_wp ! Set initial values .. ts = 0.0_nag_wp CALL exact(ts,u,npde,x,npts) u(neqn-1) = u(1) - u(2) u(neqn) = u(neqn-2) + u(neqn-3) laopt = 'S' ind = 0 itask = 1 algopt(1:30) = 0.0_nag_wp ! Theta integration algopt(1) = 1.0_nag_wp ! Sparse matrix algebra parameters algopt(29) = 0.1_nag_wp algopt(30) = 1.1_nag_wp tout = 0.5_nag_wp ! ifail: behaviour on error exit ! =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft ifail = 0 CALL d03plf(npde,ts,tout,pdedef,nmflx1,bndry1,u,npts,x,ncode,odedef, & nxi,xi,neqn,rtol,atol,itol,norm,laopt,algopt,rsave,lrsave,isave, & lisave,itask,itrace,ind,ifail) WRITE (nout,99995) npts, atol, rtol ! Set output points .. nop = 0 DO i = 1, npts, 20 nop = nop + 1 xout(nop) = x(i) END DO WRITE (nout,99996) ts WRITE (nout,99999) nop = 0 DO i = 1, npts, 20 nop = nop + 1 uout(1,nop) = u(npde*(i-1)+1) uout(2,nop) = u(npde*(i-1)+2) END DO ! Check against exact solution .. CALL exact(tout,ue,npde,xout,nop) WRITE (nout,99998) (xout(i),uout(1,i),ue(1,i),uout(2,i),ue(2,i),i=1, & nop) WRITE (nout,99997) WRITE (nout,99994) isave(1), isave(2), isave(3), isave(5) RETURN 99999 FORMAT (8X,'X',8X,'Approx U1',3X,'Exact U1',4X,'Approx U2',3X, & 'Exact U2'/) 99998 FORMAT (5F12.4) 99997 FORMAT (E11.4,4E14.4) 99996 FORMAT (' T = ',F6.3) 99995 FORMAT (/' NPTS = ',I4,' ATOL = ',E10.3,' RTOL = ',E10.3/) 99994 FORMAT (' Number of integration steps in time = ',I6/' Number ', & 'of function evaluations = ',I6/' Number of Jacobian ', & 'evaluations =',I6/' Number of iterations = ',I6) END SUBROUTINE ex1 SUBROUTINE ex2 ! .. Use Statements .. USE nag_library, ONLY : d03pek, d03plf, d03plp USE d03plfe_mod, ONLY : bndry2, el0, er0, gamma, itrace, nag_wp, & ncode2, nin, nmflx2, nout, npde2, nxi2, rl0, & rr0, uvinit ! .. Implicit None Statement .. IMPLICIT NONE ! .. Local Scalars .. REAL (KIND=nag_wp) :: d, p, tout, ts, v INTEGER :: i, ifail, ind, it, itask, itol, & k, lenode, lisave, lrsave, mlu, & ncode, neqn, npde, npts, nskip, & nwkres, nxi, outpts CHARACTER (1) :: laopt, norm ! .. Local Arrays .. REAL (KIND=nag_wp) :: algopt(30), atol(1), rtol(1), & xi(1) REAL (KIND=nag_wp), ALLOCATABLE :: rsave(:), u(:,:), ue(:,:), x(:) INTEGER, ALLOCATABLE :: isave(:) ! .. Intrinsic Functions .. INTRINSIC real ! .. Executable Statements .. WRITE (nout,*) WRITE (nout,*) WRITE (nout,*) 'Example 2' WRITE (nout,*) READ (nin,*) READ (nin,*) npts, nskip npde = npde2 ncode = ncode2 nxi = nxi2 READ (nin,*) el0, er0, gamma, rl0, rr0 nwkres = npde*(2*npts+3*npde+32) + 7*npts + 4 mlu = 3*npde - 1 neqn = npde*npts + ncode lenode = 9*neqn + 50 lisave = neqn + 24 lrsave = (3*mlu+1)*neqn + nwkres + lenode outpts = (npts-2*nskip-1)/nskip ALLOCATE (rsave(lrsave),u(npde,npts),x(npts),isave(lisave), & ue(npde,outpts)) READ (nin,*) itol READ (nin,*) norm READ (nin,*) atol(1), rtol(1) WRITE (nout,99995) gamma, el0, er0, rl0, rr0 WRITE (nout,99997) npts, atol, rtol ! Initialise mesh DO i = 1, npts x(i) = real(i-1,kind=nag_wp)/real(npts-1,kind=nag_wp) END DO ! Initial values of variables CALL uvinit(npde,npts,x,u) xi(1) = 0.0_nag_wp laopt = 'B' ind = 0 itask = 1 algopt(1:30) = 0.0_nag_wp ! Theta integration algopt(1) = 2.0_nag_wp algopt(6) = 2.0_nag_wp algopt(7) = 2.0_nag_wp ! Max. time step algopt(13) = 0.5E-2_nag_wp ts = 0.0_nag_wp WRITE (nout,99999) DO it = 1, 2 tout = real(it,kind=nag_wp)*0.1_nag_wp ! ifail: behaviour on error exit ! =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft ifail = 0 CALL d03plf(npde,ts,tout,d03plp,nmflx2,bndry2,u,npts,x,ncode,d03pek, & nxi,xi,neqn,rtol,atol,itol,norm,laopt,algopt,rsave,lrsave,isave, & lisave,itask,itrace,ind,ifail) WRITE (nout,99998) ts ! Read exact data (to 4 d.p.) at output points .. READ (nin,*) DO i = 1, outpts READ (nin,*) ue(1,i), ue(2,i), ue(3,i) END DO ! Calculate density, velocity and pressure .. k = 0 DO i = 2*nskip + 1, npts - nskip, nskip d = u(1,i) v = u(2,i)/d p = d*(gamma-1.0_nag_wp)*(u(3,i)/d-0.5_nag_wp*v**2) k = k + 1 WRITE (nout,99994) x(i), d, ue(1,k), v, ue(2,k), p, ue(3,k) END DO END DO WRITE (nout,99996) isave(1), isave(2), isave(3), isave(5) RETURN 99999 FORMAT (4X,'X',5X,'APPROX D',1X,'EXACT D',2X,'APPROX V',1X,'EXAC', & 'T V',2X,'APPROX P',1X,'EXACT P') 99998 FORMAT (/' T = ',F6.3/) 99997 FORMAT (/' NPTS = ',I4,' ATOL = ',E10.3,' RTOL = ',E10.3/) 99996 FORMAT (/' Number of integration steps in time = ',I6/' Number ', & 'of function evaluations = ',I6/' Number of Jacobian ', & 'evaluations =',I6/' Number of iterations = ',I6) 99995 FORMAT (/' GAMMA =',F6.3,' EL0 =',F6.3,' ER0 =',F6.3,' RL0 =',F6.3, & ' RR0 =',F6.3) 99994 FORMAT (1X,F7.4,6(2X,F7.4)) END SUBROUTINE ex2