Keyword: 対流, 拡散, 偏微分方程式
概要
本サンプルは保存型のソース項を用いた対流・拡散偏微分方程式を解くC言語によるサンプルプログラムです。 本サンプルは以下に示される2つの例の対流・拡散偏微分方程式を解いて結果を出力します。
※本サンプルはnAG Cライブラリに含まれる関数 nag_pde_parab_1d_cd_ode() のExampleコードです。本サンプル及び関数の詳細情報は nag_pde_parab_1d_cd_ode のマニュアルページをご参照ください。
ご相談やお問い合わせはこちらまで
出力結果
(本関数の詳細はnag_pde_parab_1d_cd_ode のマニュアルページを参照)1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65
この出力例をダウンロード |
nag_pde_parab_1d_cd_ode (d03plc) Example Program Results Example 1 Method parameters: Number of mesh points used = 201 Relative tolerance used = 2.500e-04 Absolute tolerance used = 1.000e-05 Integration Results: Global error is less than 100 times the local error tolerance. Integration Statistics: Number of time steps (nearest 50) = 150 Number of function evaluations (nearest 100) = 1400 Number of Jacobian evaluations (nearest 20) = 20 Number of iterations (nearest 100) = 400 Example 2 Problem parameters and initial conditions: gamma = 1.400 e(x<0.5,0) = 2.500 e(x>0.5,0) = 0.250 rho(x<0.5,0) = 1.000 rho(x>0.5,0) = 0.125 Method parameters: Number of mesh points used = 141 Relative tolerance used = 5.000e-04 Absolute tolerance used = 5.000e-03 Solution t x d v p 0.100 0.0000 1.0000 0.0000 1.0000 0.1000 1.0000 -0.0000 1.0000 0.2000 1.0000 -0.0000 1.0000 0.3000 1.0000 -0.0000 1.0000 0.4000 0.8668 0.1665 0.8188 0.5000 0.4299 0.9182 0.3071 0.6000 0.2969 0.9274 0.3028 0.7000 0.1250 0.0000 0.1000 0.8000 0.1250 -0.0000 0.1000 0.9000 0.1250 -0.0000 0.1000 1.0000 0.1250 0.0000 0.1000 0.200 0.0000 1.0000 0.0000 1.0000 0.1000 1.0000 -0.0000 1.0000 0.2000 1.0000 -0.0000 1.0000 0.3000 0.8718 0.1601 0.8253 0.4000 0.6113 0.5543 0.5022 0.5000 0.4245 0.9314 0.3014 0.6000 0.4259 0.9277 0.3030 0.7000 0.2772 0.9272 0.3031 0.8000 0.2657 0.9276 0.3032 0.9000 0.1250 -0.0000 0.1000 1.0000 0.1250 0.0000 0.1000 Integration Statistics: Number of time steps (nearest 50) = 150 Number of function evaluations (nearest 50) = 400 Number of Jacobian evaluations (nearest 1) = 1 Number of iterations (nearest 1) = 2
- 4~24行目には例1の計算結果が出力されています。
- 6行目にはメッシュ点の数、絶対誤差許容値、相対誤差許容値が出力されています。
- 8行目には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行目には反復数が出力されています。
ソースコード
※本サンプルソースコードはnAG数値計算ライブラリ(Windows, Linux, MAC等に対応)の関数を呼び出します。
サンプルのコンパイル及び実行方法
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642
このソースコードをダウンロード |
/* nag_pde_parab_1d_cd_ode (d03plc) Example Program. * * CLL6I261D/CLL6I261DL Version. * * Copyright 2017 Numerical Algorithms Group. * * Mark 26.1, 2017. */ #include <stdio.h> #include <nag.h> #include <nag_stdlib.h> #include <nagd03.h> #include <nagx01.h> #include <math.h> /* Structure to communicate with user-supplied function arguments */ struct user { double elo, ero, gamma, rlo, rro; }; #ifdef __cplusplus extern "C" { #endif static void nAG_CALL pdedef(Integer, double, double, const double[], const double[], Integer, const double[], const double[], double[], double[], double[], double[], Integer *, Nag_Comm *); static void nAG_CALL bndry1(Integer, Integer, double, const double[], const double[], Integer, const double[], const double[], Integer, double[], Integer *, Nag_Comm *); static void nAG_CALL bndry2(Integer, Integer, double, const double[], const double[], Integer, const double[], const double[], Integer, double[], Integer *, Nag_Comm *); static void nAG_CALL nmflx1(Integer, double, double, Integer, const double[], const double[], const double[], double[], Integer *, Nag_Comm *, Nag_D03_Save *); static void nAG_CALL nmflx2(Integer, double, double, Integer, const double[], const double[], const double[], double[], Integer *, Nag_Comm *, Nag_D03_Save *); static void nAG_CALL odedef(Integer, double, Integer, const double[], const double[], Integer, const double[], const double[], const double[], const double[], double[], Integer *, Nag_Comm *); #ifdef __cplusplus } #endif static void init1(double, double *, Integer, double *, Integer, Integer); static void init2(Integer, Integer, double *, double *, Nag_Comm *); static void exact(double, double *, Integer, const double *, Integer); static int ex1(void); static int ex2(void); #define P(I, J) p[npde*((J) -1)+(I) -1] #define UCP(I, J) ucp[npde*((J) -1)+(I) -1] #define U(I, J) u[npde*((J) -1)+(I) -1] #define UE(I, J) ue[npde*((J) -1)+(I) -1] int main(void) { Integer exit_status_ex1 = 0; Integer exit_status_ex2 = 0; printf("nag_pde_parab_1d_cd_ode (d03plc) Example Program Results\n"); exit_status_ex1 = ex1(); exit_status_ex2 = ex2(); return (exit_status_ex1 == 0 && exit_status_ex2 == 0) ? 0 : 1; } int ex1(void) { /* Constants */ const Integer npde = 2, npts = 201, nv = 2, nxi = 2; const Integer neqn = npde*npts + nv; static double ruser1[4] = { -1.0, -1.0, -1.0, -1.0 }; /* Scalars */ double tout, ts, errmax, lerr, lwgt; Integer exit_status = 0, i, ind, itask, itol, itrace, j; Integer nsteps, nfuncs, njacs, niters, ierrmax; Integer nwkres, lenode, lisave, lrsave; /* Arrays */ double *algopt = 0, *atol = 0, *rsave = 0, *rtol = 0, *u = 0, *ue = 0; double *x = 0, *xi = 0; Integer *isave = 0; /* Nag Types */ NagError fail; Nag_Comm comm; Nag_D03_Save saved; INIT_FAIL(fail); /* For communication with user-supplied functions: */ comm.user = ruser1; printf("\n\nExample 1\n\n"); /* Allocate memory */ lisave = 25*neqn + 24; nwkres = npde*(2*npts+6*nxi+3*npde+26) + nxi + nv + 7*npts + 2; lenode = 11*neqn + 50; lrsave = 4*neqn + 11*neqn/2 + 1 + nwkres + lenode; lisave = lisave*4; lrsave = lrsave*4; if (!(algopt = nAG_ALLOC(30, double)) || !(atol = nAG_ALLOC(1, double)) || !(rsave = nAG_ALLOC(lrsave, double)) || !(rtol = nAG_ALLOC(1, double)) || !(u = nAG_ALLOC(neqn, double)) || !(ue = nAG_ALLOC(npde * npts, double)) || !(x = nAG_ALLOC(npts, double)) || !(xi = nAG_ALLOC(nxi, double)) || !(isave = nAG_ALLOC(lisave, Integer))) { printf("Allocation failure\n"); exit_status = 1; goto END; } itrace = 0; itol = 1; atol[0] = 1e-5; rtol[0] = 2.5e-4; printf(" Method parameters:\n"); printf(" Number of mesh points used = %4ld\n", npts); printf(" Relative tolerance used = %12.3e\n", rtol[0]); printf(" Absolute tolerance used = %12.3e\n\n", atol[0]); /* Initialize mesh */ for (i = 0; i < npts; ++i) x[i] = i / (npts - 1.0); xi[0] = 0.0; xi[1] = 1.0; /* Set initial values */ ts = 0.0; init1(ts, u, npde, x, npts, nv); ind = 0; itask = 1; for (i = 0; i < 30; ++i) algopt[i] = 0.0; /* BDF integration */ algopt[0] = 1.0; /* Sparse matrix algebra parameters */ algopt[28] = 0.1; algopt[29] = 1.1; tout = 0.5; /* nag_pde_parab_1d_cd_ode (d03plc). * General system of convection-diffusion PDEs with source * terms in conservative form, coupled DAEs, method of * lines, upwind scheme using numerical flux function based * on Riemann solver, one space variable */ nag_pde_parab_1d_cd_ode(npde, &ts, tout, pdedef, nmflx1, bndry1, u, npts, x, nv, odedef, nxi, xi, neqn, rtol, atol, itol, Nag_OneNorm, Nag_LinAlgSparse, algopt, rsave, lrsave, isave, lisave, itask, itrace, 0, &ind, &comm, &saved, &fail); if (fail.code != NE_NOERROR) { printf("Error from nag_pde_parab_1d_cd_ode (d03plc).\n%s\n", fail.message); exit_status = 1; goto END; } /* Check against exact solution */ exact(tout, ue, npde, &x[0], npts); errmax = 0.0; for (i=1; i<npts; i++) { lerr = 0.0; for (j=0; j<npde; j++) { lwgt = rtol[0]*fabs(ue[i*npde+j]) + atol[0]; lerr += fabs(u[i*npde+j]-ue[i*npde+j])/lwgt; } lerr = lerr/(double) npde; errmax = MAX(errmax,lerr); } ierrmax = 100*(Integer)(errmax/100.0) + 100; printf("\n Integration Results:\n"); printf(" Global error is less than %4ld" " times the local error tolerance.\n", ierrmax); /* Print integration statistics (reasonably rounded) */ nsteps = 50*((isave[0]+25)/50); nfuncs = 100*((isave[1]+50)/100); njacs = 20*((isave[2]+10)/20); niters = 100*((isave[4]+50)/100); printf("\n Integration Statistics:\n"); printf(" %-30s (nearest %3d) = %6ld\n", "Number of time steps", 50, nsteps); printf(" %-30s (nearest %3d) = %6ld\n", "Number of function evaluations", 100, nfuncs); printf(" %-30s (nearest %3d) = %6ld\n", "Number of Jacobian evaluations", 20, njacs); printf(" %-30s (nearest %3d) = %6ld\n", "Number of iterations", 100, niters); END: nAG_FREE(algopt); nAG_FREE(atol); nAG_FREE(rsave); nAG_FREE(rtol); nAG_FREE(u); nAG_FREE(ue); nAG_FREE(x); nAG_FREE(xi); nAG_FREE(isave); return exit_status; } static void nAG_CALL pdedef(Integer npde, double t, double x, const double u[], const double ux[], Integer nv, const double v[], const double vdot[], double p[], double c[], double d[], double s[], Integer *ires, Nag_Comm *comm) { Integer i, j; if (comm->user[2] == -1.0) { /* printf("(User-supplied callback pdedef, first invocation.)\n"); */ comm->user[2] = 0.0; } for (i = 1; i <= npde; ++i) { c[i - 1] = 1.0; d[i - 1] = 0.0; s[i - 1] = 0.0; for (j = 1; j <= npde; ++j) { if (i == j) { P(i, j) = 1.0; } else { P(i, j) = 0.0; } } } return; } static void nAG_CALL bndry1(Integer npde, Integer npts, double t, const double x[], const double u[], Integer nv, const double v[], const double vdot[], Integer ibnd, double g[], Integer *ires, Nag_Comm *comm) { double dudx; double *ue = 0; if (comm->user[0] == -1.0) { /* printf("(User-supplied callback bndry1, first invocation.)\n"); */ comm->user[0] = 0.0; } /* Allocate memory */ if (!(ue = nAG_ALLOC(npde, double))) { printf("Allocation failure\n"); goto END; } if (ibnd == 0) { exact(t, ue, npde, &x[0], 1); g[0] = U(1, 1) + U(2, 1) - UE(1, 1) - UE(2, 1); dudx = (U(1, 2) - U(2, 2) - U(1, 1) + U(2, 1)) / (x[1] - x[0]); g[1] = vdot[0] - dudx; } else { exact(t, ue, npde, &x[npts - 1], 1); g[0] = U(1, npts) - U(2, npts) - UE(1, 1) + UE(2, 1); dudx = (U(1, npts) + U(2, npts) - U(1, npts - 1) - U(2, npts - 1)) / (x[npts - 1] - x[npts - 2]); g[1] = vdot[1] + 3.0 * dudx; } END: nAG_FREE(ue); return; } static void nAG_CALL nmflx1(Integer npde, double t, double x, Integer nv, const double v[], const double uleft[], const double uright[], double flux[], Integer *ires, Nag_Comm *comm, Nag_D03_Save *saved) { if (comm->user[1] == -1.0) { /* printf("(User-supplied callback nmflx1, first invocation.)\n"); */ comm->user[1] = 0.0; } flux[0] = 0.5 * (3.0 * uleft[0] - uright[0] + 3.0 * uleft[1] + uright[1]); flux[1] = 0.5 * (3.0 * uleft[0] + uright[0] + 3.0 * uleft[1] - uright[1]); return; } static void nAG_CALL odedef(Integer npde, double t, Integer nv, const double v[], const double vdot[], Integer nxi, const double xi[], const double ucp[], const double ucpx[], const double ucpt[], double r[], Integer *ires, Nag_Comm *comm) { if (comm->user[3] == -1.0) { /* printf("(User-supplied callback odedef, first invocation.)\n"); */ comm->user[3] = 0.0; } if (*ires == -1) { r[0] = 0.0; r[1] = 0.0; } else { r[0] = v[0] - UCP(1, 1) + UCP(2, 1); r[1] = v[1] - UCP(1, 2) - UCP(2, 2); } return; } static void exact(double t, double *u, Integer npde, const double *x, Integer npts) { /* Exact solution (for comparison and b.c. purposes) */ double f, g, x1, x2; Integer i; for (i = 0; i < npts; ++i) { x1 = x[i] - 3.0 * t; x2 = x[i] + t; f = exp(nag_pi * x1) * sin(2.0 * nag_pi * x1); g = exp(-2.0 * nag_pi * x2) * cos(2.0 * nag_pi * x2); u[npde*i] = f + g; u[npde*i+1] = f - g; } return; } static void init1(double t, double *u, Integer npde, double *x, Integer npts, Integer nv) { /* Initial solution */ double f, g, x1, x2; Integer i, neqn; neqn = npde * npts + nv; for (i = 0; i < npts; ++i) { x1 = x[i] - 3.0 * t; x2 = x[i] + t; f = exp(nag_pi * x1) * sin(2.0 * nag_pi * x1); g = exp(-2.0 * nag_pi * x2) * cos(2.0 * nag_pi * x2); u[npde*i] = f + g; u[npde*i+1] = f - g; } u[neqn - 2] = u[0] - u[1]; u[neqn - 1] = u[neqn - 3] + u[neqn - 4]; return; } int ex2(void) { const Integer npde = 3, npts = 141, nv = 0, nxi = 0; const Integer neqn = npde * npts + nv, lisave = neqn + 24; const Integer lrsave = 16392; static double ruser[2] = { -1.0, -1.0 }; double d, p, tout, ts, v; Integer exit_status = 0, i, ind, it, itask, itol, itrace; Integer nsteps, nfuncs, njacs, niters; double *algopt = 0, *atol = 0, *rsave = 0, *rtol = 0, *u = 0; double *x = 0, *xi = 0; Integer *isave = 0; NagError fail; Nag_Comm comm; Nag_D03_Save saved; struct user data; INIT_FAIL(fail); /* For communication with user-supplied functions: */ comm.user = ruser; printf("\n\nExample 2\n\n"); /* Allocate memory */ if (!(algopt = nAG_ALLOC(30, double)) || !(atol = nAG_ALLOC(1, double)) || !(rsave = nAG_ALLOC(lrsave, double)) || !(rtol = nAG_ALLOC(1, double)) || !(u = nAG_ALLOC(npde * npts, double)) || !(x = nAG_ALLOC(npts, double)) || !(xi = nAG_ALLOC(1, double)) || !(isave = nAG_ALLOC(447, Integer))) { printf("Allocation failure\n"); exit_status = -1; goto END; } /* Problem parameters */ data.elo = 2.5; data.ero = 0.25; data.gamma = 1.4; data.rlo = 1.0; data.rro = 0.125; comm.p = (Pointer) &data; itrace = 0; itol = 1; atol[0] = 0.005; rtol[0] = 5e-4; printf(" Problem parameters and initial conditions:\n"); printf(" gamma = %5.3f\n", data.gamma); printf(" e(x<0.5,0) = %5.3f", data.elo); printf(" e(x>0.5,0) = %5.3f\n", data.ero); printf(" rho(x<0.5,0) = %5.3f", data.rlo); printf(" rho(x>0.5,0) = %5.3f\n\n", data.rro); printf(" Method parameters:\n"); printf(" Number of mesh points used = %4ld\n", npts); printf(" Relative tolerance used = %12.3e\n", rtol[0]); printf(" Absolute tolerance used = %12.3e\n\n", atol[0]); /* Initialize mesh */ for (i = 0; i < npts; ++i) x[i] = i / (npts - 1.0); /* Initial values of variables */ init2(npde, npts, x, u, &comm); xi[0] = 0.0; ind = 0; itask = 1; for (i = 0; i < 30; ++i) algopt[i] = 0.0; /* Theta integration */ algopt[0] = 2.0; algopt[5] = 2.0; algopt[6] = 2.0; /* Max. time step */ algopt[12] = 0.005; ts = 0.0; printf(" Solution\n%4s%9s%9s%9s%9s\n", "t", "x", "d", "v", "p"); for (it = 0; it < 2; ++it) { tout = 0.1 * (it + 1); /* nag_pde_parab_1d_cd_ode (d03plc), see above. */ nag_pde_parab_1d_cd_ode(npde, &ts, tout, NULLFN, nmflx2, bndry2, u, npts, x, nv, NULLFN, nxi, xi, neqn, rtol, atol, itol, Nag_TwoNorm, Nag_LinAlgBand, algopt, rsave, lrsave, isave, lisave, itask, itrace, 0, &ind, &comm, &saved, &fail); if (fail.code != NE_NOERROR) { printf("Error from nag_pde_parab_1d_cd_ode (d03plc).\n%s\n", fail.message); exit_status = 1; goto END; } /* Calculate density, velocity and pressure */ for (i = 1; i <= npts; i += 14) { d = U(1, i); v = U(2, i) / d; p = d*(data.gamma - 1.0)*(U(3, i)/d - 0.5*v*v); if (i==1) { printf("%6.3f %7.4f %7.4f %7.4f %7.4f\n", ts, x[i-1], d, v, p); } else { printf("%6s %7.4f %7.4f %7.4f %7.4f\n", "", x[i-1], d, v, p); } } printf("\n"); } /* Print integration statistics (reasonably rounded) */ nsteps = 50*((isave[0]+25)/50); nfuncs = 50*((isave[1]+25)/50); njacs = isave[2]; niters = isave[4]; printf("\n Integration Statistics:\n"); printf(" %-30s (nearest %3ld) = %6ld\n", "Number of time steps", 50, nsteps); printf(" %-30s (nearest %3ld) = %6ld\n", "Number of function evaluations", 50, nfuncs); printf(" %-30s (nearest %3ld) = %6ld\n", "Number of Jacobian evaluations", 1, njacs); printf(" %-30s (nearest %3ld) = %6ld\n", "Number of iterations", 1, niters); END: nAG_FREE(algopt); nAG_FREE(atol); nAG_FREE(rsave); nAG_FREE(rtol); nAG_FREE(u); nAG_FREE(x); nAG_FREE(xi); nAG_FREE(isave); return exit_status; } static void init2(Integer npde, Integer npts, double *x, double *u, Nag_Comm *comm) { Integer i, j; struct user *data = (struct user *) comm->p; j = 0; for (i = 0; i < npts; ++i) { if (x[i] < 0.5) { u[j] = data->rlo; u[j + 1] = 0.0; u[j + 2] = data->elo; } else if (x[i] == 0.5) { u[j] = 0.5 * (data->rlo + data->rro); u[j + 1] = 0.0; u[j + 2] = 0.5 * (data->elo + data->ero); } else { u[j] = data->rro; u[j + 1] = 0.0; u[j + 2] = data->ero; } j += 3; } return; } static void nAG_CALL bndry2(Integer npde, Integer npts, double t, const double x[], const double u[], Integer nv, const double v[], const double vdot[], Integer ibnd, double g[], Integer *ires, Nag_Comm *comm) { struct user *data = (struct user *) comm->p; if (comm->user[0] == -1.0) { /* printf("(User-supplied callback bndry2, first invocation.)\n"); */ comm->user[0] = 0.0; } if (ibnd == 0) { g[0] = U(1, 1) - data->rlo; g[1] = U(2, 1); g[2] = U(3, 1) - data->elo; } else { g[0] = U(1, npts) - data->rro; g[1] = U(2, npts); g[2] = U(3, npts) - data->ero; } return; } static void nAG_CALL nmflx2(Integer npde, double t, double x, Integer nv, const double v[], const double uleft[], const double uright[], double flux[], Integer *ires, Nag_Comm *comm, Nag_D03_Save *saved) { char solver; NagError fail; struct user *data = (struct user *) comm->p; if (comm->user[1] == -1.0) { /* printf("(User-supplied callback nmflx2, first invocation.)\n"); */ comm->user[1] = 0.0; } INIT_FAIL(fail); solver = 'R'; if (solver == 'R') { /* ROE SCHEME */ /* nag_pde_parab_1d_euler_roe (d03puc). * Roe's approximate Riemann solver for Euler equations in * conservative form, for use with nag_pde_parab_1d_cd * (d03pfc), nag_pde_parab_1d_cd_ode (d03plc) and * nag_pde_parab_1d_cd_ode_remesh (d03psc) */ nag_pde_parab_1d_euler_roe(uleft, uright, data->gamma, flux, saved, &fail); } else { /* OSHER SCHEME */ /* nag_pde_parab_1d_euler_osher (d03pvc). * Osher's approximate Riemann solver for Euler equations in * conservative form, for use with nag_pde_parab_1d_cd * (d03pfc), nag_pde_parab_1d_cd_ode (d03plc) and * nag_pde_parab_1d_cd_ode_remesh (d03psc) */ nag_pde_parab_1d_euler_osher(uleft, uright, data->gamma, Nag_OsherPhysical, flux, saved, &fail); } if (fail.code != NE_NOERROR) { printf("Error from nag_pde_parab_1d_euler_osher (d03pvc).\n%s\n", fail.message); } return; }