/usr/share/doc/libsundials-serial-dev/examples/cvodes/serial/cvsFoodWeb_ASAp_kry.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 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 | /*
* -----------------------------------------------------------------
* $Revision: 1.4 $
* $Date: 2011/11/23 23:53:02 $
* -----------------------------------------------------------------
* Programmer(s): Radu Serban @ LLNL
* -----------------------------------------------------------------
* This program solves a stiff ODE system that arises from a system
* of partial differential equations. 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 the following:
*
* 1 2 ns
* c = (c , c , ..., c )
*
* and the PDEs are as follows:
*
* i i i
* dc /dt = d(i)*(c + c ) + f (x,y,c) (i=1,...,ns)
* xx yy i
*
* where
*
* i ns j
* f (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) = -a (all i)
* a(i,j) = -g (i <= np, j > np)
* a(i,j) = e (i > np, j <= np)
* b(i) = b*(1 + alpha*x*y) (i <= np)
* b(i) = -b*(1 + alpha*x*y) (i > np)
* d(i) = Dprey (i <= np)
* d(i) = Dpred (i > np)
*
* The spatial domain is the unit square. The final time is 10.
* The boundary conditions are: normal derivative = 0.
* A polynomial in x and y is used to set the initial conditions.
*
* The PDEs are discretized by central differencing on an MX by
* MY mesh. The resulting ODE system is stiff.
*
* The ODE system is solved by CVODES using Newton iteration and
* the CVSPGMR linear solver (scaled preconditioned GMRES).
*
* The preconditioner matrix used is the product of two matrices:
* (1) A matrix, only defined implicitly, based on a fixed number
* of Gauss-Seidel iterations using the diffusion terms only.
* (2) A block-diagonal matrix based on the partial derivatives
* of the interaction terms f only, using block-grouping (computing
* only a subset of the ns by ns blocks).
*
* Additionally, CVODES can integrate backwards in time the
* the semi-discrete form of the adjoint PDE:
* d(lambda)/dt = - D^T ( lambda_xx + lambda_yy )
* - F_c^T lambda
* with homogeneous Neumann boundary conditions and final conditions
* lambda(x,y,t=t_final) = - g_c^T(t_final)
* whose solution at t = 0 represents the sensitivity of
* int_x int _y g(t_final,c) dx dy dt
* with respect to the initial conditions of the original problem.
*
* In this example,
* g(t,c) = c(ISPEC), with ISPEC defined below.
* -----------------------------------------------------------------
* Reference: Peter N. Brown and Alan C. Hindmarsh, Reduced Storage
* Matrix Methods in Stiff ODE Systems, J. Appl. Math. & Comp., 31
* (1989), pp. 40-91. Also available as Lawrence Livermore National
* Laboratory Report UCRL-95088, Rev. 1, June 1987.
* -----------------------------------------------------------------
*/
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <cvodes/cvodes.h>
#include <cvodes/cvodes_spgmr.h>
#include <nvector/nvector_serial.h>
#include <sundials/sundials_dense.h>
#include <sundials/sundials_types.h>
#include <sundials/sundials_math.h>
#define ZERO RCONST(0.0)
#define ONE RCONST(1.0)
#define TWO RCONST(2.0)
/* Problem Specification Constants */
#define AA ONE /* AA = a */
#define EE RCONST(1.0e4) /* EE = e */
#define GG RCONST(0.5e-6) /* GG = g */
#define BB ONE /* BB = b */
#define DPREY ONE
#define DPRED RCONST(0.5)
#define ALPH ONE
#define NP 3
#define NS (2*NP)
/* Method Constants */
#define MX 20
#define MY 20
#define MXNS (MX*NS)
#define AX ONE
#define AY ONE
#define DX (AX/(realtype)(MX-1))
#define DY (AY/(realtype)(MY-1))
#define MP NS
#define MQ (MX*MY)
#define MXMP (MX*MP)
#define NGX 2
#define NGY 2
#define NGRP (NGX*NGY)
#define ITMAX 5
/* CVodeInit Constants */
#define NEQ (NS*MX*MY)
#define T0 ZERO
#define RTOL RCONST(1.0e-5)
#define ATOL RCONST(1.0e-5)
/* Output Constants */
#define TOUT RCONST(10.0)
/* Note: The value for species i at mesh point (j,k) is stored in */
/* component number (i-1) + j*NS + k*NS*MX of an N_Vector, */
/* where 1 <= i <= NS, 0 <= j < MX, 0 <= k < MY. */
/* Structure for user data */
typedef struct {
realtype **P[NGRP];
long int *pivot[NGRP];
int ns, mxns, mp, mq, mx, my, ngrp, ngx, ngy, mxmp;
int jgx[NGX+1], jgy[NGY+1], jigx[MX], jigy[MY];
int jxr[NGX], jyr[NGY];
realtype acoef[NS][NS], bcoef[NS], diff[NS];
realtype cox[NS], coy[NS], dx, dy, srur;
realtype fsave[NEQ];
realtype fBsave[NEQ];
N_Vector rewt;
void *cvode_mem;
int indexB;
} *WebData;
/* Adjoint calculation constants */
/* g = int_x int_y c(ISPEC) dy dx at t = Tfinal */
#define NSTEPS 80 /* check points every NSTEPS steps */
#define ISPEC 6 /* species # in objective */
/* Prototypes for user-supplied functions */
static int f(realtype t, N_Vector y, N_Vector ydot, void *user_data);
static int Precond(realtype t, N_Vector c, N_Vector fc,
booleantype jok, booleantype *jcurPtr,
realtype gamma, void *user_data,
N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3);
static int PSolve(realtype t, N_Vector c, N_Vector fc,
N_Vector r, N_Vector z,
realtype gamma, realtype delta,
int lr, void *user_data, N_Vector vtemp);
static int fB(realtype t, N_Vector c, N_Vector cB,
N_Vector cBdot, void *user_data);
static int PrecondB(realtype t, N_Vector c,
N_Vector cB, N_Vector fcB, booleantype jok,
booleantype *jcurPtr, realtype gamma,
void *user_data,
N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3);
static int PSolveB(realtype t, N_Vector c,
N_Vector cB, N_Vector fcB,
N_Vector r, N_Vector z,
realtype gamma, realtype delta,
int lr, void *user_data, N_Vector vtemp);
/* Prototypes for private functions */
static WebData AllocUserData(void);
static void InitUserData(WebData wdata);
static void SetGroups(int m, int ng, int jg[], int jig[], int jr[]);
static void CInit(N_Vector c, WebData wdata);
static void CbInit(N_Vector c, int is, WebData wdata);
static void PrintOutput(N_Vector c, int ns, int mxns, WebData wdata);
static void FreeUserData(WebData wdata);
static void WebRates(realtype x, realtype y, realtype t, realtype c[], realtype rate[],
WebData wdata);
static void WebRatesB(realtype x, realtype y, realtype t, realtype c[], realtype cB[],
realtype rate[], realtype rateB[], WebData wdata);
static void fblock (realtype t, realtype cdata[], int jx, int jy, realtype cdotdata[],
WebData wdata);
static void GSIter(realtype gamma, N_Vector z, N_Vector x, WebData wdata);
static realtype doubleIntgr(N_Vector c, int i, WebData wdata);
static int check_flag(void *flagvalue, char *funcname, int opt);
/* Small Vector Kernels */
static void v_inc_by_prod(realtype u[], realtype v[], realtype w[], int n);
static void v_sum_prods(realtype u[], realtype p[], realtype q[], realtype v[],
realtype w[], int n);
static void v_prod(realtype u[], realtype v[], realtype w[], int n);
static void v_zero(realtype u[], int n);
/*
*--------------------------------------------------------------------
* MAIN PROGRAM
*--------------------------------------------------------------------
*/
int main(int argc, char *argv[])
{
realtype abstol=ATOL, reltol=RTOL, t;
N_Vector c;
WebData wdata;
void *cvode_mem;
int flag, ncheck;
int indexB;
realtype reltolB=RTOL, abstolB=ATOL;
N_Vector cB;
c = cB = NULL;
wdata = NULL;
cvode_mem = NULL;
/* Allocate and initialize user data */
wdata = AllocUserData();
if(check_flag((void *)wdata, "AllocUserData", 2)) return(1);
InitUserData(wdata);
/* Set-up forward problem */
/* Initializations */
c = N_VNew_Serial(NEQ);
if(check_flag((void *)c, "N_VNew_Serial", 0)) return(1);
CInit(c, wdata);
/* Call CVodeCreate/CVodeInit for forward run */
printf("\nCreate and allocate CVODES memory for forward run\n");
cvode_mem = CVodeCreate(CV_BDF, CV_NEWTON);
if(check_flag((void *)cvode_mem, "CVodeCreate", 0)) return(1);
wdata->cvode_mem = cvode_mem; /* Used in Precond */
flag = CVodeSetUserData(cvode_mem, wdata);
if(check_flag(&flag, "CVodeSetUserData", 1)) return(1);
flag = CVodeInit(cvode_mem, f, T0, c);
if(check_flag(&flag, "CVodeInit", 1)) return(1);
flag = CVodeSStolerances(cvode_mem, reltol, abstol);
if(check_flag(&flag, "CVodeSStolerances", 1)) return(1);
/* Call CVSpgmr for forward run */
flag = CVSpgmr(cvode_mem, PREC_LEFT, 0);
if(check_flag(&flag, "CVSpgmr", 1)) return(1);
flag = CVSpilsSetPreconditioner(cvode_mem, Precond, PSolve);
if(check_flag(&flag, "CVSpilsSetPreconditioner", 1)) return(1);
/* Set-up adjoint calculations */
printf("\nAllocate global memory\n");
flag = CVodeAdjInit(cvode_mem, NSTEPS, CV_HERMITE);
if(check_flag(&flag, "CVodeAdjInit", 1)) return(1);
/* Perform forward run */
printf("\nForward integration\n");
flag = CVodeF(cvode_mem, TOUT, c, &t, CV_NORMAL, &ncheck);
if(check_flag(&flag, "CVodeF", 1)) return(1);
printf("\nncheck = %d\n", ncheck);
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf("\n g = int_x int_y c%d(Tfinal,x,y) dx dy = %Lf \n\n",
ISPEC, doubleIntgr(c,ISPEC,wdata));
#else
printf("\n g = int_x int_y c%d(Tfinal,x,y) dx dy = %f \n\n",
ISPEC, doubleIntgr(c,ISPEC,wdata));
#endif
/* Set-up backward problem */
/* Allocate cB */
cB = N_VNew_Serial(NEQ);
if(check_flag((void *)cB, "N_VNew_Serial", 0)) return(1);
/* Initialize cB = 0 */
CbInit(cB, ISPEC, wdata);
/* Create and allocate CVODES memory for backward run */
printf("\nCreate and allocate CVODES memory for backward run\n");
flag = CVodeCreateB(cvode_mem, CV_BDF, CV_NEWTON, &indexB);
if(check_flag(&flag, "CVodeCreateB", 1)) return(1);
flag = CVodeSetUserDataB(cvode_mem, indexB, wdata);
if(check_flag(&flag, "CVodeSetUserDataB", 1)) return(1);
flag = CVodeInitB(cvode_mem, indexB, fB, TOUT, cB);
if(check_flag(&flag, "CVodeInitB", 1)) return(1);
flag = CVodeSStolerancesB(cvode_mem, indexB, reltolB, abstolB);
if(check_flag(&flag, "CVodeSStolerancesB", 1)) return(1);
wdata->indexB = indexB;
/* Call CVSpgmr */
flag = CVSpgmrB(cvode_mem, indexB, PREC_LEFT, 0);
if(check_flag(&flag, "CVSpgmrB", 1)) return(1);
flag = CVSpilsSetPreconditionerB(cvode_mem, indexB, PrecondB, PSolveB);
if(check_flag(&flag, "CVSpilsSetPreconditionerB", 1)) return(1);
/* Perform backward integration */
printf("\nBackward integration\n");
flag = CVodeB(cvode_mem, T0, CV_NORMAL);
if(check_flag(&flag, "CVodeB", 1)) return(1);
flag = CVodeGetB(cvode_mem, indexB, &t, cB);
if(check_flag(&flag, "CVodeGetB", 1)) return(1);
PrintOutput(cB, NS, MXNS, wdata);
/* Free all memory */
CVodeFree(&cvode_mem);
N_VDestroy_Serial(c);
N_VDestroy_Serial(cB);
FreeUserData(wdata);
return(0);
}
/*
*--------------------------------------------------------------------
* FUNCTIONS CALLED BY CVODES
*--------------------------------------------------------------------
*/
/*
* This routine computes the right-hand side of the ODE system and
* returns it in cdot. The interaction rates are computed by calls to WebRates,
* and these are saved in fsave for use in preconditioning.
*/
static int f(realtype t, N_Vector c, N_Vector cdot, void *user_data)
{
int i, ic, ici, idxl, idxu, idyl, idyu, iyoff, jx, jy, ns, mxns;
realtype dcxli, dcxui, dcyli, dcyui, x, y, *cox, *coy, *fsave, dx, dy;
realtype *cdata, *cdotdata;
WebData wdata;
wdata = (WebData) user_data;
cdata = NV_DATA_S(c);
cdotdata = NV_DATA_S(cdot);
mxns = wdata->mxns;
ns = wdata->ns;
fsave = wdata->fsave;
cox = wdata->cox;
coy = wdata->coy;
mxns = wdata->mxns;
dx = wdata->dx;
dy = wdata->dy;
for (jy = 0; jy < MY; jy++) {
y = jy*dy;
iyoff = mxns*jy;
idyu = (jy == MY-1) ? -mxns : mxns;
idyl = (jy == 0) ? -mxns : mxns;
for (jx = 0; jx < MX; jx++) {
x = jx*dx;
ic = iyoff + ns*jx;
/* Get interaction rates at one point (x,y). */
WebRates(x, y, t, cdata+ic, fsave+ic, wdata);
idxu = (jx == MX-1) ? -ns : ns;
idxl = (jx == 0) ? -ns : ns;
for (i = 1; i <= ns; i++) {
ici = ic + i-1;
/* Do differencing in y. */
dcyli = cdata[ici] - cdata[ici-idyl];
dcyui = cdata[ici+idyu] - cdata[ici];
/* Do differencing in x. */
dcxli = cdata[ici] - cdata[ici-idxl];
dcxui = cdata[ici+idxu] - cdata[ici];
/* Collect terms and load cdot elements. */
cdotdata[ici] = coy[i-1]*(dcyui - dcyli) +
cox[i-1]*(dcxui - dcxli) +
fsave[ici];
}
}
}
return(0);
}
/*
* This routine generates the block-diagonal part of the Jacobian
* corresponding to the interaction rates, multiplies by -gamma, adds
* the identity matrix, and calls denseGETRF to do the LU decomposition of
* each diagonal block. The computation of the diagonal blocks uses
* the preset block and grouping information. One block per group is
* computed. The Jacobian elements are generated by difference
* quotients using calls to the routine fblock.
*
* This routine can be regarded as a prototype for the general case
* of a block-diagonal preconditioner. The blocks are of size mp, and
* there are ngrp=ngx*ngy blocks computed in the block-grouping scheme.
*/
static int Precond(realtype t, N_Vector c, N_Vector fc,
booleantype jok, booleantype *jcurPtr,
realtype gamma, void *user_data,
N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
{
realtype ***P;
long int **pivot, ier;
int i, if0, if00, ig, igx, igy, j, jj, jx, jy;
int *jxr, *jyr, ngrp, ngx, ngy, mxmp, flag;
long int mp;
realtype uround, fac, r, r0, save, srur;
realtype *f1, *fsave, *cdata, *rewtdata;
WebData wdata;
N_Vector rewt;
wdata = (WebData) user_data;
rewt = wdata->rewt;
flag = CVodeGetErrWeights(wdata->cvode_mem, rewt);
if(check_flag(&flag, "CVodeGetErrWeights", 1)) return(1);
cdata = NV_DATA_S(c);
rewtdata = NV_DATA_S(rewt);
uround = UNIT_ROUNDOFF;
P = wdata->P;
pivot = wdata->pivot;
jxr = wdata->jxr;
jyr = wdata->jyr;
mp = wdata->mp;
srur = wdata->srur;
ngrp = wdata->ngrp;
ngx = wdata->ngx;
ngy = wdata->ngy;
mxmp = wdata->mxmp;
fsave = wdata->fsave;
/* Make mp calls to fblock to approximate each diagonal block of Jacobian.
Here, fsave contains the base value of the rate vector and
r0 is a minimum increment factor for the difference quotient. */
f1 = NV_DATA_S(vtemp1);
fac = N_VWrmsNorm (fc, rewt);
r0 = RCONST(1000.0)*ABS(gamma)*uround*NEQ*fac;
if (r0 == ZERO) r0 = ONE;
for (igy = 0; igy < ngy; igy++) {
jy = jyr[igy];
if00 = jy*mxmp;
for (igx = 0; igx < ngx; igx++) {
jx = jxr[igx];
if0 = if00 + jx*mp;
ig = igx + igy*ngx;
/* Generate ig-th diagonal block */
for (j = 0; j < mp; j++) {
/* Generate the jth column as a difference quotient */
jj = if0 + j;
save = cdata[jj];
r = MAX(srur*ABS(save),r0/rewtdata[jj]);
cdata[jj] += r;
fac = -gamma/r;
fblock (t, cdata, jx, jy, f1, wdata);
for (i = 0; i < mp; i++) {
P[ig][j][i] = (f1[i] - fsave[if0+i])*fac;
}
cdata[jj] = save;
}
}
}
/* Add identity matrix and do LU decompositions on blocks. */
for (ig = 0; ig < ngrp; ig++) {
denseAddIdentity(P[ig], mp);
ier = denseGETRF(P[ig], mp, mp, pivot[ig]);
if (ier != 0) return(1);
}
*jcurPtr = TRUE;
return(0);
}
/*
* This routine applies two inverse preconditioner matrices
* to the vector r, using the interaction-only block-diagonal Jacobian
* with block-grouping, denoted Jr, and Gauss-Seidel applied to the
* diffusion contribution to the Jacobian, denoted Jd.
* It first calls GSIter for a Gauss-Seidel approximation to
* ((I - gamma*Jd)-inverse)*r, and stores the result in z.
* Then it computes ((I - gamma*Jr)-inverse)*z, using LU factors of the
* blocks in P, and pivot information in pivot, and returns the result in z.
*/
static int PSolve(realtype t, N_Vector c, N_Vector fc,
N_Vector r, N_Vector z,
realtype gamma, realtype delta,
int lr, void *user_data, N_Vector vtemp)
{
realtype ***P;
long int **pivot;
int jx, jy, igx, igy, iv, ig, *jigx, *jigy, mx, my, ngx, mp;
WebData wdata;
wdata = (WebData) user_data;
N_VScale(ONE, r, z);
/* call GSIter for Gauss-Seidel iterations */
GSIter(gamma, z, vtemp, wdata);
/* Do backsolves for inverse of block-diagonal preconditioner factor */
P = wdata->P;
pivot = wdata->pivot;
mx = wdata->mx;
my = wdata->my;
ngx = wdata->ngx;
mp = wdata->mp;
jigx = wdata->jigx;
jigy = wdata->jigy;
iv = 0;
for (jy = 0; jy < my; jy++) {
igy = jigy[jy];
for (jx = 0; jx < mx; jx++) {
igx = jigx[jx];
ig = igx + igy*ngx;
denseGETRS(P[ig], mp, pivot[ig], &(NV_DATA_S(z)[iv]));
iv += mp;
}
}
return(0);
}
/*
* This routine computes the right-hand side of the adjoint ODE system and
* returns it in cBdot. The interaction rates are computed by calls to WebRates,
* and these are saved in fsave for use in preconditioning. The adjoint
* interaction rates are computed by calls to WebRatesB.
*/
static int fB(realtype t, N_Vector c, N_Vector cB,
N_Vector cBdot, void *user_data)
{
int i, ic, ici, idxl, idxu, idyl, idyu, iyoff, jx, jy, ns, mxns;
realtype dcxli, dcxui, dcyli, dcyui, x, y, *cox, *coy, *fsave, *fBsave, dx, dy;
realtype *cdata, *cBdata, *cBdotdata;
WebData wdata;
wdata = (WebData) user_data;
cdata = NV_DATA_S(c);
cBdata = NV_DATA_S(cB);
cBdotdata = NV_DATA_S(cBdot);
mxns = wdata->mxns;
ns = wdata->ns;
fsave = wdata->fsave;
fBsave = wdata->fBsave;
cox = wdata->cox;
coy = wdata->coy;
mxns = wdata->mxns;
dx = wdata->dx;
dy = wdata->dy;
for (jy = 0; jy < MY; jy++) {
y = jy*dy;
iyoff = mxns*jy;
idyu = (jy == MY-1) ? -mxns : mxns;
idyl = (jy == 0) ? -mxns : mxns;
for (jx = 0; jx < MX; jx++) {
x = jx*dx;
ic = iyoff + ns*jx;
/* Get interaction rates at one point (x,y). */
WebRatesB(x, y, t, cdata+ic, cBdata+ic, fsave+ic, fBsave+ic, wdata);
idxu = (jx == MX-1) ? -ns : ns;
idxl = (jx == 0) ? -ns : ns;
for (i = 1; i <= ns; i++) {
ici = ic + i-1;
/* Do differencing in y. */
dcyli = cBdata[ici] - cBdata[ici-idyl];
dcyui = cBdata[ici+idyu] - cBdata[ici];
/* Do differencing in x. */
dcxli = cBdata[ici] - cBdata[ici-idxl];
dcxui = cBdata[ici+idxu] - cBdata[ici];
/* Collect terms and load cdot elements. */
cBdotdata[ici] = - coy[i-1]*(dcyui - dcyli)
- cox[i-1]*(dcxui - dcxli)
- fBsave[ici];
}
}
}
return(0);
}
/*
* Preconditioner setup function for the backward problem
*/
static int PrecondB(realtype t, N_Vector c,
N_Vector cB, N_Vector fcB, booleantype jok,
booleantype *jcurPtr, realtype gamma,
void *user_data,
N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
{
realtype ***P;
long int **pivot, ier;
int i, if0, if00, ig, igx, igy, j, jj, jx, jy;
int *jxr, *jyr, mp, ngrp, ngx, ngy, mxmp, flag;
realtype uround, fac, r, r0, save, srur;
realtype *f1, *fsave, *cdata, *rewtdata;
void *cvode_mem;
WebData wdata;
N_Vector rewt;
wdata = (WebData) user_data;
cvode_mem = CVodeGetAdjCVodeBmem(wdata->cvode_mem, wdata->indexB);
if(check_flag((void *)cvode_mem, "CVadjGetCVodeBmem", 0)) return(1);
rewt = wdata->rewt;
flag = CVodeGetErrWeights(cvode_mem, rewt);
if(check_flag(&flag, "CVodeGetErrWeights", 1)) return(1);
cdata = NV_DATA_S(c);
rewtdata = NV_DATA_S(rewt);
uround = UNIT_ROUNDOFF;
P = wdata->P;
pivot = wdata->pivot;
jxr = wdata->jxr;
jyr = wdata->jyr;
mp = wdata->mp;
srur = wdata->srur;
ngrp = wdata->ngrp;
ngx = wdata->ngx;
ngy = wdata->ngy;
mxmp = wdata->mxmp;
fsave = wdata->fsave;
/* Make mp calls to fblock to approximate each diagonal block of Jacobian.
Here, fsave contains the base value of the rate vector and
r0 is a minimum increment factor for the difference quotient. */
f1 = NV_DATA_S(vtemp1);
fac = N_VWrmsNorm (fcB, rewt);
r0 = RCONST(1000.0)*ABS(gamma)*uround*NEQ*fac;
if (r0 == ZERO) r0 = ONE;
for (igy = 0; igy < ngy; igy++) {
jy = jyr[igy];
if00 = jy*mxmp;
for (igx = 0; igx < ngx; igx++) {
jx = jxr[igx];
if0 = if00 + jx*mp;
ig = igx + igy*ngx;
/* Generate ig-th diagonal block */
for (j = 0; j < mp; j++) {
/* Generate the jth column as a difference quotient */
jj = if0 + j;
save = cdata[jj];
r = MAX(srur*ABS(save),r0/rewtdata[jj]);
cdata[jj] += r;
fac = gamma/r;
fblock (t, cdata, jx, jy, f1, wdata);
for (i = 0; i < mp; i++) {
P[ig][i][j] = (f1[i] - fsave[if0+i])*fac;
}
cdata[jj] = save;
}
}
}
/* Add identity matrix and do LU decompositions on blocks. */
for (ig = 0; ig < ngrp; ig++) {
denseAddIdentity(P[ig], mp);
ier = denseGETRF(P[ig], mp, mp, pivot[ig]);
if (ier != 0) return(1);
}
*jcurPtr = TRUE;
return(0);
}
/*
* Preconditioner solve function for the backward problem
*/
static int PSolveB(realtype t, N_Vector c,
N_Vector cB, N_Vector fcB,
N_Vector r, N_Vector z,
realtype gamma, realtype delta,
int lr, void *user_data, N_Vector vtemp)
{
realtype ***P;
long int **pivot;
int jx, jy, igx, igy, iv, ig, *jigx, *jigy, mx, my, ngx;
long int mp;
WebData wdata;
wdata = (WebData) user_data;
N_VScale(ONE, r, z);
/* call GSIter for Gauss-Seidel iterations (same routine but with gamma=-gamma) */
GSIter(-gamma, z, vtemp, wdata);
/* Do backsolves for inverse of block-diagonal preconditioner factor */
P = wdata->P;
pivot = wdata->pivot;
mx = wdata->mx;
my = wdata->my;
ngx = wdata->ngx;
mp = wdata->mp;
jigx = wdata->jigx;
jigy = wdata->jigy;
iv = 0;
for (jy = 0; jy < my; jy++) {
igy = jigy[jy];
for (jx = 0; jx < mx; jx++) {
igx = jigx[jx];
ig = igx + igy*ngx;
denseGETRS(P[ig], mp, pivot[ig], &(NV_DATA_S(z)[iv]));
iv += mp;
}
}
return(0);
}
/*
*--------------------------------------------------------------------
* PRIVATE FUNCTIONS
*--------------------------------------------------------------------
*/
/*
* Allocate space for user data structure
*/
static WebData AllocUserData(void)
{
int i, ngrp = NGRP;
long int ns = NS;
WebData wdata;
wdata = (WebData) malloc(sizeof *wdata);
for(i=0; i < ngrp; i++) {
(wdata->P)[i] = newDenseMat(ns, ns);
(wdata->pivot)[i] = newLintArray(ns);
}
wdata->rewt = N_VNew_Serial(NEQ);
return(wdata);
}
/*
* Initialize user data structure
*/
static void InitUserData(WebData wdata)
{
int i, j, ns;
realtype *bcoef, *diff, *cox, *coy, dx, dy;
realtype (*acoef)[NS];
acoef = wdata->acoef;
bcoef = wdata->bcoef;
diff = wdata->diff;
cox = wdata->cox;
coy = wdata->coy;
ns = wdata->ns = NS;
for (j = 0; j < NS; j++) { for (i = 0; i < NS; i++) acoef[i][j] = ZERO; }
for (j = 0; j < NP; j++) {
for (i = 0; i < NP; i++) {
acoef[NP+i][j] = EE;
acoef[i][NP+j] = -GG;
}
acoef[j][j] = -AA;
acoef[NP+j][NP+j] = -AA;
bcoef[j] = BB;
bcoef[NP+j] = -BB;
diff[j] = DPREY;
diff[NP+j] = DPRED;
}
/* Set remaining problem parameters */
wdata->mxns = MXNS;
dx = wdata->dx = DX;
dy = wdata->dy = DY;
for (i = 0; i < ns; i++) {
cox[i] = diff[i]/SQR(dx);
coy[i] = diff[i]/SQR(dy);
}
/* Set remaining method parameters */
wdata->mp = MP;
wdata->mq = MQ;
wdata->mx = MX;
wdata->my = MY;
wdata->srur = SQRT(UNIT_ROUNDOFF);
wdata->mxmp = MXMP;
wdata->ngrp = NGRP;
wdata->ngx = NGX;
wdata->ngy = NGY;
SetGroups(MX, NGX, wdata->jgx, wdata->jigx, wdata->jxr);
SetGroups(MY, NGY, wdata->jgy, wdata->jigy, wdata->jyr);
}
/*
* This routine sets arrays jg, jig, and jr describing
* a uniform partition of (0,1,2,...,m-1) into ng groups.
* The arrays set are:
* jg = length ng+1 array of group boundaries.
* Group ig has indices j = jg[ig],...,jg[ig+1]-1.
* jig = length m array of group indices vs node index.
* Node index j is in group jig[j].
* jr = length ng array of indices representing the groups.
* The index for group ig is j = jr[ig].
*/
static void SetGroups(int m, int ng, int jg[], int jig[], int jr[])
{
int ig, j, len1, mper, ngm1;
mper = m/ng; /* does integer division */
for (ig=0; ig < ng; ig++) jg[ig] = ig*mper;
jg[ng] = m;
ngm1 = ng - 1;
len1 = ngm1*mper;
for (j = 0; j < len1; j++) jig[j] = j/mper;
for (j = len1; j < m; j++) jig[j] = ngm1;
for (ig = 0; ig < ngm1; ig++) jr[ig] = ((2*ig+1)*mper-1)/2;
jr[ngm1] = (ngm1*mper+m-1)/2;
}
/*
* This routine computes and loads the vector of initial values.
*/
static void CInit(N_Vector c, WebData wdata)
{
int i, ici, ioff, iyoff, jx, jy, ns, mxns;
realtype argx, argy, x, y, dx, dy, x_factor, y_factor, *cdata;
cdata = NV_DATA_S(c);
ns = wdata->ns;
mxns = wdata->mxns;
dx = wdata->dx;
dy = wdata->dy;
x_factor = RCONST(4.0)/SQR(AX);
y_factor = RCONST(4.0)/SQR(AY);
for (jy = 0; jy < MY; jy++) {
y = jy*dy;
argy = SQR(y_factor*y*(AY-y));
iyoff = mxns*jy;
for (jx = 0; jx < MX; jx++) {
x = jx*dx;
argx = SQR(x_factor*x*(AX-x));
ioff = iyoff + ns*jx;
for (i = 1; i <= ns; i++) {
ici = ioff + i-1;
cdata[ici] = RCONST(10.0) + i*argx*argy;
}
}
}
}
/*
* This function computes and loads the final values for the adjoint variables
*/
static void CbInit(N_Vector c, int is, WebData wdata)
{
int i, ici, ioff, iyoff, jx, jy, ns, mxns;
realtype *cdata;
realtype gu[NS];
cdata = NV_DATA_S(c);
ns = wdata->ns;
mxns = wdata->mxns;
for ( i = 1; i <= ns; i++ ) gu[i-1] = ZERO;
gu[ISPEC-1] = ONE;
for (jy = 0; jy < MY; jy++) {
iyoff = mxns*jy;
for (jx = 0; jx < MX; jx++) {
ioff = iyoff + ns*jx;
for (i = 1; i <= ns; i++) {
ici = ioff + i-1;
cdata[ici] = gu[i-1];
}
}
}
}
/*
* This routine computes the interaction rates for the species
* c_1, ... ,c_ns (stored in c[0],...,c[ns-1]), at one spatial point
* and at time t.
*/
static void WebRates(realtype x, realtype y, realtype t, realtype c[],
realtype rate[], WebData wdata)
{
int i, j, ns;
realtype fac, *bcoef;
realtype (*acoef)[NS];
ns = wdata->ns;
acoef = wdata->acoef;
bcoef = wdata->bcoef;
for (i = 0; i < ns; i++)
rate[i] = ZERO;
for (j = 0; j < ns; j++)
for (i = 0; i < ns; i++)
rate[i] += c[j] * acoef[i][j];
fac = ONE + ALPH*x*y;
for (i = 0; i < ns; i++)
rate[i] = c[i]*(bcoef[i]*fac + rate[i]);
}
/*
* This routine computes the interaction rates for the backward problem
*/
static void WebRatesB(realtype x, realtype y, realtype t, realtype c[], realtype cB[],
realtype rate[], realtype rateB[], WebData wdata)
{
int i, j, ns;
realtype fac, *bcoef;
realtype (*acoef)[NS];
ns = wdata->ns;
acoef = wdata->acoef;
bcoef = wdata->bcoef;
fac = ONE + ALPH*x*y;
for (i = 0; i < ns; i++)
rate[i] = bcoef[i]*fac;
for (j = 0; j < ns; j++)
for (i = 0; i < ns; i++)
rate[i] += acoef[i][j]*c[j];
for (i = 0; i < ns; i++) {
rateB[i] = cB[i]*rate[i];
rate[i] = c[i]*rate[i];
}
for (j = 0; j < ns; j++)
for (i = 0; i < ns; i++)
rateB[i] += acoef[j][i]*c[j]*cB[j];
}
/*
* This routine computes one block of the interaction terms of the
* system, namely block (jx,jy), for use in preconditioning.
* Here jx and jy count from 0.
*/
static void fblock(realtype t, realtype cdata[], int jx, int jy,
realtype cdotdata[], WebData wdata)
{
int iblok, ic;
realtype x, y;
iblok = jx + jy*(wdata->mx);
y = jy*(wdata->dy);
x = jx*(wdata->dx);
ic = (wdata->ns)*(iblok);
WebRates(x, y, t, cdata+ic, cdotdata, wdata);
}
/*
* This routine performs ITMAX=5 Gauss-Seidel iterations to compute an
* approximation to (P-inverse)*z, where P = I - gamma*Jd, and
* Jd represents the diffusion contributions to the Jacobian.
* The answer is stored in z on return, and x is a temporary vector.
* The dimensions below assume a global constant NS >= ns.
* Some inner loops of length ns are implemented with the small
* vector kernels v_sum_prods, v_prod, v_inc_by_prod.
*/
static void GSIter(realtype gamma, N_Vector z, N_Vector x, WebData wdata)
{
int i, ic, iter, iyoff, jx, jy, ns, mxns, mx, my, x_loc, y_loc;
realtype beta[NS], beta2[NS], cof1[NS], gam[NS], gam2[NS];
realtype temp, *cox, *coy, *xd, *zd;
xd = NV_DATA_S(x);
zd = NV_DATA_S(z);
ns = wdata->ns;
mx = wdata->mx;
my = wdata->my;
mxns = wdata->mxns;
cox = wdata->cox;
coy = wdata->coy;
/* Write matrix as P = D - L - U.
Load local arrays beta, beta2, gam, gam2, and cof1. */
for (i = 0; i < ns; i++) {
temp = ONE/(ONE + TWO*gamma*(cox[i] + coy[i]));
beta[i] = gamma*cox[i]*temp;
beta2[i] = TWO*beta[i];
gam[i] = gamma*coy[i]*temp;
gam2[i] = TWO*gam[i];
cof1[i] = temp;
}
/* Begin iteration loop.
Load vector x with (D-inverse)*z for first iteration. */
for (jy = 0; jy < my; jy++) {
iyoff = mxns*jy;
for (jx = 0; jx < mx; jx++) {
ic = iyoff + ns*jx;
v_prod(xd+ic, cof1, zd+ic, ns); /* x[ic+i] = cof1[i]z[ic+i] */
}
}
N_VConst(ZERO, z);
/* Looping point for iterations. */
for (iter=1; iter <= ITMAX; iter++) {
/* Calculate (D-inverse)*U*x if not the first iteration. */
if (iter > 1) {
for (jy=0; jy < my; jy++) {
iyoff = mxns*jy;
for (jx=0; jx < mx; jx++) { /* order of loops matters */
ic = iyoff + ns*jx;
x_loc = (jx == 0) ? 0 : ((jx == mx-1) ? 2 : 1);
y_loc = (jy == 0) ? 0 : ((jy == my-1) ? 2 : 1);
switch (3*y_loc+x_loc) {
case 0 : /* jx == 0, jy == 0 */
/* x[ic+i] = beta2[i]x[ic+ns+i] + gam2[i]x[ic+mxns+i] */
v_sum_prods(xd+ic, beta2, xd+ic+ns, gam2, xd+ic+mxns, ns);
break;
case 1 : /* 1 <= jx <= mx-2, jy == 0 */
/* x[ic+i] = beta[i]x[ic+ns+i] + gam2[i]x[ic+mxns+i] */
v_sum_prods(xd+ic, beta, xd+ic+ns, gam2, xd+ic+mxns, ns);
break;
case 2 : /* jx == mx-1, jy == 0 */
/* x[ic+i] = gam2[i]x[ic+mxns+i] */
v_prod(xd+ic, gam2, xd+ic+mxns, ns);
break;
case 3 : /* jx == 0, 1 <= jy <= my-2 */
/* x[ic+i] = beta2[i]x[ic+ns+i] + gam[i]x[ic+mxns+i] */
v_sum_prods(xd+ic, beta2, xd+ic+ns, gam, xd+ic+mxns, ns);
break;
case 4 : /* 1 <= jx <= mx-2, 1 <= jy <= my-2 */
/* x[ic+i] = beta[i]x[ic+ns+i] + gam[i]x[ic+mxns+i] */
v_sum_prods(xd+ic, beta, xd+ic+ns, gam, xd+ic+mxns, ns);
break;
case 5 : /* jx == mx-1, 1 <= jy <= my-2 */
/* x[ic+i] = gam[i]x[ic+mxns+i] */
v_prod(xd+ic, gam, xd+ic+mxns, ns);
break;
case 6 : /* jx == 0, jy == my-1 */
/* x[ic+i] = beta2[i]x[ic+ns+i] */
v_prod(xd+ic, beta2, xd+ic+ns, ns);
break;
case 7 : /* 1 <= jx <= mx-2, jy == my-1 */
/* x[ic+i] = beta[i]x[ic+ns+i] */
v_prod(xd+ic, beta, xd+ic+ns, ns);
break;
case 8 : /* jx == mx-1, jy == my-1 */
/* x[ic+i] = ZERO */
v_zero(xd+ic, ns);
break;
}
}
}
} /* end if (iter > 1) */
/* Overwrite x with [(I - (D-inverse)*L)-inverse]*x. */
for (jy=0; jy < my; jy++) {
iyoff = mxns*jy;
for (jx=0; jx < mx; jx++) { /* order of loops matters */
ic = iyoff + ns*jx;
x_loc = (jx == 0) ? 0 : ((jx == mx-1) ? 2 : 1);
y_loc = (jy == 0) ? 0 : ((jy == my-1) ? 2 : 1);
switch (3*y_loc+x_loc) {
case 0 : /* jx == 0, jy == 0 */
break;
case 1 : /* 1 <= jx <= mx-2, jy == 0 */
/* x[ic+i] += beta[i]x[ic-ns+i] */
v_inc_by_prod(xd+ic, beta, xd+ic-ns, ns);
break;
case 2 : /* jx == mx-1, jy == 0 */
/* x[ic+i] += beta2[i]x[ic-ns+i] */
v_inc_by_prod(xd+ic, beta2, xd+ic-ns, ns);
break;
case 3 : /* jx == 0, 1 <= jy <= my-2 */
/* x[ic+i] += gam[i]x[ic-mxns+i] */
v_inc_by_prod(xd+ic, gam, xd+ic-mxns, ns);
break;
case 4 : /* 1 <= jx <= mx-2, 1 <= jy <= my-2 */
/* x[ic+i] += beta[i]x[ic-ns+i] + gam[i]x[ic-mxns+i] */
v_inc_by_prod(xd+ic, beta, xd+ic-ns, ns);
v_inc_by_prod(xd+ic, gam, xd+ic-mxns, ns);
break;
case 5 : /* jx == mx-1, 1 <= jy <= my-2 */
/* x[ic+i] += beta2[i]x[ic-ns+i] + gam[i]x[ic-mxns+i] */
v_inc_by_prod(xd+ic, beta2, xd+ic-ns, ns);
v_inc_by_prod(xd+ic, gam, xd+ic-mxns, ns);
break;
case 6 : /* jx == 0, jy == my-1 */
/* x[ic+i] += gam2[i]x[ic-mxns+i] */
v_inc_by_prod(xd+ic, gam2, xd+ic-mxns, ns);
break;
case 7 : /* 1 <= jx <= mx-2, jy == my-1 */
/* x[ic+i] += beta[i]x[ic-ns+i] + gam2[i]x[ic-mxns+i] */
v_inc_by_prod(xd+ic, beta, xd+ic-ns, ns);
v_inc_by_prod(xd+ic, gam2, xd+ic-mxns, ns);
break;
case 8 : /* jx == mx-1, jy == my-1 */
/* x[ic+i] += beta2[i]x[ic-ns+i] + gam2[i]x[ic-mxns+i] */
v_inc_by_prod(xd+ic, beta2, xd+ic-ns, ns);
v_inc_by_prod(xd+ic, gam2, xd+ic-mxns, ns);
break;
}
}
}
/* Add increment x to z : z <- z+x */
N_VLinearSum(ONE, z, ONE, x, z);
}
}
static void v_inc_by_prod(realtype u[], realtype v[], realtype w[], int n)
{
int i;
for (i=0; i < n; i++) u[i] += v[i]*w[i];
}
static void v_sum_prods(realtype u[], realtype p[], realtype q[],
realtype v[], realtype w[], int n)
{
int i;
for (i=0; i < n; i++) u[i] = p[i]*q[i] + v[i]*w[i];
}
static void v_prod(realtype u[], realtype v[], realtype w[], int n)
{
int i;
for (i=0; i < n; i++) u[i] = v[i]*w[i];
}
static void v_zero(realtype u[], int n)
{
int i;
for (i=0; i < n; i++) u[i] = ZERO;
}
/*
* Print maximum sensitivity of G for each species
*/
static void PrintOutput(N_Vector cB, int ns, int mxns, WebData wdata)
{
int i, jx, jy;
realtype *cdata, cij, cmax, x, y;
x = y = ZERO;
cdata = NV_DATA_S(cB);
for (i=1; i <= ns; i++) {
cmax = ZERO;
for (jy=MY-1; jy >= 0; jy--) {
for (jx=0; jx < MX; jx++) {
cij = cdata[(i-1) + jx*ns + jy*mxns];
if (ABS(cij) > cmax) {
cmax = cij;
x = jx * wdata->dx;
y = jy * wdata->dy;
}
}
}
printf("\nMaximum sensitivity with respect to I.C. of species %d\n", i);
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf(" mu max = %Le\n",cmax);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
printf(" mu max = %le\n",cmax);
#else
printf(" mu max = %e\n",cmax);
#endif
printf("at\n");
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf(" x = %Le\n y = %Le\n", x, y);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
printf(" x = %le\n y = %le\n", x, y);
#else
printf(" x = %e\n y = %e\n", x, y);
#endif
}
}
/*
* Compute double space integral
*/
static realtype doubleIntgr(N_Vector c, int i, WebData wdata)
{
realtype *cdata;
int ns, mx, my, mxns;
realtype dx, dy;
realtype intgr_xy, intgr_x;
int jx, jy;
cdata = NV_DATA_S(c);
ns = wdata->ns;
mx = wdata->mx;
my = wdata->my;
mxns = wdata->mxns;
dx = wdata->dx;
dy = wdata->dy;
jy = 0;
intgr_x = cdata[(i-1)+jy*mxns];
for (jx = 1; jx < mx-1; jx++) {
intgr_x += TWO*cdata[(i-1) + jx*ns + jy*mxns];
}
intgr_x += cdata[(i-1)+(mx-1)*ns+jy*mxns];
intgr_x *= RCONST(0.5)*dx;
intgr_xy = intgr_x;
for (jy = 1; jy < my-1; jy++) {
intgr_x = cdata[(i-1)+jy*mxns];
for (jx = 1; jx < mx-1; jx++) {
intgr_x += TWO*cdata[(i-1) + jx*ns + jy*mxns];
}
intgr_x += cdata[(i-1)+(mx-1)*ns+jy*mxns];
intgr_x *= RCONST(0.5)*dx;
intgr_xy += TWO*intgr_x;
}
jy = my-1;
intgr_x = cdata[(i-1)+jy*mxns];
for (jx = 1; jx < mx-1; jx++) {
intgr_x += TWO*cdata[(i-1) + jx*ns + jy*mxns];
}
intgr_x += cdata[(i-1)+(mx-1)*ns+jy*mxns];
intgr_x *= RCONST(0.5)*dx;
intgr_xy += intgr_x;
intgr_xy *= RCONST(0.5)*dy;
return(intgr_xy);
}
/*
* Free space allocated for the user data structure
*/
static void FreeUserData(WebData wdata)
{
int i, ngrp;
ngrp = wdata->ngrp;
for(i=0; i < ngrp; i++) {
destroyMat((wdata->P)[i]);
destroyArray((wdata->pivot)[i]);
}
N_VDestroy_Serial(wdata->rewt);
free(wdata);
}
/*
* 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 *errflag;
/* Check if SUNDIALS function returned NULL pointer - no memory allocated */
if (opt == 0 && flagvalue == NULL) {
fprintf(stderr, "\nSUNDIALS_ERROR: %s() failed - returned NULL pointer\n\n",
funcname);
return(1); }
/* Check if flag < 0 */
else if (opt == 1) {
errflag = (int *) flagvalue;
if (*errflag < 0) {
fprintf(stderr, "\nSUNDIALS_ERROR: %s() failed with flag = %d\n\n",
funcname, *errflag);
return(1); }}
/* Check if function returned NULL pointer - no memory allocated */
else if (opt == 2 && flagvalue == NULL) {
fprintf(stderr, "\nMEMORY_ERROR: %s() failed - returned NULL pointer\n\n",
funcname);
return(1); }
return(0);
}
|