/usr/share/doc/libsundials-serial-dev/examples/ida/serial/idaHeat2D_bnd.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 | /*
* -----------------------------------------------------------------
* $Revision: 1.2 $
* $Date: 2009/09/30 23:25:59 $
* -----------------------------------------------------------------
* Programmer(s): Allan Taylor, Alan Hindmarsh and
* Radu Serban @ LLNL
* -----------------------------------------------------------------
* Example problem for IDA: 2D heat equation, serial, banded.
*
* This example solves a discretized 2D heat equation problem.
* This version uses the band solver IDABand, and IDACalcIC.
*
* The DAE system solved is a spatial discretization of the PDE
* du/dt = d^2u/dx^2 + d^2u/dy^2
* on the unit square. The boundary condition is u = 0 on all edges.
* Initial conditions are given by u = 16 x (1 - x) y (1 - y).
* The PDE is treated with central differences on a uniform M x M
* grid. The values of u at the interior points satisfy ODEs, and
* equations u = 0 at the boundaries are appended, to form a DAE
* system of size N = M^2. Here M = 10.
*
* The system is solved with IDA using the banded linear system
* solver, half-bandwidths equal to M, and default
* difference-quotient Jacobian. For purposes of illustration,
* IDACalcIC is called to compute correct values at the boundary,
* given incorrect values as input initial guesses. The constraints
* u >= 0 are posed for all components. Output is taken at
* t = 0, .01, .02, .04, ..., 10.24. (Output at t = 0 is for
* IDACalcIC cost statistics only.)
* -----------------------------------------------------------------
*/
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <ida/ida.h>
#include <ida/ida_band.h>
#include <nvector/nvector_serial.h>
#include <sundials/sundials_types.h>
/* Problem Constants */
#define NOUT 11
#define MGRID 10
#define NEQ MGRID*MGRID
#define ZERO RCONST(0.0)
#define ONE RCONST(1.0)
#define TWO RCONST(2.0)
#define BVAL RCONST(0.1)
/* Type: UserData */
typedef struct {
long int mm;
realtype dx;
realtype coeff;
} *UserData;
/* Prototypes of functions called by IDA */
int heatres(realtype tres, N_Vector uu, N_Vector up, N_Vector resval, void *user_data);
/* Prototypes of private functions */
static void PrintHeader(realtype rtol, realtype atol);
static void PrintOutput(void *mem, realtype t, N_Vector u);
static int SetInitialProfile(UserData data, N_Vector uu, N_Vector up,
N_Vector id, N_Vector res);
static int check_flag(void *flagvalue, char *funcname, int opt);
/*
*--------------------------------------------------------------------
* MAIN PROGRAM
*--------------------------------------------------------------------
*/
int main(void)
{
void *mem;
UserData data;
N_Vector uu, up, constraints, id, res;
int ier, iout;
long int mu, ml, netf, ncfn;
realtype rtol, atol, t0, t1, tout, tret;
mem = NULL;
data = NULL;
uu = up = constraints = id = res = NULL;
/* Create vectors uu, up, res, constraints, id. */
uu = N_VNew_Serial(NEQ);
if(check_flag((void *)uu, "N_VNew_Serial", 0)) return(1);
up = N_VNew_Serial(NEQ);
if(check_flag((void *)up, "N_VNew_Serial", 0)) return(1);
res = N_VNew_Serial(NEQ);
if(check_flag((void *)res, "N_VNew_Serial", 0)) return(1);
constraints = N_VNew_Serial(NEQ);
if(check_flag((void *)constraints, "N_VNew_Serial", 0)) return(1);
id = N_VNew_Serial(NEQ);
if(check_flag((void *)id, "N_VNew_Serial", 0)) return(1);
/* Create and load problem data block. */
data = (UserData) malloc(sizeof *data);
if(check_flag((void *)data, "malloc", 2)) return(1);
data->mm = MGRID;
data->dx = ONE/(MGRID - ONE);
data->coeff = ONE/( (data->dx) * (data->dx) );
/* Initialize uu, up, id. */
SetInitialProfile(data, uu, up, id, res);
/* Set constraints to all 1's for nonnegative solution values. */
N_VConst(ONE, constraints);
/* Set remaining input parameters. */
t0 = ZERO;
t1 = RCONST(0.01);
rtol = ZERO;
atol = RCONST(1.0e-3);
/* Call IDACreate and IDAMalloc to initialize solution */
mem = IDACreate();
if(check_flag((void *)mem, "IDACreate", 0)) return(1);
ier = IDASetUserData(mem, data);
if(check_flag(&ier, "IDASetUserData", 1)) return(1);
ier = IDASetId(mem, id);
if(check_flag(&ier, "IDASetId", 1)) return(1);
ier = IDASetConstraints(mem, constraints);
if(check_flag(&ier, "IDASetConstraints", 1)) return(1);
N_VDestroy_Serial(constraints);
ier = IDAInit(mem, heatres, t0, uu, up);
if(check_flag(&ier, "IDAInit", 1)) return(1);
ier = IDASStolerances(mem, rtol, atol);
if(check_flag(&ier, "IDASStolerances", 1)) return(1);
/* Call IDABand to specify the linear solver. */
mu = MGRID; ml = MGRID;
ier = IDABand(mem, NEQ, mu, ml);
if(check_flag(&ier, "IDABand", 1)) return(1);
/* Call IDACalcIC to correct the initial values. */
ier = IDACalcIC(mem, IDA_YA_YDP_INIT, t1);
if(check_flag(&ier, "IDACalcIC", 1)) return(1);
/* Print output heading. */
PrintHeader(rtol, atol);
PrintOutput(mem, t0, uu);
/* Loop over output times, call IDASolve, and print results. */
for (tout = t1, iout = 1; iout <= NOUT; iout++, tout *= TWO) {
ier = IDASolve(mem, tout, &tret, uu, up, IDA_NORMAL);
if(check_flag(&ier, "IDASolve", 1)) return(1);
PrintOutput(mem, tret, uu);
}
/* Print remaining counters and free memory. */
ier = IDAGetNumErrTestFails(mem, &netf);
check_flag(&ier, "IDAGetNumErrTestFails", 1);
ier = IDAGetNumNonlinSolvConvFails(mem, &ncfn);
check_flag(&ier, "IDAGetNumNonlinSolvConvFails", 1);
printf("\n netf = %ld, ncfn = %ld \n", netf, ncfn);
IDAFree(&mem);
N_VDestroy_Serial(uu);
N_VDestroy_Serial(up);
N_VDestroy_Serial(id);
N_VDestroy_Serial(res);
free(data);
return(0);
}
/*
*--------------------------------------------------------------------
* FUNCTIONS CALLED BY KINSOL
*--------------------------------------------------------------------
*/
/*
* heatres: heat equation system residual function
* This uses 5-point central differencing on the interior points, and
* includes algebraic equations for the boundary values.
* So for each interior point, the residual component has the form
* res_i = u'_i - (central difference)_i
* while for each boundary point, it is res_i = u_i.
*/
int heatres(realtype tres, N_Vector uu, N_Vector up, N_Vector resval,
void *user_data)
{
long int mm, i, j, offset, loc;
realtype *uv, *upv, *resv, coeff;
UserData data;
uv = NV_DATA_S(uu); upv = NV_DATA_S(up); resv = NV_DATA_S(resval);
data = (UserData)user_data;
mm = data->mm;
coeff = data->coeff;
/* Initialize resval to uu, to take care of boundary equations. */
N_VScale(ONE, uu, resval);
/* Loop over interior points; set res = up - (central difference). */
for (j = 1; j < mm-1; j++) {
offset = mm*j;
for (i = 1; i < mm-1; i++) {
loc = offset + i;
resv[loc] = upv[loc] - coeff *
(uv[loc-1] + uv[loc+1] + uv[loc-mm] + uv[loc+mm] - RCONST(4.0)*uv[loc]);
}
}
return(0);
}
/*
*--------------------------------------------------------------------
* PRIVATE FUNCTIONS
*--------------------------------------------------------------------
*/
/*
* SetInitialProfile: routine to initialize u, up, and id vectors.
*/
static int SetInitialProfile(UserData data, N_Vector uu, N_Vector up,
N_Vector id, N_Vector res)
{
realtype xfact, yfact, *udata, *updata, *iddata;
long int mm, mm1, i, j, offset, loc;
mm = data->mm;
mm1 = mm - 1;
udata = NV_DATA_S(uu);
updata = NV_DATA_S(up);
iddata = NV_DATA_S(id);
/* Initialize id to 1's. */
N_VConst(ONE, id);
/* Initialize uu on all grid points. */
for (j = 0; j < mm; j++) {
yfact = data->dx * j;
offset = mm*j;
for (i = 0;i < mm; i++) {
xfact = data->dx * i;
loc = offset + i;
udata[loc] = RCONST(16.0) * xfact * (ONE - xfact) * yfact * (ONE - yfact);
}
}
/* Initialize up vector to 0. */
N_VConst(ZERO, up);
/* heatres sets res to negative of ODE RHS values at interior points. */
heatres(ZERO, uu, up, res, data);
/* Copy -res into up to get correct interior initial up values. */
N_VScale(-ONE, res, up);
/* Finally, set values of u, up, and id at boundary points. */
for (j = 0; j < mm; j++) {
offset = mm*j;
for (i = 0;i < mm; i++) {
loc = offset + i;
if (j == 0 || j == mm1 || i == 0 || i == mm1 ) {
udata[loc] = BVAL; updata[loc] = ZERO; iddata[loc] = ZERO; }
}
}
return(0);
}
/*
* Print first lines of output (problem description)
*/
static void PrintHeader(realtype rtol, realtype atol)
{
printf("\nidaHeat2D_bnd: Heat equation, serial example problem for IDA\n");
printf(" Discretized heat equation on 2D unit square.\n");
printf(" Zero boundary conditions,");
printf(" polynomial initial conditions.\n");
printf(" Mesh dimensions: %d x %d", MGRID, MGRID);
printf(" Total system size: %d\n\n", NEQ);
#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("Constraints set to force all solution components >= 0. \n");
printf("Linear solver: IDABAND, banded direct solver \n");
printf(" difference quotient Jacobian, half-bandwidths = %d \n",MGRID);
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf("IDACalcIC called with input boundary values = %Lg \n",BVAL);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
printf("IDACalcIC called with input boundary values = %lg \n",BVAL);
#else
printf("IDACalcIC called with input boundary values = %g \n",BVAL);
#endif
/* Print output table heading and initial line of table. */
printf("\n Output Summary (umax = max-norm of solution) \n\n");
printf(" time umax k nst nni nje nre nreLS h \n" );
printf(" . . . . . . . . . . . . . . . . . . . . . \n");
}
/*
* Print Output
*/
static void PrintOutput(void *mem, realtype t, N_Vector uu)
{
int ier;
realtype umax, hused;
long int nst, nni, nje, nre, nreLS;
int kused;
umax = N_VMaxNorm(uu);
ier = IDAGetLastOrder(mem, &kused);
check_flag(&ier, "IDAGetLastOrder", 1);
ier = IDAGetNumSteps(mem, &nst);
check_flag(&ier, "IDAGetNumSteps", 1);
ier = IDAGetNumNonlinSolvIters(mem, &nni);
check_flag(&ier, "IDAGetNumNonlinSolvIters", 1);
ier = IDAGetNumResEvals(mem, &nre);
check_flag(&ier, "IDAGetNumResEvals", 1);
ier = IDAGetLastStep(mem, &hused);
check_flag(&ier, "IDAGetLastStep", 1);
ier = IDADlsGetNumJacEvals(mem, &nje);
check_flag(&ier, "IDADlsGetNumJacEvals", 1);
ier = IDADlsGetNumResEvals(mem, &nreLS);
check_flag(&ier, "IDADlsGetNumResEvals", 1);
#if defined(SUNDIALS_EXTENDED_PRECISION)
printf(" %5.2Lf %13.5Le %d %3ld %3ld %3ld %4ld %4ld %9.2Le \n",
t, umax, kused, nst, nni, nje, nre, nreLS, hused);
#elif defined(SUNDIALS_DOUBLE_PRECISION)
printf(" %5.2f %13.5le %d %3ld %3ld %3ld %4ld %4ld %9.2le \n",
t, umax, kused, nst, nni, nje, nre, nreLS, hused);
#else
printf(" %5.2f %13.5e %d %3ld %3ld %3ld %4ld %4ld %9.2e \n",
t, umax, kused, nst, nni, nje, nre, nreLS, hused);
#endif
}
/*
* 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);
} else if (opt == 1) {
/* Check if flag < 0 */
errflag = (int *) flagvalue;
if (*errflag < 0) {
fprintf(stderr,
"\nSUNDIALS_ERROR: %s() failed with flag = %d\n\n",
funcname, *errflag);
return(1);
}
} else if (opt == 2 && flagvalue == NULL) {
/* Check if function returned NULL pointer - no memory allocated */
fprintf(stderr,
"\nMEMORY_ERROR: %s() failed - returned NULL pointer\n\n",
funcname);
return(1);
}
return(0);
}
|