/usr/share/doc/libsundials-serial-dev/examples/cvodes/parallel/cvsDiurnal_kry_p.c is in libsundials-serial-dev 2.5.0-3+b3.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
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 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 | /*
* -----------------------------------------------------------------
* $Revision: 1.4 $
* $Date: 2010/12/01 23:00:48 $
* -----------------------------------------------------------------
* Programmer(s): S. D. Cohen, A. C. Hindmarsh, M. R. Wittman, and
* Radu Serban @ LLNL
* -----------------------------------------------------------------
* Example problem:
*
* An ODE system is generated from the following 2-species diurnal
* kinetics advection-diffusion PDE system in 2 space dimensions:
*
* dc(i)/dt = Kh*(d/dx)^2 c(i) + V*dc(i)/dx + (d/dy)(Kv(y)*dc(i)/dy)
* + Ri(c1,c2,t) for i = 1,2, where
* R1(c1,c2,t) = -q1*c1*c3 - q2*c1*c2 + 2*q3(t)*c3 + q4(t)*c2 ,
* R2(c1,c2,t) = q1*c1*c3 - q2*c1*c2 - q4(t)*c2 ,
* Kv(y) = Kv0*exp(y/5) ,
* Kh, V, Kv0, q1, q2, and c3 are constants, and q3(t) and q4(t)
* vary diurnally. The problem is posed on the square
* 0 <= x <= 20, 30 <= y <= 50 (all in km),
* with homogeneous Neumann boundary conditions, and for time t in
* 0 <= t <= 86400 sec (1 day).
* The PDE system is treated by central differences on a uniform
* mesh, with simple polynomial initial profiles.
*
* The problem is solved by CVODES on NPE processors, treated
* as a rectangular process grid of size NPEX by NPEY, with
* NPE = NPEX*NPEY. Each processor contains a subgrid of size MXSUB
* by MYSUB of the (x,y) mesh. Thus the actual mesh sizes are
* MX = MXSUB*NPEX and MY = MYSUB*NPEY, and the ODE system size is
* neq = 2*MX*MY.
*
* The solution is done with the BDF/GMRES method (i.e. using the
* CVSPGMR linear solver) and the block-diagonal part of the
* Newton matrix as a left preconditioner. A copy of the
* block-diagonal part of the Jacobian is saved and conditionally
* reused within the preconditioner routine.
*
* Performance data and sampled solution values are printed at
* selected output times, and all performance counters are printed
* on completion.
*
* This version uses MPI for user routines.
*
* Execution: mpirun -np N cvsDiurnal_kry_p with N = NPEX*NPEY
* (see constants below).
* -----------------------------------------------------------------
*/
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <cvodes/cvodes.h> /* prototypes for CVODE fcts. */
#include <cvodes/cvodes_spgmr.h> /* prototypes and constants for CVSPGMR solver */
#include <nvector/nvector_parallel.h> /* definition N_Vector and macro NV_DATA_P */
#include <sundials/sundials_dense.h> /* prototypes for small dense matrix fcts. */
#include <sundials/sundials_types.h> /* definitions of realtype, booleantype */
#include <sundials/sundials_math.h> /* definition of macros SQR and EXP */
#include <mpi.h> /* MPI constants and types */
/* Problem Constants */
#define NVARS 2 /* number of species */
#define KH RCONST(4.0e-6) /* horizontal diffusivity Kh */
#define VEL RCONST(0.001) /* advection velocity V */
#define KV0 RCONST(1.0e-8) /* coefficient in Kv(y) */
#define Q1 RCONST(1.63e-16) /* coefficients q1, q2, c3 */
#define Q2 RCONST(4.66e-16)
#define C3 RCONST(3.7e16)
#define A3 RCONST(22.62) /* coefficient in expression for q3(t) */
#define A4 RCONST(7.601) /* coefficient in expression for q4(t) */
#define C1_SCALE RCONST(1.0e6) /* coefficients in initial profiles */
#define C2_SCALE RCONST(1.0e12)
#define T0 RCONST(0.0) /* initial time */
#define NOUT 12 /* number of output times */
#define TWOHR RCONST(7200.0) /* number of seconds in two hours */
#define HALFDAY RCONST(4.32e4) /* number of seconds in a half day */
#define PI RCONST(3.1415926535898) /* pi */
#define XMIN RCONST(0.0) /* grid boundaries in x */
#define XMAX RCONST(20.0)
#define YMIN RCONST(30.0) /* grid boundaries in y */
#define YMAX RCONST(50.0)
#define NPEX 2 /* no. PEs in x direction of PE array */
#define NPEY 2 /* no. PEs in y direction of PE array */
/* Total no. PEs = NPEX*NPEY */
#define MXSUB 5 /* no. x points per subgrid */
#define MYSUB 5 /* no. y points per subgrid */
#define MX (NPEX*MXSUB) /* MX = number of x mesh points */
#define MY (NPEY*MYSUB) /* MY = number of y mesh points */
/* Spatial mesh is MX by MY */
/* CVodeMalloc Constants */
#define RTOL RCONST(1.0e-5) /* scalar relative tolerance */
#define FLOOR RCONST(100.0) /* value of C1 or C2 at which tolerances */
/* change from relative to absolute */
#define ATOL (RTOL*FLOOR) /* scalar absolute tolerance */
/* User-defined matrix accessor macro: IJth */
/* IJth is defined in order to write code which indexes into dense
matrices with a (row,column) pair, where 1 <= row,column <= NVARS.
IJth(a,i,j) references the (i,j)th entry of the small matrix realtype **a,
where 1 <= i,j <= NVARS. The small matrix routines in sundials_dense.h
work with matrices stored by column in a 2-dimensional array. In C,
arrays are indexed starting at 0, not 1. */
#define IJth(a,i,j) (a[j-1][i-1])
/* Type : UserData
contains problem constants, preconditioner blocks, pivot arrays,
grid constants, and processor indices, as well as data needed
for the preconditiner */
typedef struct {
realtype q4, om, dx, dy, hdco, haco, vdco;
realtype uext[NVARS*(MXSUB+2)*(MYSUB+2)];
int my_pe, isubx, isuby;
long int nvmxsub, nvmxsub2;
MPI_Comm comm;
/* For preconditioner */
realtype **P[MXSUB][MYSUB], **Jbd[MXSUB][MYSUB];
long int *pivot[MXSUB][MYSUB];
} *UserData;
/* Private Helper Functions */
static void InitUserData(int my_pe, MPI_Comm comm, UserData data);
static void FreeUserData(UserData data);
static void SetInitialProfiles(N_Vector u, UserData data);
static void PrintOutput(void *cvode_mem, int my_pe, MPI_Comm comm,
N_Vector u, realtype t);
static void PrintFinalStats(void *cvode_mem);
static void BSend(MPI_Comm comm,
int my_pe, int isubx, int isuby,
long int dsizex, long int dsizey,
realtype udata[]);
static void BRecvPost(MPI_Comm comm, MPI_Request request[],
int my_pe, int isubx, int isuby,
long int dsizex, long int dsizey,
realtype uext[], realtype buffer[]);
static void BRecvWait(MPI_Request request[],
int isubx, int isuby,
long int dsizex, realtype uext[],
realtype buffer[]);
static void ucomm(realtype t, N_Vector u, UserData data);
static void fcalc(realtype t, realtype udata[], realtype dudata[],
UserData data);
/* Functions Called by the Solver */
static int f(realtype t, N_Vector u, N_Vector udot, void *user_data);
static int Precond(realtype tn, N_Vector u, N_Vector fu,
booleantype jok, booleantype *jcurPtr,
realtype gamma, void *user_data,
N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3);
static int PSolve(realtype tn, N_Vector u, N_Vector fu,
N_Vector r, N_Vector z,
realtype gamma, realtype delta,
int lr, void *user_data, N_Vector vtemp);
/* Private function to check function return values */
static int check_flag(void *flagvalue, char *funcname, int opt, int id);
/***************************** Main Program ******************************/
int main(int argc, char *argv[])
{
realtype abstol, reltol, t, tout;
N_Vector u;
UserData data;
void *cvode_mem;
int iout, flag, my_pe, npes;
long int neq, local_N;
MPI_Comm comm;
u = NULL;
data = NULL;
cvode_mem = NULL;
/* Set problem size neq */
neq = NVARS*MX*MY;
/* Get processor number and total number of pe's */
MPI_Init(&argc, &argv);
comm = MPI_COMM_WORLD;
MPI_Comm_size(comm, &npes);
MPI_Comm_rank(comm, &my_pe);
if (npes != NPEX*NPEY) {
if (my_pe == 0)
fprintf(stderr, "\nMPI_ERROR(0): npes = %d is not equal to NPEX*NPEY = %d\n\n",
npes,NPEX*NPEY);
MPI_Finalize();
return(1);
}
/* Set local length */
local_N = NVARS*MXSUB*MYSUB;
/* Allocate and load user data block; allocate preconditioner block */
data = (UserData) malloc(sizeof *data);
if (check_flag((void *)data, "malloc", 2, my_pe)) MPI_Abort(comm, 1);
InitUserData(my_pe, comm, data);
/* Allocate u, and set initial values and tolerances */
u = N_VNew_Parallel(comm, local_N, neq);
if (check_flag((void *)u, "N_VNew", 0, my_pe)) MPI_Abort(comm, 1);
SetInitialProfiles(u, data);
abstol = ATOL; reltol = RTOL;
/* Call CVodeCreate to create the solver memory and specify the
* Backward Differentiation Formula and the use of a Newton iteration */
cvode_mem = CVodeCreate(CV_BDF, CV_NEWTON);
if (check_flag((void *)cvode_mem, "CVodeCreate", 0, my_pe)) MPI_Abort(comm, 1);
/* Set the pointer to user-defined data */
flag = CVodeSetUserData(cvode_mem, data);
if (check_flag(&flag, "CVodeSetUserData", 1, my_pe)) MPI_Abort(comm, 1);
/* Call CVodeInit to initialize the integrator memory and specify the
* user's right hand side function in u'=f(t,u), the inital time T0, and
* the initial dependent variable vector u. */
flag = CVodeInit(cvode_mem, f, T0, u);
if(check_flag(&flag, "CVodeInit", 1, my_pe)) return(1);
/* Call CVodeSStolerances to specify the scalar relative tolerance
* and scalar absolute tolerances */
flag = CVodeSStolerances(cvode_mem, reltol, abstol);
if (check_flag(&flag, "CVodeSStolerances", 1, my_pe)) return(1);
/* Call CVSpgmr to specify the linear solver CVSPGMR
with left preconditioning and the maximum Krylov dimension maxl */
flag = CVSpgmr(cvode_mem, PREC_LEFT, 0);
if (check_flag(&flag, "CVSpgmr", 1, my_pe)) MPI_Abort(comm, 1);
/* Set preconditioner setup and solve routines Precond and PSolve,
and the pointer to the user-defined block data */
flag = CVSpilsSetPreconditioner(cvode_mem, Precond, PSolve);
if (check_flag(&flag, "CVSpilsSetPreconditioner", 1, my_pe)) MPI_Abort(comm, 1);
if (my_pe == 0)
printf("\n2-species diurnal advection-diffusion problem\n\n");
/* In loop over output points, call CVode, print results, test for error */
for (iout=1, tout = TWOHR; iout <= NOUT; iout++, tout += TWOHR) {
flag = CVode(cvode_mem, tout, u, &t, CV_NORMAL);
if (check_flag(&flag, "CVode", 1, my_pe)) break;
PrintOutput(cvode_mem, my_pe, comm, u, t);
}
/* Print final statistics */
if (my_pe == 0) PrintFinalStats(cvode_mem);
/* Free memory */
N_VDestroy_Parallel(u);
FreeUserData(data);
CVodeFree(&cvode_mem);
MPI_Finalize();
return(0);
}
/*********************** Private Helper Functions ************************/
/* Load constants in data */
static void InitUserData(int my_pe, MPI_Comm comm, UserData data)
{
int isubx, isuby;
int lx, ly;
/* Set problem constants */
data->om = PI/HALFDAY;
data->dx = (XMAX-XMIN)/((realtype)(MX-1));
data->dy = (YMAX-YMIN)/((realtype)(MY-1));
data->hdco = KH/SQR(data->dx);
data->haco = VEL/(RCONST(2.0)*data->dx);
data->vdco = (RCONST(1.0)/SQR(data->dy))*KV0;
/* Set machine-related constants */
data->comm = comm;
data->my_pe = my_pe;
/* isubx and isuby are the PE grid indices corresponding to my_pe */
isuby = my_pe/NPEX;
isubx = my_pe - isuby*NPEX;
data->isubx = isubx;
data->isuby = isuby;
/* Set the sizes of a boundary x-line in u and uext */
data->nvmxsub = NVARS*MXSUB;
data->nvmxsub2 = NVARS*(MXSUB+2);
/* Preconditioner-related fields */
for (lx = 0; lx < MXSUB; lx++) {
for (ly = 0; ly < MYSUB; ly++) {
(data->P)[lx][ly] = newDenseMat(NVARS, NVARS);
(data->Jbd)[lx][ly] = newDenseMat(NVARS, NVARS);
(data->pivot)[lx][ly] = newLintArray(NVARS);
}
}
}
/* Free user data memory */
static void FreeUserData(UserData data)
{
int lx, ly;
for (lx = 0; lx < MXSUB; lx++) {
for (ly = 0; ly < MYSUB; ly++) {
destroyMat((data->P)[lx][ly]);
destroyMat((data->Jbd)[lx][ly]);
destroyArray((data->pivot)[lx][ly]);
}
}
free(data);
}
/* Set initial conditions in u */
static void SetInitialProfiles(N_Vector u, UserData data)
{
int isubx, isuby, lx, ly, jx, jy;
long int offset;
realtype dx, dy, x, y, cx, cy, xmid, ymid;
realtype *udata;
/* Set pointer to data array in vector u */
udata = NV_DATA_P(u);
/* Get mesh spacings, and subgrid indices for this PE */
dx = data->dx; dy = data->dy;
isubx = data->isubx; isuby = data->isuby;
/* Load initial profiles of c1 and c2 into local u vector.
Here lx and ly are local mesh point indices on the local subgrid,
and jx and jy are the global mesh point indices. */
offset = 0;
xmid = RCONST(0.5)*(XMIN + XMAX);
ymid = RCONST(0.5)*(YMIN + YMAX);
for (ly = 0; ly < MYSUB; ly++) {
jy = ly + isuby*MYSUB;
y = YMIN + jy*dy;
cy = SQR(RCONST(0.1)*(y - ymid));
cy = RCONST(1.0) - cy + RCONST(0.5)*SQR(cy);
for (lx = 0; lx < MXSUB; lx++) {
jx = lx + isubx*MXSUB;
x = XMIN + jx*dx;
cx = SQR(RCONST(0.1)*(x - xmid));
cx = RCONST(1.0) - cx + RCONST(0.5)*SQR(cx);
udata[offset ] = C1_SCALE*cx*cy;
udata[offset+1] = C2_SCALE*cx*cy;
offset = offset + 2;
}
}
}
/* Print current t, step count, order, stepsize, and sampled c1,c2 values */
static void PrintOutput(void *cvode_mem, int my_pe, MPI_Comm comm,
N_Vector u, realtype t)
{
int qu, flag;
realtype hu, *udata, tempu[2];
int npelast;
long int i0, i1, nst;
MPI_Status status;
npelast = NPEX*NPEY - 1;
udata = NV_DATA_P(u);
/* Send c1,c2 at top right mesh point to PE 0 */
if (my_pe == npelast) {
i0 = NVARS*MXSUB*MYSUB - 2;
i1 = i0 + 1;
if (npelast != 0)
MPI_Send(&udata[i0], 2, PVEC_REAL_MPI_TYPE, 0, 0, comm);
else {
tempu[0] = udata[i0];
tempu[1] = udata[i1];
}
}
/* On PE 0, receive c1,c2 at top right, then print performance data
and sampled solution values */
if (my_pe == 0) {
if (npelast != 0)
MPI_Recv(&tempu[0], 2, PVEC_REAL_MPI_TYPE, npelast, 0, comm, &status);
flag = CVodeGetNumSteps(cvode_mem, &nst);
check_flag(&flag, "CVodeGetNumSteps", 1, my_pe);
flag = CVodeGetLastOrder(cvode_mem, &qu);
check_flag(&flag, "CVodeGetLastOrder", 1, my_pe);
flag = CVodeGetLastStep(cvode_mem, &hu);
check_flag(&flag, "CVodeGetLastStep", 1, my_pe);
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf("t = %.2Le no. steps = %ld order = %d stepsize = %.2Le\n",
t, nst, qu, hu);
printf("At bottom left: c1, c2 = %12.3Le %12.3Le \n", udata[0], udata[1]);
printf("At top right: c1, c2 = %12.3Le %12.3Le \n\n", tempu[0], tempu[1]);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
printf("t = %.2le no. steps = %ld order = %d stepsize = %.2le\n",
t, nst, qu, hu);
printf("At bottom left: c1, c2 = %12.3le %12.3le \n", udata[0], udata[1]);
printf("At top right: c1, c2 = %12.3le %12.3le \n\n", tempu[0], tempu[1]);
#else
printf("t = %.2e no. steps = %ld order = %d stepsize = %.2e\n",
t, nst, qu, hu);
printf("At bottom left: c1, c2 = %12.3e %12.3e \n", udata[0], udata[1]);
printf("At top right: c1, c2 = %12.3e %12.3e \n\n", tempu[0], tempu[1]);
#endif
}
}
/* Print final statistics contained in iopt */
static void PrintFinalStats(void *cvode_mem)
{
long int lenrw, leniw ;
long int lenrwLS, leniwLS;
long int nst, nfe, nsetups, nni, ncfn, netf;
long int nli, npe, nps, ncfl, nfeLS;
int flag;
flag = CVodeGetWorkSpace(cvode_mem, &lenrw, &leniw);
check_flag(&flag, "CVodeGetWorkSpace", 1, 0);
flag = CVodeGetNumSteps(cvode_mem, &nst);
check_flag(&flag, "CVodeGetNumSteps", 1, 0);
flag = CVodeGetNumRhsEvals(cvode_mem, &nfe);
check_flag(&flag, "CVodeGetNumRhsEvals", 1, 0);
flag = CVodeGetNumLinSolvSetups(cvode_mem, &nsetups);
check_flag(&flag, "CVodeGetNumLinSolvSetups", 1, 0);
flag = CVodeGetNumErrTestFails(cvode_mem, &netf);
check_flag(&flag, "CVodeGetNumErrTestFails", 1, 0);
flag = CVodeGetNumNonlinSolvIters(cvode_mem, &nni);
check_flag(&flag, "CVodeGetNumNonlinSolvIters", 1, 0);
flag = CVodeGetNumNonlinSolvConvFails(cvode_mem, &ncfn);
check_flag(&flag, "CVodeGetNumNonlinSolvConvFails", 1, 0);
flag = CVSpilsGetWorkSpace(cvode_mem, &lenrwLS, &leniwLS);
check_flag(&flag, "CVSpilsGetWorkSpace", 1, 0);
flag = CVSpilsGetNumLinIters(cvode_mem, &nli);
check_flag(&flag, "CVSpilsGetNumLinIters", 1, 0);
flag = CVSpilsGetNumPrecEvals(cvode_mem, &npe);
check_flag(&flag, "CVSpilsGetNumPrecEvals", 1, 0);
flag = CVSpilsGetNumPrecSolves(cvode_mem, &nps);
check_flag(&flag, "CVSpilsGetNumPrecSolves", 1, 0);
flag = CVSpilsGetNumConvFails(cvode_mem, &ncfl);
check_flag(&flag, "CVSpilsGetNumConvFails", 1, 0);
flag = CVSpilsGetNumRhsEvals(cvode_mem, &nfeLS);
check_flag(&flag, "CVSpilsGetNumRhsEvals", 1, 0);
printf("\nFinal Statistics: \n\n");
printf("lenrw = %5ld leniw = %5ld\n", lenrw, leniw);
printf("lenrwls = %5ld leniwls = %5ld\n", lenrwLS, leniwLS);
printf("nst = %5ld\n" , nst);
printf("nfe = %5ld nfels = %5ld\n" , nfe, nfeLS);
printf("nni = %5ld nli = %5ld\n" , nni, nli);
printf("nsetups = %5ld netf = %5ld\n" , nsetups, netf);
printf("npe = %5ld nps = %5ld\n" , npe, nps);
printf("ncfn = %5ld ncfl = %5ld\n\n", ncfn, ncfl);
}
/* Routine to send boundary data to neighboring PEs */
static void BSend(MPI_Comm comm,
int my_pe, int isubx, int isuby,
long int dsizex, long int dsizey,
realtype udata[])
{
int i, ly;
long int offsetu, offsetbuf;
realtype bufleft[NVARS*MYSUB], bufright[NVARS*MYSUB];
/* If isuby > 0, send data from bottom x-line of u */
if (isuby != 0)
MPI_Send(&udata[0], dsizex, PVEC_REAL_MPI_TYPE, my_pe-NPEX, 0, comm);
/* If isuby < NPEY-1, send data from top x-line of u */
if (isuby != NPEY-1) {
offsetu = (MYSUB-1)*dsizex;
MPI_Send(&udata[offsetu], dsizex, PVEC_REAL_MPI_TYPE, my_pe+NPEX, 0, comm);
}
/* If isubx > 0, send data from left y-line of u (via bufleft) */
if (isubx != 0) {
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NVARS;
offsetu = ly*dsizex;
for (i = 0; i < NVARS; i++)
bufleft[offsetbuf+i] = udata[offsetu+i];
}
MPI_Send(&bufleft[0], dsizey, PVEC_REAL_MPI_TYPE, my_pe-1, 0, comm);
}
/* If isubx < NPEX-1, send data from right y-line of u (via bufright) */
if (isubx != NPEX-1) {
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NVARS;
offsetu = offsetbuf*MXSUB + (MXSUB-1)*NVARS;
for (i = 0; i < NVARS; i++)
bufright[offsetbuf+i] = udata[offsetu+i];
}
MPI_Send(&bufright[0], dsizey, PVEC_REAL_MPI_TYPE, my_pe+1, 0, comm);
}
}
/* Routine to start receiving boundary data from neighboring PEs.
Notes:
1) buffer should be able to hold 2*NVARS*MYSUB realtype entries, should be
passed to both the BRecvPost and BRecvWait functions, and should not
be manipulated between the two calls.
2) request should have 4 entries, and should be passed in both calls also. */
static void BRecvPost(MPI_Comm comm, MPI_Request request[],
int my_pe, int isubx, int isuby,
long int dsizex, long int dsizey,
realtype uext[], realtype buffer[])
{
long int offsetue;
/* Have bufleft and bufright use the same buffer */
realtype *bufleft = buffer, *bufright = buffer+NVARS*MYSUB;
/* If isuby > 0, receive data for bottom x-line of uext */
if (isuby != 0)
MPI_Irecv(&uext[NVARS], dsizex, PVEC_REAL_MPI_TYPE,
my_pe-NPEX, 0, comm, &request[0]);
/* If isuby < NPEY-1, receive data for top x-line of uext */
if (isuby != NPEY-1) {
offsetue = NVARS*(1 + (MYSUB+1)*(MXSUB+2));
MPI_Irecv(&uext[offsetue], dsizex, PVEC_REAL_MPI_TYPE,
my_pe+NPEX, 0, comm, &request[1]);
}
/* If isubx > 0, receive data for left y-line of uext (via bufleft) */
if (isubx != 0) {
MPI_Irecv(&bufleft[0], dsizey, PVEC_REAL_MPI_TYPE,
my_pe-1, 0, comm, &request[2]);
}
/* If isubx < NPEX-1, receive data for right y-line of uext (via bufright) */
if (isubx != NPEX-1) {
MPI_Irecv(&bufright[0], dsizey, PVEC_REAL_MPI_TYPE,
my_pe+1, 0, comm, &request[3]);
}
}
/* Routine to finish receiving boundary data from neighboring PEs.
Notes:
1) buffer should be able to hold 2*NVARS*MYSUB realtype entries, should be
passed to both the BRecvPost and BRecvWait functions, and should not
be manipulated between the two calls.
2) request should have 4 entries, and should be passed in both calls also. */
static void BRecvWait(MPI_Request request[],
int isubx, int isuby,
long int dsizex, realtype uext[],
realtype buffer[])
{
int i, ly;
long int dsizex2, offsetue, offsetbuf;
realtype *bufleft = buffer, *bufright = buffer+NVARS*MYSUB;
MPI_Status status;
dsizex2 = dsizex + 2*NVARS;
/* If isuby > 0, receive data for bottom x-line of uext */
if (isuby != 0)
MPI_Wait(&request[0],&status);
/* If isuby < NPEY-1, receive data for top x-line of uext */
if (isuby != NPEY-1)
MPI_Wait(&request[1],&status);
/* If isubx > 0, receive data for left y-line of uext (via bufleft) */
if (isubx != 0) {
MPI_Wait(&request[2],&status);
/* Copy the buffer to uext */
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NVARS;
offsetue = (ly+1)*dsizex2;
for (i = 0; i < NVARS; i++)
uext[offsetue+i] = bufleft[offsetbuf+i];
}
}
/* If isubx < NPEX-1, receive data for right y-line of uext (via bufright) */
if (isubx != NPEX-1) {
MPI_Wait(&request[3],&status);
/* Copy the buffer to uext */
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NVARS;
offsetue = (ly+2)*dsizex2 - NVARS;
for (i = 0; i < NVARS; i++)
uext[offsetue+i] = bufright[offsetbuf+i];
}
}
}
/* ucomm routine. This routine performs all communication
between processors of data needed to calculate f. */
static void ucomm(realtype t, N_Vector u, UserData data)
{
realtype *udata, *uext, buffer[2*NVARS*MYSUB];
MPI_Comm comm;
int my_pe, isubx, isuby;
long int nvmxsub, nvmysub;
MPI_Request request[4];
udata = NV_DATA_P(u);
/* Get comm, my_pe, subgrid indices, data sizes, extended array uext */
comm = data->comm; my_pe = data->my_pe;
isubx = data->isubx; isuby = data->isuby;
nvmxsub = data->nvmxsub;
nvmysub = NVARS*MYSUB;
uext = data->uext;
/* Start receiving boundary data from neighboring PEs */
BRecvPost(comm, request, my_pe, isubx, isuby, nvmxsub, nvmysub, uext, buffer);
/* Send data from boundary of local grid to neighboring PEs */
BSend(comm, my_pe, isubx, isuby, nvmxsub, nvmysub, udata);
/* Finish receiving boundary data from neighboring PEs */
BRecvWait(request, isubx, isuby, nvmxsub, uext, buffer);
}
/* fcalc routine. Compute f(t,y). This routine assumes that communication
between processors of data needed to calculate f has already been done,
and this data is in the work array uext. */
static void fcalc(realtype t, realtype udata[],
realtype dudata[], UserData data)
{
realtype *uext;
realtype q3, c1, c2, c1dn, c2dn, c1up, c2up, c1lt, c2lt;
realtype c1rt, c2rt, cydn, cyup, hord1, hord2, horad1, horad2;
realtype qq1, qq2, qq3, qq4, rkin1, rkin2, s, vertd1, vertd2, ydn, yup;
realtype q4coef, dely, verdco, hordco, horaco;
int i, lx, ly, jx, jy;
int isubx, isuby;
long int nvmxsub, nvmxsub2, offsetu, offsetue;
/* Get subgrid indices, data sizes, extended work array uext */
isubx = data->isubx; isuby = data->isuby;
nvmxsub = data->nvmxsub; nvmxsub2 = data->nvmxsub2;
uext = data->uext;
/* Copy local segment of u vector into the working extended array uext */
offsetu = 0;
offsetue = nvmxsub2 + NVARS;
for (ly = 0; ly < MYSUB; ly++) {
for (i = 0; i < nvmxsub; i++) uext[offsetue+i] = udata[offsetu+i];
offsetu = offsetu + nvmxsub;
offsetue = offsetue + nvmxsub2;
}
/* To facilitate homogeneous Neumann boundary conditions, when this is
a boundary PE, copy data from the first interior mesh line of u to uext */
/* If isuby = 0, copy x-line 2 of u to uext */
if (isuby == 0) {
for (i = 0; i < nvmxsub; i++) uext[NVARS+i] = udata[nvmxsub+i];
}
/* If isuby = NPEY-1, copy x-line MYSUB-1 of u to uext */
if (isuby == NPEY-1) {
offsetu = (MYSUB-2)*nvmxsub;
offsetue = (MYSUB+1)*nvmxsub2 + NVARS;
for (i = 0; i < nvmxsub; i++) uext[offsetue+i] = udata[offsetu+i];
}
/* If isubx = 0, copy y-line 2 of u to uext */
if (isubx == 0) {
for (ly = 0; ly < MYSUB; ly++) {
offsetu = ly*nvmxsub + NVARS;
offsetue = (ly+1)*nvmxsub2;
for (i = 0; i < NVARS; i++) uext[offsetue+i] = udata[offsetu+i];
}
}
/* If isubx = NPEX-1, copy y-line MXSUB-1 of u to uext */
if (isubx == NPEX-1) {
for (ly = 0; ly < MYSUB; ly++) {
offsetu = (ly+1)*nvmxsub - 2*NVARS;
offsetue = (ly+2)*nvmxsub2 - NVARS;
for (i = 0; i < NVARS; i++) uext[offsetue+i] = udata[offsetu+i];
}
}
/* Make local copies of problem variables, for efficiency */
dely = data->dy;
verdco = data->vdco;
hordco = data->hdco;
horaco = data->haco;
/* Set diurnal rate coefficients as functions of t, and save q4 in
data block for use by preconditioner evaluation routine */
s = sin((data->om)*t);
if (s > RCONST(0.0)) {
q3 = EXP(-A3/s);
q4coef = EXP(-A4/s);
} else {
q3 = RCONST(0.0);
q4coef = RCONST(0.0);
}
data->q4 = q4coef;
/* Loop over all grid points in local subgrid */
for (ly = 0; ly < MYSUB; ly++) {
jy = ly + isuby*MYSUB;
/* Set vertical diffusion coefficients at jy +- 1/2 */
ydn = YMIN + (jy - RCONST(0.5))*dely;
yup = ydn + dely;
cydn = verdco*EXP(RCONST(0.2)*ydn);
cyup = verdco*EXP(RCONST(0.2)*yup);
for (lx = 0; lx < MXSUB; lx++) {
jx = lx + isubx*MXSUB;
/* Extract c1 and c2, and set kinetic rate terms */
offsetue = (lx+1)*NVARS + (ly+1)*nvmxsub2;
c1 = uext[offsetue];
c2 = uext[offsetue+1];
qq1 = Q1*c1*C3;
qq2 = Q2*c1*c2;
qq3 = q3*C3;
qq4 = q4coef*c2;
rkin1 = -qq1 - qq2 + RCONST(2.0)*qq3 + qq4;
rkin2 = qq1 - qq2 - qq4;
/* Set vertical diffusion terms */
c1dn = uext[offsetue-nvmxsub2];
c2dn = uext[offsetue-nvmxsub2+1];
c1up = uext[offsetue+nvmxsub2];
c2up = uext[offsetue+nvmxsub2+1];
vertd1 = cyup*(c1up - c1) - cydn*(c1 - c1dn);
vertd2 = cyup*(c2up - c2) - cydn*(c2 - c2dn);
/* Set horizontal diffusion and advection terms */
c1lt = uext[offsetue-2];
c2lt = uext[offsetue-1];
c1rt = uext[offsetue+2];
c2rt = uext[offsetue+3];
hord1 = hordco*(c1rt - RCONST(2.0)*c1 + c1lt);
hord2 = hordco*(c2rt - RCONST(2.0)*c2 + c2lt);
horad1 = horaco*(c1rt - c1lt);
horad2 = horaco*(c2rt - c2lt);
/* Load all terms into dudata */
offsetu = lx*NVARS + ly*nvmxsub;
dudata[offsetu] = vertd1 + hord1 + horad1 + rkin1;
dudata[offsetu+1] = vertd2 + hord2 + horad2 + rkin2;
}
}
}
/***************** Functions Called by the Solver *************************/
/* f routine. Evaluate f(t,y). First call ucomm to do communication of
subgrid boundary data into uext. Then calculate f by a call to fcalc. */
static int f(realtype t, N_Vector u, N_Vector udot, void *user_data)
{
realtype *udata, *dudata;
UserData data;
udata = NV_DATA_P(u);
dudata = NV_DATA_P(udot);
data = (UserData) user_data;
/* Call ucomm to do inter-processor communication */
ucomm(t, u, data);
/* Call fcalc to calculate all right-hand sides */
fcalc(t, udata, dudata, data);
return(0);
}
/* Preconditioner setup routine. Generate and preprocess P. */
static int Precond(realtype tn, N_Vector u, N_Vector fu,
booleantype jok, booleantype *jcurPtr,
realtype gamma, void *user_data,
N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
{
realtype c1, c2, cydn, cyup, diag, ydn, yup, q4coef, dely, verdco, hordco;
realtype **(*P)[MYSUB], **(*Jbd)[MYSUB];
long int nvmxsub, offset;
long int *(*pivot)[MYSUB], ier;
int lx, ly, jx, jy, isubx, isuby;
realtype *udata, **a, **j;
UserData data;
/* Make local copies of pointers in user_data, pointer to u's data,
and PE index pair */
data = (UserData) user_data;
P = data->P;
Jbd = data->Jbd;
pivot = data->pivot;
udata = NV_DATA_P(u);
isubx = data->isubx; isuby = data->isuby;
nvmxsub = data->nvmxsub;
if (jok) {
/* jok = TRUE: Copy Jbd to P */
for (ly = 0; ly < MYSUB; ly++)
for (lx = 0; lx < MXSUB; lx++)
denseCopy(Jbd[lx][ly], P[lx][ly], NVARS, NVARS);
*jcurPtr = FALSE;
}
else {
/* jok = FALSE: Generate Jbd from scratch and copy to P */
/* Make local copies of problem variables, for efficiency */
q4coef = data->q4;
dely = data->dy;
verdco = data->vdco;
hordco = data->hdco;
/* Compute 2x2 diagonal Jacobian blocks (using q4 values
computed on the last f call). Load into P. */
for (ly = 0; ly < MYSUB; ly++) {
jy = ly + isuby*MYSUB;
ydn = YMIN + (jy - RCONST(0.5))*dely;
yup = ydn + dely;
cydn = verdco*EXP(RCONST(0.2)*ydn);
cyup = verdco*EXP(RCONST(0.2)*yup);
diag = -(cydn + cyup + RCONST(2.0)*hordco);
for (lx = 0; lx < MXSUB; lx++) {
jx = lx + isubx*MXSUB;
offset = lx*NVARS + ly*nvmxsub;
c1 = udata[offset];
c2 = udata[offset+1];
j = Jbd[lx][ly];
a = P[lx][ly];
IJth(j,1,1) = (-Q1*C3 - Q2*c2) + diag;
IJth(j,1,2) = -Q2*c1 + q4coef;
IJth(j,2,1) = Q1*C3 - Q2*c2;
IJth(j,2,2) = (-Q2*c1 - q4coef) + diag;
denseCopy(j, a, NVARS, NVARS);
}
}
*jcurPtr = TRUE;
}
/* Scale by -gamma */
for (ly = 0; ly < MYSUB; ly++)
for (lx = 0; lx < MXSUB; lx++)
denseScale(-gamma, P[lx][ly], NVARS, NVARS);
/* Add identity matrix and do LU decompositions on blocks in place */
for (lx = 0; lx < MXSUB; lx++) {
for (ly = 0; ly < MYSUB; ly++) {
denseAddIdentity(P[lx][ly], NVARS);
ier = denseGETRF(P[lx][ly], NVARS, NVARS, pivot[lx][ly]);
if (ier != 0) return(1);
}
}
return(0);
}
/* Preconditioner solve routine */
static int PSolve(realtype tn, N_Vector u, N_Vector fu,
N_Vector r, N_Vector z,
realtype gamma, realtype delta,
int lr, void *user_data, N_Vector vtemp)
{
realtype **(*P)[MYSUB];
long int nvmxsub;
long int *(*pivot)[MYSUB];
int lx, ly;
realtype *zdata, *v;
UserData data;
/* Extract the P and pivot arrays from user_data */
data = (UserData) user_data;
P = data->P;
pivot = data->pivot;
/* Solve the block-diagonal system Px = r using LU factors stored
in P and pivot data in pivot, and return the solution in z.
First copy vector r to z. */
N_VScale(RCONST(1.0), r, z);
nvmxsub = data->nvmxsub;
zdata = NV_DATA_P(z);
for (lx = 0; lx < MXSUB; lx++) {
for (ly = 0; ly < MYSUB; ly++) {
v = &(zdata[lx*NVARS + ly*nvmxsub]);
denseGETRS(P[lx][ly], NVARS, pivot[lx][ly], v);
}
}
return(0);
}
/*********************** Private Helper Function ************************/
/* Check function return value...
opt == 0 means SUNDIALS function allocates memory so check if
returned NULL pointer
opt == 1 means SUNDIALS function returns a flag so check if
flag >= 0
opt == 2 means function allocates memory so check if returned
NULL pointer */
static int check_flag(void *flagvalue, char *funcname, int opt, int id)
{
int *errflag;
/* Check if SUNDIALS function returned NULL pointer - no memory allocated */
if (opt == 0 && flagvalue == NULL) {
fprintf(stderr, "\nSUNDIALS_ERROR(%d): %s() failed - returned NULL pointer\n\n",
id, funcname);
return(1); }
/* Check if flag < 0 */
else if (opt == 1) {
errflag = (int *) flagvalue;
if (*errflag < 0) {
fprintf(stderr, "\nSUNDIALS_ERROR(%d): %s() failed with flag = %d\n\n",
id, funcname, *errflag);
return(1); }}
/* Check if function returned NULL pointer - no memory allocated */
else if (opt == 2 && flagvalue == NULL) {
fprintf(stderr, "\nMEMORY_ERROR(%d): %s() failed - returned NULL pointer\n\n",
id, funcname);
return(1); }
return(0);
}
|