/usr/share/doc/libsundials-serial-dev/examples/idas/parallel/idasFoodWeb_kry_bbd_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 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 | /*
* -----------------------------------------------------------------
* $Revision: 1.2 $
* $Date: 2010/12/01 23:06:37 $
* -----------------------------------------------------------------
* Programmer(s): Allan Taylor, Alan Hindmarsh and
* Radu Serban @ LLNL
* -----------------------------------------------------------------
* Example program for IDA: Food web, parallel, GMRES, IDABBD
* preconditioner.
*
* This example program for IDAS uses IDASPGMR as the linear solver.
* It is written for a parallel computer system and uses the
* IDABBDPRE band-block-diagonal preconditioner module for the
* IDASPGMR package. It was originally run on a Sun SPARC cluster
* and used MPICH.
*
* The mathematical problem solved in this example is a DAE system
* that arises from a system of partial differential equations after
* spatial discretization. The PDE system is a food web population
* model, with predator-prey interaction and diffusion on the unit
* square in two dimensions. The dependent variable vector is:
*
* 1 2 ns
* c = (c , c , ..., c ) , ns = 2 * np
*
* and the PDE's are as follows:
*
* i i i
* dc /dt = d(i)*(c + c ) + R (x,y,c) (i = 1,...,np)
* xx yy i
*
* i i
* 0 = d(i)*(c + c ) + R (x,y,c) (i = np+1,...,ns)
* xx yy i
*
* where the reaction terms R are:
*
* i ns j
* R (x,y,c) = c * (b(i) + sum a(i,j)*c )
* i j=1
*
* The number of species is ns = 2 * np, with the first np being
* prey and the last np being predators. The coefficients a(i,j),
* b(i), d(i) are:
*
* a(i,i) = -AA (all i)
* a(i,j) = -GG (i <= np , j > np)
* a(i,j) = EE (i > np, j <= np)
* all other a(i,j) = 0
* b(i) = BB*(1+ alpha * x*y + beta*sin(4 pi x)*sin(4 pi y)) (i <= np)
* b(i) =-BB*(1+ alpha * x*y + beta*sin(4 pi x)*sin(4 pi y)) (i > np)
* d(i) = DPREY (i <= np)
* d(i) = DPRED (i > np)
*
* Note: The above equations are written in 1-based indices,
* whereas the code has 0-based indices, being written in C.
*
* The various scalar parameters required are set using '#define'
* statements or directly in routine InitUserData. In this program,
* np = 1, ns = 2. The boundary conditions are homogeneous Neumann:
* normal derivative = 0.
*
* A polynomial in x and y is used to set the initial values of the
* first np variables (the prey variables) at each x,y location,
* while initial values for the remaining (predator) variables are
* set to a flat value, which is corrected by IDACalcIC.
*
* The PDEs are discretized by central differencing on a MX by MY
* mesh, and so the system size Neq is the product
* MX * MY * NUM_SPECIES. The system is actually implemented on
* submeshes, processor by processor, with an MXSUB by MYSUB mesh
* on each of NPEX * NPEY processors.
*
* The DAE system is solved by IDAS using the IDASPGMR linear solver,
* in conjunction with the preconditioner module IDABBDPRE. The
* preconditioner uses a 5-diagonal band-block-diagonal
* approximation (half-bandwidths = 2). Output is printed at
* t = 0, .001, .01, .1, .4, .7, 1.
* -----------------------------------------------------------------
* References:
* [1] Peter N. Brown and Alan C. Hindmarsh,
* Reduced Storage Matrix Methods in Stiff ODE systems,
* Journal of Applied Mathematics and Computation, Vol. 31
* (May 1989), pp. 40-91.
*
* [2] Peter N. Brown, Alan C. Hindmarsh, and Linda R. Petzold,
* Using Krylov Methods in the Solution of Large-Scale
* Differential-Algebraic Systems, SIAM J. Sci. Comput., 15
* (1994), pp. 1467-1488.
*
* [3] Peter N. Brown, Alan C. Hindmarsh, and Linda R. Petzold,
* Consistent Initial Condition Calculation for Differential-
* Algebraic Systems, SIAM J. Sci. Comput., 19 (1998),
* pp. 1495-1512.
* -----------------------------------------------------------------
*/
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <idas/idas.h>
#include <idas/idas_spgmr.h>
#include <idas/idas_bbdpre.h>
#include <nvector/nvector_parallel.h>
#include <sundials/sundials_dense.h>
#include <sundials/sundials_types.h>
#include <sundials/sundials_math.h>
#include <mpi.h>
/* Problem Constants */
#define NPREY 1 /* Number of prey (= number of predators). */
#define NUM_SPECIES 2*NPREY
#define PI RCONST(3.1415926535898) /* pi */
#define FOURPI (RCONST(4.0)*PI) /* 4 pi */
#define MXSUB 10 /* Number of x mesh points per processor subgrid */
#define MYSUB 10 /* Number of y mesh points per processor subgrid */
#define NPEX 2 /* Number of subgrids in the x direction */
#define NPEY 2 /* Number of subgrids in the y direction */
#define MX (MXSUB*NPEX) /* MX = number of x mesh points */
#define MY (MYSUB*NPEY) /* MY = number of y mesh points */
#define NSMXSUB (NUM_SPECIES * MXSUB)
#define NEQ (NUM_SPECIES*MX*MY) /* Number of equations in system */
#define AA RCONST(1.0) /* Coefficient in above eqns. for a */
#define EE RCONST(10000.) /* Coefficient in above eqns. for a */
#define GG RCONST(0.5e-6) /* Coefficient in above eqns. for a */
#define BB RCONST(1.0) /* Coefficient in above eqns. for b */
#define DPREY RCONST(1.0) /* Coefficient in above eqns. for d */
#define DPRED RCONST(0.05) /* Coefficient in above eqns. for d */
#define ALPHA RCONST(50.) /* Coefficient alpha in above eqns. */
#define BETA RCONST(1000.) /* Coefficient beta in above eqns. */
#define AX RCONST(1.0) /* Total range of x variable */
#define AY RCONST(1.0) /* Total range of y variable */
#define RTOL RCONST(1.e-5) /* rtol tolerance */
#define ATOL RCONST(1.e-5) /* atol tolerance */
#define ZERO RCONST(0.) /* 0. */
#define ONE RCONST(1.0) /* 1. */
#define NOUT 6
#define TMULT RCONST(10.0) /* Multiplier for tout values */
#define TADD RCONST(0.3) /* Increment for tout values */
/* User-defined vector accessor macro IJ_Vptr. */
/*
* IJ_Vptr is defined in order to express the underlying 3-d structure of the
* dependent variable vector from its underlying 1-d storage (an N_Vector).
* IJ_Vptr(vv,i,j) returns a pointer to the location in vv corresponding to
* species index is = 0, x-index ix = i, and y-index jy = j.
*/
#define IJ_Vptr(vv,i,j) (&NV_Ith_P(vv, (i)*NUM_SPECIES + (j)*NSMXSUB ))
/* Type: UserData. Contains problem constants, preconditioner data, etc. */
typedef struct {
int ns, np, thispe, npes, ixsub, jysub, npex, npey;
int mxsub, mysub, nsmxsub, nsmxsub2;
realtype dx, dy, **acoef;
realtype cox[NUM_SPECIES], coy[NUM_SPECIES], bcoef[NUM_SPECIES],
rhs[NUM_SPECIES], cext[(MXSUB+2)*(MYSUB+2)*NUM_SPECIES];
MPI_Comm comm;
N_Vector rates;
long int n_local;
} *UserData;
/* Prototypes for functions called by the IDA Solver. */
static int resweb(realtype tt,
N_Vector cc, N_Vector cp, N_Vector rr,
void *user_data);
static int reslocal(long int Nlocal, realtype tt,
N_Vector cc, N_Vector cp, N_Vector res,
void *user_data);
static int rescomm(long int Nlocal, realtype tt,
N_Vector cc, N_Vector cp,
void *user_data);
/* Prototypes for supporting functions */
static void BSend(MPI_Comm comm, int thispe, int ixsub, int jysub,
int dsizex, int dsizey, realtype carray[]);
static void BRecvPost(MPI_Comm comm, MPI_Request request[], int thispe,
int ixsub, int jysub,
int dsizex, int dsizey,
realtype cext[], realtype buffer[]);
static void BRecvWait(MPI_Request request[], int ixsub, int jysub,
int dsizex, realtype cext[], realtype buffer[]);
static void WebRates(realtype xx, realtype yy, realtype *cxy, realtype *ratesxy,
UserData webdata);
static realtype dotprod(int size, realtype *x1, realtype *x2);
/* Prototypes for private functions */
static void InitUserData(UserData webdata, int thispe, int npes,
MPI_Comm comm);
static void SetInitialProfiles(N_Vector cc, N_Vector cp, N_Vector id,
N_Vector scrtch, UserData webdata);
static void PrintHeader(int SystemSize, int maxl,
int mudq, int mldq,
int mukeep, int mlkeep,
realtype rtol, realtype atol);
static void PrintOutput(void *mem, N_Vector cc, realtype time,
UserData webdata, MPI_Comm comm);
static void PrintFinalStats(void *mem);
static int check_flag(void *flagvalue, char *funcname, int opt, int id);
/*
*--------------------------------------------------------------------
* MAIN PROGRAM
*--------------------------------------------------------------------
*/
int main(int argc, char *argv[])
{
MPI_Comm comm;
void *mem;
UserData webdata;
long int SystemSize, local_N, mudq, mldq, mukeep, mlkeep;
realtype rtol, atol, t0, tout, tret;
N_Vector cc, cp, res, id;
int thispe, npes, maxl, iout, retval;
cc = cp = res = id = NULL;
webdata = NULL;
mem = NULL;
/* Set communicator, and get processor number and total number of PE's. */
MPI_Init(&argc, &argv);
comm = MPI_COMM_WORLD;
MPI_Comm_rank(comm, &thispe);
MPI_Comm_size(comm, &npes);
if (npes != NPEX*NPEY) {
if (thispe == 0)
fprintf(stderr,
"\nMPI_ERROR(0): npes = %d not equal to NPEX*NPEY = %d\n",
npes, NPEX*NPEY);
MPI_Finalize();
return(1);
}
/* Set local length (local_N) and global length (SystemSize). */
local_N = MXSUB*MYSUB*NUM_SPECIES;
SystemSize = NEQ;
/* Set up user data block webdata. */
webdata = (UserData) malloc(sizeof *webdata);
webdata->rates = N_VNew_Parallel(comm, local_N, SystemSize);
webdata->acoef = newDenseMat(NUM_SPECIES, NUM_SPECIES);
InitUserData(webdata, thispe, npes, comm);
/* Create needed vectors, and load initial values.
The vector res is used temporarily only. */
cc = N_VNew_Parallel(comm, local_N, SystemSize);
if(check_flag((void *)cc, "N_VNew_Parallel", 0, thispe)) MPI_Abort(comm, 1);
cp = N_VNew_Parallel(comm, local_N, SystemSize);
if(check_flag((void *)cp, "N_VNew_Parallel", 0, thispe)) MPI_Abort(comm, 1);
res = N_VNew_Parallel(comm, local_N, SystemSize);
if(check_flag((void *)res, "N_VNew_Parallel", 0, thispe)) MPI_Abort(comm, 1);
id = N_VNew_Parallel(comm, local_N, SystemSize);
if(check_flag((void *)id, "N_VNew_Parallel", 0, thispe)) MPI_Abort(comm, 1);
SetInitialProfiles(cc, cp, id, res, webdata);
N_VDestroy_Parallel(res);
/* Set remaining inputs to IDAMalloc. */
t0 = ZERO;
rtol = RTOL;
atol = ATOL;
/* Call IDACreate and IDAMalloc to initialize solution */
mem = IDACreate();
if(check_flag((void *)mem, "IDACreate", 0, thispe)) MPI_Abort(comm, 1);
retval = IDASetUserData(mem, webdata);
if(check_flag(&retval, "IDASetUserData", 1, thispe)) MPI_Abort(comm, 1);
retval = IDASetId(mem, id);
if(check_flag(&retval, "IDASetId", 1, thispe)) MPI_Abort(comm, 1);
retval = IDAInit(mem, resweb, t0, cc, cp);
if(check_flag(&retval, "IDAInit", 1, thispe)) MPI_Abort(comm, 1);
retval = IDASStolerances(mem, rtol, atol);
if(check_flag(&retval, "IDASStolerances", 1, thispe)) MPI_Abort(comm, 1);
/* Call IDASpgmr to specify the IDA linear solver IDASPGMR */
maxl = 16;
retval = IDASpgmr(mem, maxl);
if(check_flag(&retval, "IDASpgmr", 1, thispe)) MPI_Abort(comm, 1);
/* Call IDABBDPrecInit to initialize the band-block-diagonal preconditioner.
The half-bandwidths for the difference quotient evaluation are exact
for the system Jacobian, but only a 5-diagonal band matrix is retained. */
mudq = mldq = NSMXSUB;
mukeep = mlkeep = 2;
retval = IDABBDPrecInit(mem, local_N, mudq, mldq, mukeep, mlkeep,
ZERO, reslocal, NULL);
if(check_flag(&retval, "IDABBDPrecInit", 1, thispe)) MPI_Abort(comm, 1);
/* Call IDACalcIC (with default options) to correct the initial values. */
tout = RCONST(0.001);
retval = IDACalcIC(mem, IDA_YA_YDP_INIT, tout);
if(check_flag(&retval, "IDACalcIC", 1, thispe)) MPI_Abort(comm, 1);
/* On PE 0, print heading, basic parameters, initial values. */
if (thispe == 0) PrintHeader(SystemSize, maxl,
mudq, mldq, mukeep, mlkeep,
rtol, atol);
PrintOutput(mem, cc, t0, webdata, comm);
/* Call IDA in tout loop, normal mode, and print selected output. */
for (iout = 1; iout <= NOUT; iout++) {
retval = IDASolve(mem, tout, &tret, cc, cp, IDA_NORMAL);
if(check_flag(&retval, "IDASolve", 1, thispe)) MPI_Abort(comm, 1);
PrintOutput(mem, cc, tret, webdata, comm);
if (iout < 3) tout *= TMULT;
else tout += TADD;
}
/* On PE 0, print final set of statistics. */
if (thispe == 0) PrintFinalStats(mem);
/* Free memory. */
N_VDestroy_Parallel(cc);
N_VDestroy_Parallel(cp);
N_VDestroy_Parallel(id);
IDAFree(&mem);
destroyMat(webdata->acoef);
N_VDestroy_Parallel(webdata->rates);
free(webdata);
MPI_Finalize();
return(0);
}
/*
*--------------------------------------------------------------------
* PRIVATE FUNCTIONS
*--------------------------------------------------------------------
*/
/*
* InitUserData: Load problem constants in webdata (of type UserData).
*/
static void InitUserData(UserData webdata, int thispe, int npes,
MPI_Comm comm)
{
int i, j, np;
realtype *a1,*a2, *a3, *a4, dx2, dy2, **acoef, *bcoef, *cox, *coy;
webdata->jysub = thispe / NPEX;
webdata->ixsub = thispe - (webdata->jysub)*NPEX;
webdata->mxsub = MXSUB;
webdata->mysub = MYSUB;
webdata->npex = NPEX;
webdata->npey = NPEY;
webdata->ns = NUM_SPECIES;
webdata->np = NPREY;
webdata->dx = AX/(MX-1);
webdata->dy = AY/(MY-1);
webdata->thispe = thispe;
webdata->npes = npes;
webdata->nsmxsub = MXSUB * NUM_SPECIES;
webdata->nsmxsub2 = (MXSUB+2)*NUM_SPECIES;
webdata->comm = comm;
webdata->n_local = MXSUB*MYSUB*NUM_SPECIES;
/* Set up the coefficients a and b plus others found in the equations. */
np = webdata->np;
dx2 = (webdata->dx)*(webdata->dx);
dy2 = (webdata->dy)*(webdata->dy);
acoef = webdata->acoef;
bcoef = webdata->bcoef;
cox = webdata->cox;
coy = webdata->coy;
for (i = 0; i < np; i++) {
a1 = &(acoef[i][np]);
a2 = &(acoef[i+np][0]);
a3 = &(acoef[i][0]);
a4 = &(acoef[i+np][np]);
/* Fill in the portion of acoef in the four quadrants, row by row. */
for (j = 0; j < np; j++) {
*a1++ = -GG;
*a2++ = EE;
*a3++ = ZERO;
*a4++ = ZERO;
}
/* Reset the diagonal elements of acoef to -AA. */
acoef[i][i] = -AA; acoef[i+np][i+np] = -AA;
/* Set coefficients for b and diffusion terms. */
bcoef[i] = BB; bcoef[i+np] = -BB;
cox[i] = DPREY/dx2; cox[i+np] = DPRED/dx2;
coy[i] = DPREY/dy2; coy[i+np] = DPRED/dy2;
}
}
/*
* SetInitialProfiles: Set initial conditions in cc, cp, and id.
* A polynomial profile is used for the prey cc values, and a constant
* (1.0e5) is loaded as the initial guess for the predator cc values.
* The id values are set to 1 for the prey and 0 for the predators.
* The prey cp values are set according to the given system, and
* the predator cp values are set to zero.
*/
static void SetInitialProfiles(N_Vector cc, N_Vector cp, N_Vector id,
N_Vector res, UserData webdata)
{
int ixsub, jysub, mxsub, mysub, nsmxsub, np, ix, jy, is;
realtype *cxy, *idxy, *cpxy, dx, dy, xx, yy, xyfactor;
ixsub = webdata->ixsub;
jysub = webdata->jysub;
mxsub = webdata->mxsub;
mysub = webdata->mxsub;
nsmxsub = webdata->nsmxsub;
dx = webdata->dx;
dy = webdata->dy;
np = webdata->np;
/* Loop over grid, load cc values and id values. */
for (jy = 0; jy < mysub; jy++) {
yy = (jy + jysub*mysub) * dy;
for (ix = 0; ix < mxsub; ix++) {
xx = (ix + ixsub*mxsub) * dx;
xyfactor = 16.*xx*(1. - xx)*yy*(1. - yy);
xyfactor *= xyfactor;
cxy = IJ_Vptr(cc,ix,jy);
idxy = IJ_Vptr(id,ix,jy);
for (is = 0; is < NUM_SPECIES; is++) {
if (is < np) { cxy[is] = RCONST(10.0) + (realtype)(is+1)*xyfactor; idxy[is] = ONE; }
else { cxy[is] = 1.0e5; idxy[is] = ZERO; }
}
}
}
/* Set c' for the prey by calling the residual function with cp = 0. */
N_VConst(ZERO, cp);
resweb(ZERO, cc, cp, res, webdata);
N_VScale(-ONE, res, cp);
/* Set c' for predators to 0. */
for (jy = 0; jy < mysub; jy++) {
for (ix = 0; ix < mxsub; ix++) {
cpxy = IJ_Vptr(cp,ix,jy);
for (is = np; is < NUM_SPECIES; is++) cpxy[is] = ZERO;
}
}
}
/*
* Print first lines of output (problem description)
* and table headerr
*/
static void PrintHeader(int SystemSize, int maxl,
int mudq, int mldq,
int mukeep, int mlkeep,
realtype rtol, realtype atol)
{
printf("\nidasFoodWeb_kry_bbd_p: Predator-prey DAE parallel example problem for IDA \n\n");
printf("Number of species ns: %d", NUM_SPECIES);
printf(" Mesh dimensions: %d x %d", MX, MY);
printf(" Total system size: %d\n",SystemSize);
printf("Subgrid dimensions: %d x %d", MXSUB, MYSUB);
printf(" Processor array: %d x %d\n", NPEX, NPEY);
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf("Tolerance parameters: rtol = %Lg atol = %Lg\n", rtol, atol);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
printf("Tolerance parameters: rtol = %lg atol = %lg\n", rtol, atol);
#else
printf("Tolerance parameters: rtol = %g atol = %g\n", rtol, atol);
#endif
printf("Linear solver: IDASPGMR Max. Krylov dimension maxl: %d\n", maxl);
printf("Preconditioner: band-block-diagonal (IDABBDPRE), with parameters\n");
printf(" mudq = %d, mldq = %d, mukeep = %d, mlkeep = %d\n",
mudq, mldq, mukeep, mlkeep);
printf("CalcIC called to correct initial predator concentrations \n\n");
printf("-----------------------------------------------------------\n");
printf(" t bottom-left top-right");
printf(" | nst k h\n");
printf("-----------------------------------------------------------\n\n");
}
/*
* PrintOutput: Print output values at output time t = tt.
* Selected run statistics are printed. Then values of c1 and c2
* are printed for the bottom left and top right grid points only.
*/
static void PrintOutput(void *mem, N_Vector cc, realtype tt,
UserData webdata, MPI_Comm comm)
{
MPI_Status status;
realtype *cdata, clast[2], hused;
long int nst;
int i, kused, flag, thispe, npelast, ilast;;
thispe = webdata->thispe;
npelast = webdata->npes - 1;
cdata = NV_DATA_P(cc);
/* Send conc. at top right mesh point from PE npes-1 to PE 0. */
if (thispe == npelast) {
ilast = NUM_SPECIES*MXSUB*MYSUB - 2;
if (npelast != 0)
MPI_Send(&cdata[ilast], 2, PVEC_REAL_MPI_TYPE, 0, 0, comm);
else { clast[0] = cdata[ilast]; clast[1] = cdata[ilast+1]; }
}
/* On PE 0, receive conc. at top right from PE npes - 1.
Then print performance data and sampled solution values. */
if (thispe == 0) {
if (npelast != 0)
MPI_Recv(&clast[0], 2, PVEC_REAL_MPI_TYPE, npelast, 0, comm, &status);
flag = IDAGetLastOrder(mem, &kused);
check_flag(&flag, "IDAGetLastOrder", 1, thispe);
flag = IDAGetNumSteps(mem, &nst);
check_flag(&flag, "IDAGetNumSteps", 1, thispe);
flag = IDAGetLastStep(mem, &hused);
check_flag(&flag, "IDAGetLastStep", 1, thispe);
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf("%8.2Le %12.4Le %12.4Le | %3ld %1d %12.4Le\n",
tt, cdata[0], clast[0], nst, kused, hused);
for (i=1;i<NUM_SPECIES;i++)
printf(" %12.4Le %12.4Le |\n",cdata[i],clast[i]);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
printf("%8.2le %12.4le %12.4le | %3ld %1d %12.4le\n",
tt, cdata[0], clast[0], nst, kused, hused);
for (i=1;i<NUM_SPECIES;i++)
printf(" %12.4le %12.4le |\n",cdata[i],clast[i]);
#else
printf("%8.2e %12.4e %12.4e | %3ld %1d %12.4e\n",
tt, cdata[0], clast[0], nst, kused, hused);
for (i=1;i<NUM_SPECIES;i++)
printf(" %12.4e %12.4e |\n",cdata[i],clast[i]);
#endif
printf("\n");
}
}
/*
* PrintFinalStats: Print final run data contained in iopt.
*/
static void PrintFinalStats(void *mem)
{
long int nst, nre, nreLS, netf, ncfn, nni, ncfl, nli, npe, nps, nge;
int flag;
flag = IDAGetNumSteps(mem, &nst);
check_flag(&flag, "IDAGetNumSteps", 1, 0);
flag = IDAGetNumResEvals(mem, &nre);
check_flag(&flag, "IDAGetNumResEvals", 1, 0);
flag = IDAGetNumErrTestFails(mem, &netf);
check_flag(&flag, "IDAGetNumErrTestFails", 1, 0);
flag = IDAGetNumNonlinSolvConvFails(mem, &ncfn);
check_flag(&flag, "IDAGetNumNonlinSolvConvFails", 1, 0);
flag = IDAGetNumNonlinSolvIters(mem, &nni);
check_flag(&flag, "IDAGetNumNonlinSolvIters", 1, 0);
flag = IDASpilsGetNumConvFails(mem, &ncfl);
check_flag(&flag, "IDASpilsGetNumConvFails", 1, 0);
flag = IDASpilsGetNumLinIters(mem, &nli);
check_flag(&flag, "IDASpilsGetNumLinIters", 1, 0);
flag = IDASpilsGetNumPrecEvals(mem, &npe);
check_flag(&flag, "IDASpilsGetNumPrecEvals", 1, 0);
flag = IDASpilsGetNumPrecSolves(mem, &nps);
check_flag(&flag, "IDASpilsGetNumPrecSolves", 1, 0);
flag = IDASpilsGetNumResEvals(mem, &nreLS);
check_flag(&flag, "IDASpilsGetNumResEvals", 1, 0);
flag = IDABBDPrecGetNumGfnEvals(mem, &nge);
check_flag(&flag, "IDABBDPrecGetNumGfnEvals", 1, 0);
printf("-----------------------------------------------------------\n");
printf("\nFinal statistics: \n\n");
printf("Number of steps = %ld\n", nst);
printf("Number of residual evaluations = %ld\n", nre+nreLS);
printf("Number of nonlinear iterations = %ld\n", nni);
printf("Number of error test failures = %ld\n", netf);
printf("Number of nonlinear conv. failures = %ld\n\n", ncfn);
printf("Number of linear iterations = %ld\n", nli);
printf("Number of linear conv. failures = %ld\n\n", ncfl);
printf("Number of preconditioner setups = %ld\n", npe);
printf("Number of preconditioner solves = %ld\n", nps);
printf("Number of local residual evals. = %ld\n", nge);
}
/*
* 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;
if (opt == 0 && flagvalue == NULL) {
/* Check if SUNDIALS function returned NULL pointer - no memory allocated */
fprintf(stderr,
"\nSUNDIALS_ERROR(%d): %s() failed - returned NULL pointer\n\n",
id, funcname);
return(1);
} else if (opt == 1) {
/* Check if flag < 0 */
errflag = (int *) flagvalue;
if (*errflag < 0) {
fprintf(stderr,
"\nSUNDIALS_ERROR(%d): %s() failed with flag = %d\n\n",
id, funcname, *errflag);
return(1);
}
} else if (opt == 2 && flagvalue == NULL) {
/* Check if function returned NULL pointer - no memory allocated */
fprintf(stderr,
"\nMEMORY_ERROR(%d): %s() failed - returned NULL pointer\n\n",
id, funcname);
return(1);
}
return(0);
}
/*
*--------------------------------------------------------------------
* FUNCTIONS CALLED BY IDA & SUPPORTING FUNCTIONS
*--------------------------------------------------------------------
*/
/*
* resweb: System residual function for predator-prey system.
* To compute the residual function F, this routine calls:
* rescomm, for needed communication, and then
* reslocal, for computation of the residuals on this processor.
*/
static int resweb(realtype tt,
N_Vector cc, N_Vector cp, N_Vector rr,
void *user_data)
{
int retval;
UserData webdata;
long int Nlocal;
webdata = (UserData) user_data;
Nlocal = webdata->n_local;
/* Call rescomm to do inter-processor communication. */
retval = rescomm(Nlocal, tt, cc, cp, user_data);
/* Call reslocal to calculate the local portion of residual vector. */
retval = reslocal(Nlocal, tt, cc, cp, rr, user_data);
return(0);
}
/*
* rescomm: Communication routine in support of resweb.
* This routine performs all inter-processor communication of components
* of the cc vector needed to calculate F, namely the components at all
* interior subgrid boundaries (ghost cell data). It loads this data
* into a work array cext (the local portion of c, extended).
* The message-passing uses blocking sends, non-blocking receives,
* and receive-waiting, in routines BRecvPost, BSend, BRecvWait.
*/
static int rescomm(long int Nlocal, realtype tt,
N_Vector cc, N_Vector cp,
void *user_data)
{
UserData webdata;
realtype *cdata, *cext, buffer[2*NUM_SPECIES*MYSUB];
int thispe, ixsub, jysub, nsmxsub, nsmysub;
MPI_Comm comm;
MPI_Request request[4];
webdata = (UserData) user_data;
cdata = NV_DATA_P(cc);
/* Get comm, thispe, subgrid indices, data sizes, extended array cext. */
comm = webdata->comm;
thispe = webdata->thispe;
ixsub = webdata->ixsub;
jysub = webdata->jysub;
cext = webdata->cext;
nsmxsub = webdata->nsmxsub;
nsmysub = (webdata->ns)*(webdata->mysub);
/* Start receiving boundary data from neighboring PEs. */
BRecvPost(comm, request, thispe, ixsub, jysub, nsmxsub, nsmysub,
cext, buffer);
/* Send data from boundary of local grid to neighboring PEs. */
BSend(comm, thispe, ixsub, jysub, nsmxsub, nsmysub, cdata);
/* Finish receiving boundary data from neighboring PEs. */
BRecvWait(request, ixsub, jysub, nsmxsub, cext, buffer);
return(0);
}
/*
* BRecvPost: Start receiving boundary data from neighboring PEs.
* (1) buffer should be able to hold 2*NUM_SPECIES*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 is also passed in both calls.
*/
static void BRecvPost(MPI_Comm comm, MPI_Request request[], int my_pe,
int ixsub, int jysub,
int dsizex, int dsizey,
realtype cext[], realtype buffer[])
{
int offsetce;
/* Have bufleft and bufright use the same buffer. */
realtype *bufleft = buffer, *bufright = buffer+NUM_SPECIES*MYSUB;
/* If jysub > 0, receive data for bottom x-line of cext. */
if (jysub != 0)
MPI_Irecv(&cext[NUM_SPECIES], dsizex, PVEC_REAL_MPI_TYPE,
my_pe-NPEX, 0, comm, &request[0]);
/* If jysub < NPEY-1, receive data for top x-line of cext. */
if (jysub != NPEY-1) {
offsetce = NUM_SPECIES*(1 + (MYSUB+1)*(MXSUB+2));
MPI_Irecv(&cext[offsetce], dsizex, PVEC_REAL_MPI_TYPE,
my_pe+NPEX, 0, comm, &request[1]);
}
/* If ixsub > 0, receive data for left y-line of cext (via bufleft). */
if (ixsub != 0) {
MPI_Irecv(&bufleft[0], dsizey, PVEC_REAL_MPI_TYPE,
my_pe-1, 0, comm, &request[2]);
}
/* If ixsub < NPEX-1, receive data for right y-line of cext (via bufright). */
if (ixsub != NPEX-1) {
MPI_Irecv(&bufright[0], dsizey, PVEC_REAL_MPI_TYPE,
my_pe+1, 0, comm, &request[3]);
}
}
/*
* BRecvWait: Finish receiving boundary data from neighboring PEs.
* (1) buffer should be able to hold 2*NUM_SPECIES*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 is also passed in both calls.
*/
static void BRecvWait(MPI_Request request[], int ixsub, int jysub,
int dsizex, realtype cext[], realtype buffer[])
{
int i;
int ly, dsizex2, offsetce, offsetbuf;
realtype *bufleft = buffer, *bufright = buffer+NUM_SPECIES*MYSUB;
MPI_Status status;
dsizex2 = dsizex + 2*NUM_SPECIES;
/* If jysub > 0, receive data for bottom x-line of cext. */
if (jysub != 0)
MPI_Wait(&request[0],&status);
/* If jysub < NPEY-1, receive data for top x-line of cext. */
if (jysub != NPEY-1)
MPI_Wait(&request[1],&status);
/* If ixsub > 0, receive data for left y-line of cext (via bufleft). */
if (ixsub != 0) {
MPI_Wait(&request[2],&status);
/* Copy the buffer to cext */
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NUM_SPECIES;
offsetce = (ly+1)*dsizex2;
for (i = 0; i < NUM_SPECIES; i++)
cext[offsetce+i] = bufleft[offsetbuf+i];
}
}
/* If ixsub < NPEX-1, receive data for right y-line of cext (via bufright). */
if (ixsub != NPEX-1) {
MPI_Wait(&request[3],&status);
/* Copy the buffer to cext */
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NUM_SPECIES;
offsetce = (ly+2)*dsizex2 - NUM_SPECIES;
for (i = 0; i < NUM_SPECIES; i++)
cext[offsetce+i] = bufright[offsetbuf+i];
}
}
}
/*
* BSend: Send boundary data to neighboring PEs.
* This routine sends components of cc from internal subgrid boundaries
* to the appropriate neighbor PEs.
*/
static void BSend(MPI_Comm comm, int my_pe, int ixsub, int jysub,
int dsizex, int dsizey, realtype cdata[])
{
int i;
int ly, offsetc, offsetbuf;
realtype bufleft[NUM_SPECIES*MYSUB], bufright[NUM_SPECIES*MYSUB];
/* If jysub > 0, send data from bottom x-line of cc. */
if (jysub != 0)
MPI_Send(&cdata[0], dsizex, PVEC_REAL_MPI_TYPE, my_pe-NPEX, 0, comm);
/* If jysub < NPEY-1, send data from top x-line of cc. */
if (jysub != NPEY-1) {
offsetc = (MYSUB-1)*dsizex;
MPI_Send(&cdata[offsetc], dsizex, PVEC_REAL_MPI_TYPE, my_pe+NPEX, 0, comm);
}
/* If ixsub > 0, send data from left y-line of cc (via bufleft). */
if (ixsub != 0) {
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NUM_SPECIES;
offsetc = ly*dsizex;
for (i = 0; i < NUM_SPECIES; i++)
bufleft[offsetbuf+i] = cdata[offsetc+i];
}
MPI_Send(&bufleft[0], dsizey, PVEC_REAL_MPI_TYPE, my_pe-1, 0, comm);
}
/* If ixsub < NPEX-1, send data from right y-line of cc (via bufright). */
if (ixsub != NPEX-1) {
for (ly = 0; ly < MYSUB; ly++) {
offsetbuf = ly*NUM_SPECIES;
offsetc = offsetbuf*MXSUB + (MXSUB-1)*NUM_SPECIES;
for (i = 0; i < NUM_SPECIES; i++)
bufright[offsetbuf+i] = cdata[offsetc+i];
}
MPI_Send(&bufright[0], dsizey, PVEC_REAL_MPI_TYPE, my_pe+1, 0, comm);
}
}
/* Define lines are for ease of readability in the following functions. */
#define mxsub (webdata->mxsub)
#define mysub (webdata->mysub)
#define npex (webdata->npex)
#define npey (webdata->npey)
#define ixsub (webdata->ixsub)
#define jysub (webdata->jysub)
#define nsmxsub (webdata->nsmxsub)
#define nsmxsub2 (webdata->nsmxsub2)
#define np (webdata->np)
#define dx (webdata->dx)
#define dy (webdata->dy)
#define cox (webdata->cox)
#define coy (webdata->coy)
#define rhs (webdata->rhs)
#define cext (webdata->cext)
#define rates (webdata->rates)
#define ns (webdata->ns)
#define acoef (webdata->acoef)
#define bcoef (webdata->bcoef)
/*
* reslocal: Compute res = F(t,cc,cp).
* This routine assumes that all inter-processor communication of data
* needed to calculate F has already been done. Components at interior
* subgrid boundaries are assumed to be in the work array cext.
* The local portion of the cc vector is first copied into cext.
* The exterior Neumann boundary conditions are explicitly handled here
* by copying data from the first interior mesh line to the ghost cell
* locations in cext. Then the reaction and diffusion terms are
* evaluated in terms of the cext array, and the residuals are formed.
* The reaction terms are saved separately in the vector webdata->rates
* for use by the preconditioner setup routine.
*/
static int reslocal(long int Nlocal, realtype tt,
N_Vector cc, N_Vector cp, N_Vector rr,
void *user_data)
{
realtype *cdata, *ratesxy, *cpxy, *resxy,
xx, yy, dcyli, dcyui, dcxli, dcxui;
int ix, jy, is, i, locc, ylocce, locce;
UserData webdata;
webdata = (UserData) user_data;
/* Get data pointers, subgrid data, array sizes, work array cext. */
cdata = NV_DATA_P(cc);
/* Copy local segment of cc vector into the working extended array cext. */
locc = 0;
locce = nsmxsub2 + NUM_SPECIES;
for (jy = 0; jy < mysub; jy++) {
for (i = 0; i < nsmxsub; i++) cext[locce+i] = cdata[locc+i];
locc = locc + nsmxsub;
locce = locce + nsmxsub2;
}
/* To facilitate homogeneous Neumann boundary conditions, when this is
a boundary PE, copy data from the first interior mesh line of cc to cext. */
/* If jysub = 0, copy x-line 2 of cc to cext. */
if (jysub == 0)
{ for (i = 0; i < nsmxsub; i++) cext[NUM_SPECIES+i] = cdata[nsmxsub+i]; }
/* If jysub = npey-1, copy x-line mysub-1 of cc to cext. */
if (jysub == npey-1) {
locc = (mysub-2)*nsmxsub;
locce = (mysub+1)*nsmxsub2 + NUM_SPECIES;
for (i = 0; i < nsmxsub; i++) cext[locce+i] = cdata[locc+i];
}
/* If ixsub = 0, copy y-line 2 of cc to cext. */
if (ixsub == 0) {
for (jy = 0; jy < mysub; jy++) {
locc = jy*nsmxsub + NUM_SPECIES;
locce = (jy+1)*nsmxsub2;
for (i = 0; i < NUM_SPECIES; i++) cext[locce+i] = cdata[locc+i];
}
}
/* If ixsub = npex-1, copy y-line mxsub-1 of cc to cext. */
if (ixsub == npex-1) {
for (jy = 0; jy < mysub; jy++) {
locc = (jy+1)*nsmxsub - 2*NUM_SPECIES;
locce = (jy+2)*nsmxsub2 - NUM_SPECIES;
for (i = 0; i < NUM_SPECIES; i++) cext[locce+i] = cdata[locc+i];
}
}
/* Loop over all grid points, setting local array rates to right-hand sides.
Then set rr values appropriately for prey/predator components of F. */
for (jy = 0; jy < mysub; jy++) {
ylocce = (jy+1)*nsmxsub2;
yy = (jy+jysub*mysub)*dy;
for (ix = 0; ix < mxsub; ix++) {
locce = ylocce + (ix+1)*NUM_SPECIES;
xx = (ix + ixsub*mxsub)*dx;
ratesxy = IJ_Vptr(rates,ix,jy);
WebRates(xx, yy, &(cext[locce]), ratesxy, webdata);
resxy = IJ_Vptr(rr,ix,jy);
cpxy = IJ_Vptr(cp,ix,jy);
for (is = 0; is < NUM_SPECIES; is++) {
dcyli = cext[locce+is] - cext[locce+is-nsmxsub2];
dcyui = cext[locce+is+nsmxsub2] - cext[locce+is];
dcxli = cext[locce+is] - cext[locce+is-NUM_SPECIES];
dcxui = cext[locce+is+NUM_SPECIES] - cext[locce+is];
rhs[is] = cox[is]*(dcxui-dcxli) + coy[is]*(dcyui-dcyli) + ratesxy[is];
if (is < np) resxy[is] = cpxy[is] - rhs[is];
else resxy[is] = - rhs[is];
}
}
}
return(0);
}
/*
* WebRates: Evaluate reaction rates at a given spatial point.
* At a given (x,y), evaluate the array of ns reaction terms R.
*/
static void WebRates(realtype xx, realtype yy, realtype *cxy, realtype *ratesxy,
UserData webdata)
{
int is;
realtype fac;
for (is = 0; is < NUM_SPECIES; is++)
ratesxy[is] = dotprod(NUM_SPECIES, cxy, acoef[is]);
fac = ONE + ALPHA*xx*yy + BETA*sin(FOURPI*xx)*sin(FOURPI*yy);
for (is = 0; is < NUM_SPECIES; is++)
ratesxy[is] = cxy[is]*( bcoef[is]*fac + ratesxy[is] );
}
/*
* dotprod: dot product routine for realtype arrays, for use by WebRates.
*/
static realtype dotprod(int size, realtype *x1, realtype *x2)
{
int i;
realtype *xx1, *xx2, temp = ZERO;
xx1 = x1;
xx2 = x2;
for (i = 0; i < size; i++)
temp += (*xx1++) * (*xx2++);
return(temp);
}
|