1 /* 2 Note: 3 To performance a convergence study, a reference solution can be generated with a small stepsize and stored into a binary file, e.g. -ts_type ssp -ts_dt 1e-5 -ts_max_steps 30000 -ts_view_solution binary:reference.bin 4 Errors can be computed in the following runs with -simulation -f reference.bin 5 6 Multirate RK options: 7 -ts_rk_dtratio is the ratio between larger time step size and small time step size 8 -ts_rk_multirate_type has three choices: none (for single step RK) 9 combined (for multirate method and user just need to provide the same RHS with the single step RK) 10 partitioned (for multiraet method and user need to provide two RHS, one is for fast components and the orther is for slow components 11 */ 12 13 static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" 14 " advection - Variable coefficient scalar advection\n" 15 " u_t + (a*u)_x = 0\n" 16 " a is a piecewise function with two values say a0 and a1 (both a0 and a1 are constant).\n" 17 " in this toy problem we assume a=a0 when x<0 and a=a1 when x>0 with a0<a1 which has the same discontinuous point with initial condition.\n" 18 " we don't provide an exact solution, so you need to generate reference solution to do convergnece staudy,\n" 19 " more precisely, you need firstly generate a reference to a binary file say file.bin, then on commend line,\n" 20 " you should type -simulation -f file.bin.\n" 21 " you can choose the number of grids by -da_grid_x.\n" 22 " you can choose the value of a by -physics_advect_a1 and -physics_advect_a2.\n"; 23 24 #include <petscts.h> 25 #include <petscdm.h> 26 #include <petscdmda.h> 27 #include <petscdraw.h> 28 #include <petsc/private/tsimpl.h> 29 30 #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */ 31 32 #include "finitevolume1d.h" 33 34 static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) { 35 PetscReal range = xmax - xmin; 36 return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range); 37 } 38 39 /* --------------------------------- Advection ----------------------------------- */ 40 41 typedef struct { 42 PetscReal a[2]; /* advective velocity */ 43 } AdvectCtx; 44 45 static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed, PetscReal x, PetscReal xmin, PetscReal xmax) { 46 AdvectCtx *ctx = (AdvectCtx *)vctx; 47 PetscReal *speed; 48 49 PetscFunctionBeginUser; 50 speed = ctx->a; 51 if (x == 0 || x == xmin || x == xmax) 52 flux[0] = PetscMax(0, (speed[0] + speed[1]) / 2.0) * uL[0] + PetscMin(0, (speed[0] + speed[1]) / 2.0) * uR[0]; /* if condition need to be changed base on different problem, '0' is the discontinuous point of a */ 53 else if (x < 0) flux[0] = PetscMax(0, speed[0]) * uL[0] + PetscMin(0, speed[0]) * uR[0]; /* else if condition need to be changed based on diferent problem, 'x = 0' is discontinuous point of a */ 54 else flux[0] = PetscMax(0, speed[1]) * uL[0] + PetscMin(0, speed[1]) * uR[0]; 55 *maxspeed = *speed; 56 PetscFunctionReturn(0); 57 } 58 59 static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds, PetscReal x) { 60 AdvectCtx *ctx = (AdvectCtx *)vctx; 61 62 PetscFunctionBeginUser; 63 X[0] = 1.; 64 Xi[0] = 1.; 65 if (x < 0) speeds[0] = ctx->a[0]; /* x is discontinuous point of a */ 66 else speeds[0] = ctx->a[1]; 67 PetscFunctionReturn(0); 68 } 69 70 static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) { 71 AdvectCtx *ctx = (AdvectCtx *)vctx; 72 PetscReal *a = ctx->a, x0; 73 74 PetscFunctionBeginUser; 75 if (x < 0) { /* x is cell center */ 76 switch (bctype) { 77 case FVBC_OUTFLOW: x0 = x - a[0] * t; break; 78 case FVBC_PERIODIC: x0 = RangeMod(x - a[0] * t, xmin, xmax); break; 79 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType"); 80 } 81 switch (initial) { 82 case 0: u[0] = (x0 < 0) ? 1 : -1; break; 83 case 1: u[0] = (x0 < 0) ? -1 : 1; break; 84 case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 85 case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break; 86 case 4: u[0] = PetscAbs(x0); break; 87 case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break; 88 case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break; 89 case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break; 90 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition"); 91 } 92 } else { 93 switch (bctype) { 94 case FVBC_OUTFLOW: x0 = x - a[1] * t; break; 95 case FVBC_PERIODIC: x0 = RangeMod(x - a[1] * t, xmin, xmax); break; 96 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType"); 97 } 98 switch (initial) { 99 case 0: u[0] = (x0 < 0) ? 1 : -1; break; 100 case 1: u[0] = (x0 < 0) ? -1 : 1; break; 101 case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 102 case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break; 103 case 4: u[0] = PetscAbs(x0); break; 104 case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break; 105 case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break; 106 case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break; 107 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition"); 108 } 109 } 110 PetscFunctionReturn(0); 111 } 112 113 static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) { 114 AdvectCtx *user; 115 116 PetscFunctionBeginUser; 117 PetscCall(PetscNew(&user)); 118 ctx->physics.sample = PhysicsSample_Advect; 119 ctx->physics.riemann = PhysicsRiemann_Advect; 120 ctx->physics.characteristic = PhysicsCharacteristic_Advect; 121 ctx->physics.destroy = PhysicsDestroy_SimpleFree; 122 ctx->physics.user = user; 123 ctx->physics.dof = 1; 124 PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0])); 125 user->a[0] = 1; 126 user->a[1] = 1; 127 PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", ""); 128 { 129 PetscCall(PetscOptionsReal("-physics_advect_a1", "Speed1", "", user->a[0], &user->a[0], NULL)); 130 PetscCall(PetscOptionsReal("-physics_advect_a2", "Speed2", "", user->a[1], &user->a[1], NULL)); 131 } 132 PetscOptionsEnd(); 133 PetscFunctionReturn(0); 134 } 135 136 /* --------------------------------- Finite Volume Solver for slow parts ----------------------------------- */ 137 138 PetscErrorCode FVRHSFunctionslow(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 139 FVCtx *ctx = (FVCtx *)vctx; 140 PetscInt i, j, k, Mx, dof, xs, xm, len_slow; 141 PetscReal hx, cfl_idt = 0; 142 PetscScalar *x, *f, *slope; 143 Vec Xloc; 144 DM da; 145 146 PetscFunctionBeginUser; 147 PetscCall(TSGetDM(ts, &da)); 148 PetscCall(DMGetLocalVector(da, &Xloc)); 149 PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 150 hx = (ctx->xmax - ctx->xmin) / Mx; 151 PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 152 PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 153 PetscCall(VecZeroEntries(F)); 154 PetscCall(DMDAVecGetArray(da, Xloc, &x)); 155 PetscCall(VecGetArray(F, &f)); 156 PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 157 PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 158 PetscCall(ISGetSize(ctx->iss, &len_slow)); 159 160 if (ctx->bctype == FVBC_OUTFLOW) { 161 for (i = xs - 2; i < 0; i++) { 162 for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 163 } 164 for (i = Mx; i < xs + xm + 2; i++) { 165 for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 166 } 167 } 168 for (i = xs - 1; i < xs + xm + 1; i++) { 169 struct _LimitInfo info; 170 PetscScalar *cjmpL, *cjmpR; 171 if (i < len_slow + 1) { 172 /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 173 PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 174 /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 175 PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 176 cjmpL = &ctx->cjmpLR[0]; 177 cjmpR = &ctx->cjmpLR[dof]; 178 for (j = 0; j < dof; j++) { 179 PetscScalar jmpL, jmpR; 180 jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 181 jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 182 for (k = 0; k < dof; k++) { 183 cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 184 cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 185 } 186 } 187 /* Apply limiter to the left and right characteristic jumps */ 188 info.m = dof; 189 info.hx = hx; 190 (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 191 for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 192 for (j = 0; j < dof; j++) { 193 PetscScalar tmp = 0; 194 for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 195 slope[i * dof + j] = tmp; 196 } 197 } 198 } 199 200 for (i = xs; i < xs + xm + 1; i++) { 201 PetscReal maxspeed; 202 PetscScalar *uL, *uR; 203 if (i < len_slow) { /* slow parts can be changed based on a */ 204 uL = &ctx->uLR[0]; 205 uR = &ctx->uLR[dof]; 206 for (j = 0; j < dof; j++) { 207 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 208 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 209 } 210 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 211 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 212 if (i > xs) { 213 for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx; 214 } 215 if (i < xs + xm) { 216 for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hx; 217 } 218 } 219 if (i == len_slow) { /* slow parts can be changed based on a */ 220 uL = &ctx->uLR[0]; 221 uR = &ctx->uLR[dof]; 222 for (j = 0; j < dof; j++) { 223 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 224 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 225 } 226 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 227 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 228 if (i > xs) { 229 for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx; 230 } 231 } 232 } 233 PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 234 PetscCall(VecRestoreArray(F, &f)); 235 PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 236 PetscCall(DMRestoreLocalVector(da, &Xloc)); 237 PetscFunctionReturn(0); 238 } 239 240 /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 241 242 PetscErrorCode FVRHSFunctionfast(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 243 FVCtx *ctx = (FVCtx *)vctx; 244 PetscInt i, j, k, Mx, dof, xs, xm, len_slow; 245 PetscReal hx, cfl_idt = 0; 246 PetscScalar *x, *f, *slope; 247 Vec Xloc; 248 DM da; 249 250 PetscFunctionBeginUser; 251 PetscCall(TSGetDM(ts, &da)); 252 PetscCall(DMGetLocalVector(da, &Xloc)); 253 PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 254 hx = (ctx->xmax - ctx->xmin) / Mx; 255 PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 256 PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 257 PetscCall(VecZeroEntries(F)); 258 PetscCall(DMDAVecGetArray(da, Xloc, &x)); 259 PetscCall(VecGetArray(F, &f)); 260 PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 261 PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 262 PetscCall(ISGetSize(ctx->iss, &len_slow)); 263 264 if (ctx->bctype == FVBC_OUTFLOW) { 265 for (i = xs - 2; i < 0; i++) { 266 for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 267 } 268 for (i = Mx; i < xs + xm + 2; i++) { 269 for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 270 } 271 } 272 for (i = xs - 1; i < xs + xm + 1; i++) { 273 struct _LimitInfo info; 274 PetscScalar *cjmpL, *cjmpR; 275 if (i > len_slow - 2) { 276 PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 277 PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 278 cjmpL = &ctx->cjmpLR[0]; 279 cjmpR = &ctx->cjmpLR[dof]; 280 for (j = 0; j < dof; j++) { 281 PetscScalar jmpL, jmpR; 282 jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 283 jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 284 for (k = 0; k < dof; k++) { 285 cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 286 cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 287 } 288 } 289 /* Apply limiter to the left and right characteristic jumps */ 290 info.m = dof; 291 info.hx = hx; 292 (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 293 for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 294 for (j = 0; j < dof; j++) { 295 PetscScalar tmp = 0; 296 for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 297 slope[i * dof + j] = tmp; 298 } 299 } 300 } 301 302 for (i = xs; i < xs + xm + 1; i++) { 303 PetscReal maxspeed; 304 PetscScalar *uL, *uR; 305 if (i > len_slow) { 306 uL = &ctx->uLR[0]; 307 uR = &ctx->uLR[dof]; 308 for (j = 0; j < dof; j++) { 309 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 310 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 311 } 312 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 313 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 314 if (i > xs) { 315 for (j = 0; j < dof; j++) f[(i - len_slow - 1) * dof + j] -= ctx->flux[j] / hx; 316 } 317 if (i < xs + xm) { 318 for (j = 0; j < dof; j++) f[(i - len_slow) * dof + j] += ctx->flux[j] / hx; 319 } 320 } 321 if (i == len_slow) { 322 uL = &ctx->uLR[0]; 323 uR = &ctx->uLR[dof]; 324 for (j = 0; j < dof; j++) { 325 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 326 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 327 } 328 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 329 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 330 if (i < xs + xm) { 331 for (j = 0; j < dof; j++) f[(i - len_slow) * dof + j] += ctx->flux[j] / hx; 332 } 333 } 334 } 335 PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 336 PetscCall(VecRestoreArray(F, &f)); 337 PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 338 PetscCall(DMRestoreLocalVector(da, &Xloc)); 339 PetscFunctionReturn(0); 340 } 341 342 PetscErrorCode FVRHSFunctionslow2(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 343 FVCtx *ctx = (FVCtx *)vctx; 344 PetscInt i, j, k, Mx, dof, xs, xm, len_slow1, len_slow2; 345 PetscReal hx, cfl_idt = 0; 346 PetscScalar *x, *f, *slope; 347 Vec Xloc; 348 DM da; 349 350 PetscFunctionBeginUser; 351 PetscCall(TSGetDM(ts, &da)); 352 PetscCall(DMGetLocalVector(da, &Xloc)); 353 PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 354 hx = (ctx->xmax - ctx->xmin) / Mx; 355 PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 356 PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 357 PetscCall(VecZeroEntries(F)); 358 PetscCall(DMDAVecGetArray(da, Xloc, &x)); 359 PetscCall(VecGetArray(F, &f)); 360 PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 361 PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 362 PetscCall(ISGetSize(ctx->iss, &len_slow1)); 363 PetscCall(ISGetSize(ctx->iss2, &len_slow2)); 364 if (ctx->bctype == FVBC_OUTFLOW) { 365 for (i = xs - 2; i < 0; i++) { 366 for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 367 } 368 for (i = Mx; i < xs + xm + 2; i++) { 369 for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 370 } 371 } 372 for (i = xs - 1; i < xs + xm + 1; i++) { 373 struct _LimitInfo info; 374 PetscScalar *cjmpL, *cjmpR; 375 if (i < len_slow1 + len_slow2 + 1 && i > len_slow1 - 2) { 376 PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 377 PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 378 cjmpL = &ctx->cjmpLR[0]; 379 cjmpR = &ctx->cjmpLR[dof]; 380 for (j = 0; j < dof; j++) { 381 PetscScalar jmpL, jmpR; 382 jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 383 jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 384 for (k = 0; k < dof; k++) { 385 cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 386 cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 387 } 388 } 389 /* Apply limiter to the left and right characteristic jumps */ 390 info.m = dof; 391 info.hx = hx; 392 (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 393 for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 394 for (j = 0; j < dof; j++) { 395 PetscScalar tmp = 0; 396 for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 397 slope[i * dof + j] = tmp; 398 } 399 } 400 } 401 402 for (i = xs; i < xs + xm + 1; i++) { 403 PetscReal maxspeed; 404 PetscScalar *uL, *uR; 405 if (i < len_slow1 + len_slow2 && i > len_slow1) { /* slow parts can be changed based on a */ 406 uL = &ctx->uLR[0]; 407 uR = &ctx->uLR[dof]; 408 for (j = 0; j < dof; j++) { 409 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 410 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 411 } 412 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 413 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 414 if (i > xs) { 415 for (j = 0; j < dof; j++) f[(i - len_slow1 - 1) * dof + j] -= ctx->flux[j] / hx; 416 } 417 if (i < xs + xm) { 418 for (j = 0; j < dof; j++) f[(i - len_slow1) * dof + j] += ctx->flux[j] / hx; 419 } 420 } 421 if (i == len_slow1 + len_slow2) { /* slow parts can be changed based on a */ 422 uL = &ctx->uLR[0]; 423 uR = &ctx->uLR[dof]; 424 for (j = 0; j < dof; j++) { 425 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 426 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 427 } 428 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 429 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 430 if (i > xs) { 431 for (j = 0; j < dof; j++) f[(i - len_slow1 - 1) * dof + j] -= ctx->flux[j] / hx; 432 } 433 } 434 if (i == len_slow1) { /* slow parts can be changed based on a */ 435 uL = &ctx->uLR[0]; 436 uR = &ctx->uLR[dof]; 437 for (j = 0; j < dof; j++) { 438 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 439 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 440 } 441 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 442 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 443 if (i < xs + xm) { 444 for (j = 0; j < dof; j++) f[(i - len_slow1) * dof + j] += ctx->flux[j] / hx; 445 } 446 } 447 } 448 449 PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 450 PetscCall(VecRestoreArray(F, &f)); 451 PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 452 PetscCall(DMRestoreLocalVector(da, &Xloc)); 453 PetscFunctionReturn(0); 454 } 455 456 PetscErrorCode FVRHSFunctionfast2(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 457 FVCtx *ctx = (FVCtx *)vctx; 458 PetscInt i, j, k, Mx, dof, xs, xm, len_slow1, len_slow2; 459 PetscReal hx, cfl_idt = 0; 460 PetscScalar *x, *f, *slope; 461 Vec Xloc; 462 DM da; 463 464 PetscFunctionBeginUser; 465 PetscCall(TSGetDM(ts, &da)); 466 PetscCall(DMGetLocalVector(da, &Xloc)); 467 PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 468 hx = (ctx->xmax - ctx->xmin) / Mx; 469 PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 470 PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 471 PetscCall(VecZeroEntries(F)); 472 PetscCall(DMDAVecGetArray(da, Xloc, &x)); 473 PetscCall(VecGetArray(F, &f)); 474 PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 475 PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 476 PetscCall(ISGetSize(ctx->iss, &len_slow1)); 477 PetscCall(ISGetSize(ctx->iss2, &len_slow2)); 478 479 if (ctx->bctype == FVBC_OUTFLOW) { 480 for (i = xs - 2; i < 0; i++) { 481 for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 482 } 483 for (i = Mx; i < xs + xm + 2; i++) { 484 for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 485 } 486 } 487 for (i = xs - 1; i < xs + xm + 1; i++) { 488 struct _LimitInfo info; 489 PetscScalar *cjmpL, *cjmpR; 490 if (i > len_slow1 + len_slow2 - 2) { 491 PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 492 PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 493 cjmpL = &ctx->cjmpLR[0]; 494 cjmpR = &ctx->cjmpLR[dof]; 495 for (j = 0; j < dof; j++) { 496 PetscScalar jmpL, jmpR; 497 jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 498 jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 499 for (k = 0; k < dof; k++) { 500 cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 501 cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 502 } 503 } 504 /* Apply limiter to the left and right characteristic jumps */ 505 info.m = dof; 506 info.hx = hx; 507 (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 508 for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 509 for (j = 0; j < dof; j++) { 510 PetscScalar tmp = 0; 511 for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 512 slope[i * dof + j] = tmp; 513 } 514 } 515 } 516 517 for (i = xs; i < xs + xm + 1; i++) { 518 PetscReal maxspeed; 519 PetscScalar *uL, *uR; 520 if (i > len_slow1 + len_slow2) { 521 uL = &ctx->uLR[0]; 522 uR = &ctx->uLR[dof]; 523 for (j = 0; j < dof; j++) { 524 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 525 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 526 } 527 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 528 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 529 if (i > xs) { 530 for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2 - 1) * dof + j] -= ctx->flux[j] / hx; 531 } 532 if (i < xs + xm) { 533 for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2) * dof + j] += ctx->flux[j] / hx; 534 } 535 } 536 if (i == len_slow1 + len_slow2) { 537 uL = &ctx->uLR[0]; 538 uR = &ctx->uLR[dof]; 539 for (j = 0; j < dof; j++) { 540 uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 541 uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 542 } 543 PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 544 cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 545 if (i < xs + xm) { 546 for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2) * dof + j] += ctx->flux[j] / hx; 547 } 548 } 549 } 550 PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 551 PetscCall(VecRestoreArray(F, &f)); 552 PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 553 PetscCall(DMRestoreLocalVector(da, &Xloc)); 554 PetscFunctionReturn(0); 555 } 556 557 int main(int argc, char *argv[]) { 558 char lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m"; 559 PetscFunctionList limiters = 0, physics = 0; 560 MPI_Comm comm; 561 TS ts; 562 DM da; 563 Vec X, X0, R; 564 FVCtx ctx; 565 PetscInt i, k, dof, xs, xm, Mx, draw = 0, *index_slow, *index_fast, islow = 0, ifast = 0; 566 PetscBool view_final = PETSC_FALSE; 567 PetscReal ptime; 568 569 PetscFunctionBeginUser; 570 PetscCall(PetscInitialize(&argc, &argv, 0, help)); 571 comm = PETSC_COMM_WORLD; 572 PetscCall(PetscMemzero(&ctx, sizeof(ctx))); 573 574 /* Register limiters to be available on the command line */ 575 PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit_Upwind)); 576 PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit_LaxWendroff)); 577 PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit_BeamWarming)); 578 PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit_Fromm)); 579 PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit_Minmod)); 580 PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit_Superbee)); 581 PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit_MC)); 582 PetscCall(PetscFunctionListAdd(&limiters, "vanleer", Limit_VanLeer)); 583 PetscCall(PetscFunctionListAdd(&limiters, "vanalbada", Limit_VanAlbada)); 584 PetscCall(PetscFunctionListAdd(&limiters, "vanalbadatvd", Limit_VanAlbadaTVD)); 585 PetscCall(PetscFunctionListAdd(&limiters, "koren", Limit_Koren)); 586 PetscCall(PetscFunctionListAdd(&limiters, "korensym", Limit_KorenSym)); 587 PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit_Koren3)); 588 PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon2", Limit_CadaTorrilhon2)); 589 PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r0p1", Limit_CadaTorrilhon3R0p1)); 590 PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r1", Limit_CadaTorrilhon3R1)); 591 PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r10", Limit_CadaTorrilhon3R10)); 592 PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r100", Limit_CadaTorrilhon3R100)); 593 594 /* Register physical models to be available on the command line */ 595 PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect)); 596 597 ctx.comm = comm; 598 ctx.cfl = 0.9; 599 ctx.bctype = FVBC_PERIODIC; 600 ctx.xmin = -1.0; 601 ctx.xmax = 1.0; 602 PetscOptionsBegin(comm, NULL, "Finite Volume solver options", ""); 603 PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL)); 604 PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL)); 605 PetscCall(PetscOptionsFList("-limit", "Name of flux limiter to use", "", limiters, lname, lname, sizeof(lname), NULL)); 606 PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL)); 607 PetscCall(PetscOptionsString("-view_final", "Write final solution in ASCII MATLAB format to given file name", "", final_fname, final_fname, sizeof(final_fname), &view_final)); 608 PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL)); 609 PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL)); 610 PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL)); 611 PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL)); 612 PetscCall(PetscOptionsBool("-recursive_split", "Split the domain recursively", "", ctx.recursive, &ctx.recursive, NULL)); 613 PetscOptionsEnd(); 614 615 /* Choose the limiter from the list of registered limiters */ 616 PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit)); 617 PetscCheck(ctx.limit, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname); 618 619 /* Choose the physics from the list of registered models */ 620 { 621 PetscErrorCode (*r)(FVCtx *); 622 PetscCall(PetscFunctionListFind(physics, physname, &r)); 623 PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname); 624 /* Create the physics, will set the number of fields and their names */ 625 PetscCall((*r)(&ctx)); 626 } 627 628 /* Create a DMDA to manage the parallel grid */ 629 PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics.dof, 2, NULL, &da)); 630 PetscCall(DMSetFromOptions(da)); 631 PetscCall(DMSetUp(da)); 632 /* Inform the DMDA of the field names provided by the physics. */ 633 /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 634 for (i = 0; i < ctx.physics.dof; i++) { PetscCall(DMDASetFieldName(da, i, ctx.physics.fieldname[i])); } 635 PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 636 PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 637 638 /* Set coordinates of cell centers */ 639 PetscCall(DMDASetUniformCoordinates(da, ctx.xmin + 0.5 * (ctx.xmax - ctx.xmin) / Mx, ctx.xmax + 0.5 * (ctx.xmax - ctx.xmin) / Mx, 0, 0, 0, 0)); 640 641 /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 642 PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope)); 643 PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds)); 644 645 /* Create a vector to store the solution and to save the initial state */ 646 PetscCall(DMCreateGlobalVector(da, &X)); 647 PetscCall(VecDuplicate(X, &X0)); 648 PetscCall(VecDuplicate(X, &R)); 649 650 /*-------------------------------- create index for slow parts and fast parts ----------------------------------------*/ 651 PetscCall(PetscMalloc1(xm * dof, &index_slow)); 652 PetscCall(PetscMalloc1(xm * dof, &index_fast)); 653 for (i = xs; i < xs + xm; i++) { 654 if (ctx.xmin + i * (ctx.xmax - ctx.xmin) / (PetscReal)Mx + 0.5 * (ctx.xmax - ctx.xmin) / (PetscReal)Mx < 0) 655 for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k; 656 else 657 for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k; 658 } /* this step need to be changed based on discontinuous point of a */ 659 PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss)); 660 PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf)); 661 662 /* Create a time-stepping object */ 663 PetscCall(TSCreate(comm, &ts)); 664 PetscCall(TSSetDM(ts, da)); 665 PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction, &ctx)); 666 PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss)); 667 PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf)); 668 PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow, &ctx)); 669 PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast, &ctx)); 670 671 if (ctx.recursive) { 672 TS subts; 673 islow = 0; 674 ifast = 0; 675 for (i = xs; i < xs + xm; i++) { 676 PetscReal coord = ctx.xmin + i * (ctx.xmax - ctx.xmin) / (PetscReal)Mx + 0.5 * (ctx.xmax - ctx.xmin) / (PetscReal)Mx; 677 if (coord >= 0 && coord < ctx.xmin + (ctx.xmax - ctx.xmin) * 3 / 4.) 678 for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k; 679 if (coord >= ctx.xmin + (ctx.xmax - ctx.xmin) * 3 / 4.) 680 for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k; 681 } /* this step need to be changed based on discontinuous point of a */ 682 PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss2)); 683 PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf2)); 684 685 PetscCall(TSRHSSplitGetSubTS(ts, "fast", &subts)); 686 PetscCall(TSRHSSplitSetIS(subts, "slow", ctx.iss2)); 687 PetscCall(TSRHSSplitSetIS(subts, "fast", ctx.isf2)); 688 PetscCall(TSRHSSplitSetRHSFunction(subts, "slow", NULL, FVRHSFunctionslow2, &ctx)); 689 PetscCall(TSRHSSplitSetRHSFunction(subts, "fast", NULL, FVRHSFunctionfast2, &ctx)); 690 } 691 692 /*PetscCall(TSSetType(ts,TSSSP));*/ 693 PetscCall(TSSetType(ts, TSMPRK)); 694 PetscCall(TSSetMaxTime(ts, 10)); 695 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 696 697 /* Compute initial conditions and starting time step */ 698 PetscCall(FVSample(&ctx, da, 0, X0)); 699 PetscCall(FVRHSFunction(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */ 700 PetscCall(VecCopy(X0, X)); /* The function value was not used so we set X=X0 again */ 701 PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt)); 702 PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 703 PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 704 { 705 PetscInt steps; 706 PetscScalar mass_initial, mass_final, mass_difference; 707 708 PetscCall(TSSolve(ts, X)); 709 PetscCall(TSGetSolveTime(ts, &ptime)); 710 PetscCall(TSGetStepNumber(ts, &steps)); 711 PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps)); 712 /* calculate the total mass at initial time and final time */ 713 mass_initial = 0.0; 714 mass_final = 0.0; 715 PetscCall(VecSum(X0, &mass_initial)); 716 PetscCall(VecSum(X, &mass_final)); 717 mass_difference = (ctx.xmax - ctx.xmin) / (PetscScalar)Mx * (mass_final - mass_initial); 718 PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_difference)); 719 if (ctx.simulation) { 720 PetscViewer fd; 721 char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 722 Vec XR; 723 PetscReal nrm1, nrmsup; 724 PetscBool flg; 725 726 PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg)); 727 PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option"); 728 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd)); 729 PetscCall(VecDuplicate(X0, &XR)); 730 PetscCall(VecLoad(XR, fd)); 731 PetscCall(PetscViewerDestroy(&fd)); 732 PetscCall(VecAYPX(XR, -1, X)); 733 PetscCall(VecNorm(XR, NORM_1, &nrm1)); 734 PetscCall(VecNorm(XR, NORM_INFINITY, &nrmsup)); 735 nrm1 /= Mx; 736 PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g ||x-x_e||_sup %g\n", (double)nrm1, (double)nrmsup)); 737 PetscCall(VecDestroy(&XR)); 738 } 739 } 740 741 PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 742 if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD)); 743 if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD)); 744 if (draw & 0x4) { 745 Vec Y; 746 PetscCall(VecDuplicate(X, &Y)); 747 PetscCall(FVSample(&ctx, da, ptime, Y)); 748 PetscCall(VecAYPX(Y, -1, X)); 749 PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD)); 750 PetscCall(VecDestroy(&Y)); 751 } 752 753 if (view_final) { 754 PetscViewer viewer; 755 PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer)); 756 PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB)); 757 PetscCall(VecView(X, viewer)); 758 PetscCall(PetscViewerPopFormat(viewer)); 759 PetscCall(PetscViewerDestroy(&viewer)); 760 } 761 762 /* Clean up */ 763 PetscCall((*ctx.physics.destroy)(ctx.physics.user)); 764 for (i = 0; i < ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i])); 765 PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope)); 766 PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds)); 767 PetscCall(ISDestroy(&ctx.iss)); 768 PetscCall(ISDestroy(&ctx.isf)); 769 PetscCall(ISDestroy(&ctx.iss2)); 770 PetscCall(ISDestroy(&ctx.isf2)); 771 PetscCall(VecDestroy(&X)); 772 PetscCall(VecDestroy(&X0)); 773 PetscCall(VecDestroy(&R)); 774 PetscCall(DMDestroy(&da)); 775 PetscCall(TSDestroy(&ts)); 776 PetscCall(PetscFree(index_slow)); 777 PetscCall(PetscFree(index_fast)); 778 PetscCall(PetscFunctionListDestroy(&limiters)); 779 PetscCall(PetscFunctionListDestroy(&physics)); 780 PetscCall(PetscFinalize()); 781 return 0; 782 } 783 784 /*TEST 785 build: 786 requires: !complex 787 depends: finitevolume1d.c 788 789 test: 790 args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 791 792 test: 793 suffix: 2 794 args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 795 output_file: output/ex5_1.out 796 797 test: 798 suffix: 3 799 args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0 800 801 test: 802 suffix: 4 803 args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 804 output_file: output/ex5_3.out 805 806 test: 807 suffix: 5 808 args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 -recursive_split 809 810 test: 811 suffix: 6 812 args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 -recursive_split 813 TEST*/ 814