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
|
/*-----------------------------------------------------------------
* Programmer(s): Daniel R. Reynolds @ SMU
*---------------------------------------------------------------
* SUNDIALS Copyright Start
* Copyright (c) 2002-2022, Lawrence Livermore National Security
* and Southern Methodist University.
* All rights reserved.
*
* See the top-level LICENSE and NOTICE files for details.
*
* SPDX-License-Identifier: BSD-3-Clause
* SUNDIALS Copyright End
*---------------------------------------------------------------
* Example problem:
*
* The following test simulates a brusselator problem from chemical
* kinetics. This is an ODE system with 3 components, Y = [u,v,w],
* satisfying the equations,
* du/dt = a - (w+1)*u + v*u^2
* dv/dt = w*u - v*u^2
* dw/dt = (b-w)/ep - w*u
* for t in the interval [0.0, 10.0], with initial conditions
* Y0 = [u0,v0,w0].
*
* We have 3 different testing scenarios:
*
* Test 1: u0=3.9, v0=1.1, w0=2.8, a=1.2, b=2.5, ep=1.0e-5
* Here, all three components exhibit a rapid transient change
* during the first 0.2 time units, followed by a slow and
* smooth evolution.
*
* Test 2: u0=1.2, v0=3.1, w0=3, a=1, b=3.5, ep=5.0e-6
* Here, w experiences a fast initial transient, jumping 0.5
* within a few steps. All values proceed smoothly until
* around t=6.5, when both u and v undergo a sharp transition,
* with u increaseing from around 0.5 to 5 and v decreasing
* from around 6 to 1 in less than 0.5 time units. After this
* transition, both u and v continue to evolve somewhat
* rapidly for another 1.4 time units, and finish off smoothly.
*
* Test 3: u0=3, v0=3, w0=3.5, a=0.5, b=3, ep=5.0e-4
* Here, all components undergo very rapid initial transients
* during the first 0.3 time units, and all then proceed very
* smoothly for the remainder of the simulation.
*
* This file is hard-coded to use test 3.
*
* This program solves the problem with the ARK method, using an
* accelerated fixed-point iteration for the nonlinear solver.
*
* 100 outputs are printed at equal intervals, and run statistics
* are printed at the end.
*-----------------------------------------------------------------*/
/* Header files */
#include <stdio.h>
#include <math.h>
#include <arkode/arkode_arkstep.h> /* prototypes for ARKStep fcts., consts */
#include <nvector/nvector_serial.h> /* serial N_Vector types, fcts., macros */
#include <sunnonlinsol/sunnonlinsol_fixedpoint.h> /* access to FP nonlinear solver */
#include <sundials/sundials_types.h> /* def. of type 'realtype' */
#if defined(SUNDIALS_EXTENDED_PRECISION)
#define GSYM "Lg"
#define ESYM "Le"
#define FSYM "Lf"
#else
#define GSYM "g"
#define ESYM "e"
#define FSYM "f"
#endif
/* User-supplied Functions Called by the Solver */
static int fe(realtype t, N_Vector y, N_Vector ydot, void *user_data);
static int fi(realtype t, N_Vector y, N_Vector ydot, void *user_data);
/* Private function to check function return values */
static int check_flag(void *flagvalue, const char *funcname, int opt);
/* Main Program */
int main(int argc, char *argv[])
{
/* general problem parameters */
realtype T0 = RCONST(0.0); /* initial time */
realtype Tf = RCONST(10.0); /* final time */
realtype dTout = RCONST(1.0); /* time between outputs */
sunindextype NEQ = 3; /* number of dependent vars. */
int Nt = (int) ceil(Tf/dTout); /* number of output times */
int test = 3; /* test problem to run */
realtype reltol = 1.0e-6; /* tolerances */
realtype abstol = 1.0e-10;
int fp_m = 3; /* dimension of acceleration subspace */
int maxcor = 10; /* maximum # of nonlinear iterations/step */
realtype a, b, ep, u0, v0, w0;
realtype rdata[3];
int monitor = 0; /* turn on/off monitoring */
/* general problem variables */
int flag; /* reusable error-checking flag */
N_Vector y = NULL; /* empty vector for storing solution */
SUNNonlinearSolver NLS = NULL; /* empty nonlinear solver object */
void *arkode_mem = NULL; /* empty ARKode memory structure */
FILE *UFID, *INFOFID;
realtype t, tout;
int iout;
long int nst, nst_a, nfe, nfi, nni, ncfn, netf;
/* Create SUNDIALS context */
SUNContext ctx = NULL;
flag = SUNContext_Create(NULL, &ctx);
if (check_flag(&flag, "SUNContext_Create", 1)) return 1;
/* read inputs */
if (argc == 2) {
monitor = atoi(argv[1]);
}
/* set up the test problem according to the desired test */
if (test == 1) {
u0 = RCONST(3.9);
v0 = RCONST(1.1);
w0 = RCONST(2.8);
a = RCONST(1.2);
b = RCONST(2.5);
ep = RCONST(1.0e-5);
} else if (test == 3) {
u0 = RCONST(3.0);
v0 = RCONST(3.0);
w0 = RCONST(3.5);
a = RCONST(0.5);
b = RCONST(3.0);
ep = RCONST(5.0e-4);
} else {
u0 = RCONST(1.2);
v0 = RCONST(3.1);
w0 = RCONST(3.0);
a = RCONST(1.0);
b = RCONST(3.5);
ep = RCONST(5.0e-6);
}
/* Initial problem output */
printf("\nBrusselator ODE test problem, fixed-point solver:\n");
printf(" initial conditions: u0 = %"GSYM", v0 = %"GSYM", w0 = %"GSYM"\n",u0,v0,w0);
printf(" problem parameters: a = %"GSYM", b = %"GSYM", ep = %"GSYM"\n",a,b,ep);
printf(" reltol = %.1"ESYM", abstol = %.1"ESYM"\n\n",reltol,abstol);
/* Open up info output file */
INFOFID = NULL;
if (monitor) INFOFID = fopen("ark_brusselator_fp-info.txt","w");
/* Initialize data structures */
rdata[0] = a; /* set user data */
rdata[1] = b;
rdata[2] = ep;
y = N_VNew_Serial(NEQ, ctx); /* Create serial vector for solution */
if (check_flag((void *)y, "N_VNew_Serial", 0)) return 1;
NV_Ith_S(y,0) = u0; /* Set initial conditions */
NV_Ith_S(y,1) = v0;
NV_Ith_S(y,2) = w0;
/* Call ARKStepCreate to initialize the ARK timestepper module and
specify the right-hand side functions in y'=fe(t,y)+fi(t,y),
the inital time T0, and the initial dependent variable vector y. */
arkode_mem = ARKStepCreate(fe, fi, T0, y, ctx);
if (check_flag((void *)arkode_mem, "ARKStepCreate", 0)) return 1;
/* Initialize fixed-point nonlinear solver and attach to ARKStep */
NLS = SUNNonlinSol_FixedPoint(y, fp_m, ctx);
if (check_flag((void *)NLS, "SUNNonlinSol_FixedPoint", 0)) return 1;
if (monitor) {
flag = SUNNonlinSolSetPrintLevel_FixedPoint(NLS, 1);
if (check_flag(&flag, "SUNNonlinSolSetPrintLevel_Newton", 1)) return(1);
flag = SUNNonlinSolSetInfoFile_FixedPoint(NLS, INFOFID);
if (check_flag(&flag, "SUNNonlinSolSetPrintLevel_Newton", 1)) return(1);
}
flag = ARKStepSetNonlinearSolver(arkode_mem, NLS);
if (check_flag(&flag, "ARKStepSetNonlinearSolver", 1)) return 1;
/* Set routines */
flag = ARKStepSetUserData(arkode_mem, (void *) rdata); /* Pass rdata to user functions */
if (check_flag(&flag, "ARKStepSetUserData", 1)) return 1;
flag = ARKStepSStolerances(arkode_mem, reltol, abstol); /* Specify tolerances */
if (check_flag(&flag, "ARKStepSStolerances", 1)) return 1;
flag = ARKStepSetMaxNonlinIters(arkode_mem, maxcor); /* Increase default iterations */
if (check_flag(&flag, "ARKStepSetMaxNonlinIters", 1)) return 1;
/* Open output stream for results, output comment line */
UFID = fopen("solution.txt","w");
fprintf(UFID,"# t u v w\n");
/* output initial condition to disk */
fprintf(UFID," %.16"ESYM" %.16"ESYM" %.16"ESYM" %.16"ESYM"\n",
T0, NV_Ith_S(y,0), NV_Ith_S(y,1), NV_Ith_S(y,2));
/* Main time-stepping loop: calls ARKStepEvolve to perform the integration, then
prints results. Stops when the final time has been reached */
t = T0;
tout = T0+dTout;
printf(" t u v w\n");
printf(" ----------------------------------------------\n");
for (iout=0; iout<Nt; iout++) {
flag = ARKStepEvolve(arkode_mem, tout, y, &t, ARK_NORMAL); /* call integrator */
if (check_flag(&flag, "ARKStepEvolve", 1)) break;
printf(" %10.6"FSYM" %10.6"FSYM" %10.6"FSYM" %10.6"FSYM"\n", /* access/print solution */
t, NV_Ith_S(y,0), NV_Ith_S(y,1), NV_Ith_S(y,2));
fprintf(UFID," %.16"ESYM" %.16"ESYM" %.16"ESYM" %.16"ESYM"\n",
t, NV_Ith_S(y,0), NV_Ith_S(y,1), NV_Ith_S(y,2));
if (flag >= 0) { /* successful solve: update time */
tout += dTout;
tout = (tout > Tf) ? Tf : tout;
} else { /* unsuccessful solve: break */
fprintf(stderr,"Solver failure, stopping integration\n");
break;
}
}
printf(" ----------------------------------------------\n");
fclose(UFID);
if (monitor) fclose(INFOFID);
/* Print some final statistics */
flag = ARKStepGetNumSteps(arkode_mem, &nst);
check_flag(&flag, "ARKStepGetNumSteps", 1);
flag = ARKStepGetNumStepAttempts(arkode_mem, &nst_a);
check_flag(&flag, "ARKStepGetNumStepAttempts", 1);
flag = ARKStepGetNumRhsEvals(arkode_mem, &nfe, &nfi);
check_flag(&flag, "ARKStepGetNumRhsEvals", 1);
flag = ARKStepGetNumErrTestFails(arkode_mem, &netf);
check_flag(&flag, "ARKStepGetNumErrTestFails", 1);
flag = ARKStepGetNumNonlinSolvIters(arkode_mem, &nni);
check_flag(&flag, "ARKStepGetNumNonlinSolvIters", 1);
flag = ARKStepGetNumNonlinSolvConvFails(arkode_mem, &ncfn);
check_flag(&flag, "ARKStepGetNumNonlinSolvConvFails", 1);
printf("\nFinal Solver Statistics:\n");
printf(" Internal solver steps = %li (attempted = %li)\n", nst, nst_a);
printf(" Total RHS evals: Fe = %li, Fi = %li\n", nfe, nfi);
printf(" Total number of fixed-point iterations = %li\n", nni);
printf(" Total number of nonlinear solver convergence failures = %li\n", ncfn);
printf(" Total number of error test failures = %li\n\n", netf);
/* Clean up and return with successful completion */
N_VDestroy(y); /* Free y vector */
ARKStepFree(&arkode_mem); /* Free integrator memory */
SUNNonlinSolFree(NLS); /* Free NLS object */
SUNContext_Free(&ctx); /* Free context */
return 0;
}
/*-------------------------------
* Functions called by the solver
*-------------------------------*/
/* fi routine to compute the implicit portion of the ODE RHS. */
static int fi(realtype t, N_Vector y, N_Vector ydot, void *user_data)
{
realtype *rdata = (realtype *) user_data; /* cast user_data to realtype */
realtype b = rdata[1]; /* access data entries */
realtype ep = rdata[2];
realtype w = NV_Ith_S(y,2); /* access solution values */
/* fill in the RHS function */
NV_Ith_S(ydot,0) = 0.0;
NV_Ith_S(ydot,1) = 0.0;
NV_Ith_S(ydot,2) = (b-w)/ep;
return 0; /* Return with success */
}
/* fe routine to compute the explicit portion of the ODE RHS. */
static int fe(realtype t, N_Vector y, N_Vector ydot, void *user_data)
{
realtype *rdata = (realtype *) user_data; /* cast user_data to realtype */
realtype a = rdata[0]; /* access data entries */
realtype u = NV_Ith_S(y,0); /* access solution values */
realtype v = NV_Ith_S(y,1);
realtype w = NV_Ith_S(y,2);
/* fill in the RHS function */
NV_Ith_S(ydot,0) = a - (w+1.0)*u + v*u*u;
NV_Ith_S(ydot,1) = w*u - v*u*u;
NV_Ith_S(ydot,2) = -w*u;
return 0; /* Return with success */
}
/*-------------------------------
* Private helper functions
*-------------------------------*/
/* 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, const 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;
}
/*---- end of file ----*/
|