1 static char help[] = "Second Order TVD Finite Volume Example.\n"; 2 /*F 3 4 We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops 5 over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values, 6 \begin{equation} 7 f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R) 8 \end{equation} 9 and then update the cell values given the cell volume. 10 \begin{eqnarray} 11 f^L_i &-=& \frac{f_i}{vol^L} \\ 12 f^R_i &+=& \frac{f_i}{vol^R} 13 \end{eqnarray} 14 15 As an example, we can consider the shallow water wave equation, 16 \begin{eqnarray} 17 h_t + \nabla\cdot \left( uh \right) &=& 0 \\ 18 (uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0 19 \end{eqnarray} 20 where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity. 21 22 A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function, 23 \begin{eqnarray} 24 f^{L,R}_h &=& uh^{L,R} \cdot \hat n \\ 25 f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\ 26 c^{L,R} &=& \sqrt{g h^{L,R}} \\ 27 s &=& \max\left( \left|\frac{uh^L \cdot \hat n}{h^L}\right| + c^L, \left|\frac{uh^R \cdot \hat n}{h^R}\right| + c^R \right) \\ 28 f_i &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right) 29 \end{eqnarray} 30 where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux. 31 32 The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial 33 over a neighborhood of the given element. 34 35 The mesh is read in from an ExodusII file, usually generated by Cubit. 36 37 The example also shows how to handle AMR in a time-dependent TS solver. 38 F*/ 39 #include <petscdmplex.h> 40 #include <petscdmforest.h> 41 #include <petscds.h> 42 #include <petscts.h> 43 44 #include "ex11.h" 45 46 static PetscFunctionList PhysicsList, PhysicsRiemannList_SW, PhysicsRiemannList_Euler; 47 48 /* 'User' implements a discretization of a continuous model. */ 49 typedef struct _n_User *User; 50 typedef PetscErrorCode (*SolutionFunction)(Model, PetscReal, const PetscReal *, PetscScalar *, void *); 51 typedef PetscErrorCode (*SetUpBCFunction)(DM, PetscDS, Physics); 52 typedef PetscErrorCode (*FunctionalFunction)(Model, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); 53 typedef PetscErrorCode (*SetupFields)(Physics, PetscSection); 54 static PetscErrorCode ModelSolutionSetDefault(Model, SolutionFunction, void *); 55 static PetscErrorCode ModelFunctionalRegister(Model, const char *, PetscInt *, FunctionalFunction, void *); 56 static PetscErrorCode OutputVTK(DM, const char *, PetscViewer *); 57 58 typedef struct _n_FunctionalLink *FunctionalLink; 59 struct _n_FunctionalLink { 60 char *name; 61 FunctionalFunction func; 62 void *ctx; 63 PetscInt offset; 64 FunctionalLink next; 65 }; 66 67 struct _n_Model { 68 MPI_Comm comm; /* Does not do collective communication, but some error conditions can be collective */ 69 Physics physics; 70 FunctionalLink functionalRegistry; 71 PetscInt maxComputed; 72 PetscInt numMonitored; 73 FunctionalLink *functionalMonitored; 74 PetscInt numCall; 75 FunctionalLink *functionalCall; 76 SolutionFunction solution; 77 SetUpBCFunction setupbc; 78 void *solutionctx; 79 PetscReal maxspeed; /* estimate of global maximum speed (for CFL calculation) */ 80 PetscReal bounds[2 * DIM]; 81 PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *); 82 void *errorCtx; 83 PetscErrorCode (*setupCEED)(DM, Physics); 84 }; 85 86 struct _n_User { 87 PetscInt vtkInterval; /* For monitor */ 88 char outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */ 89 PetscInt monitorStepOffset; 90 Model model; 91 PetscBool vtkmon; 92 }; 93 94 #ifdef PETSC_HAVE_LIBCEED 95 // Free a plain data context that was allocated using PETSc; returning libCEED error codes 96 static int FreeContextPetsc(void *data) 97 { 98 if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 99 return CEED_ERROR_SUCCESS; 100 } 101 #endif 102 103 /******************* Advect ********************/ 104 typedef enum { 105 ADVECT_SOL_TILTED, 106 ADVECT_SOL_BUMP, 107 ADVECT_SOL_BUMP_CAVITY 108 } AdvectSolType; 109 static const char *const AdvectSolTypes[] = {"TILTED", "BUMP", "BUMP_CAVITY", "AdvectSolType", "ADVECT_SOL_", 0}; 110 typedef enum { 111 ADVECT_SOL_BUMP_CONE, 112 ADVECT_SOL_BUMP_COS 113 } AdvectSolBumpType; 114 static const char *const AdvectSolBumpTypes[] = {"CONE", "COS", "AdvectSolBumpType", "ADVECT_SOL_BUMP_", 0}; 115 116 typedef struct { 117 PetscReal wind[DIM]; 118 } Physics_Advect_Tilted; 119 typedef struct { 120 PetscReal center[DIM]; 121 PetscReal radius; 122 AdvectSolBumpType type; 123 } Physics_Advect_Bump; 124 125 typedef struct { 126 PetscReal inflowState; 127 AdvectSolType soltype; 128 union 129 { 130 Physics_Advect_Tilted tilted; 131 Physics_Advect_Bump bump; 132 } sol; 133 struct { 134 PetscInt Solution; 135 PetscInt Error; 136 } functional; 137 } Physics_Advect; 138 139 static const struct FieldDescription PhysicsFields_Advect[] = { 140 {"U", 1}, 141 {NULL, 0} 142 }; 143 144 static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 145 { 146 Physics phys = (Physics)ctx; 147 Physics_Advect *advect = (Physics_Advect *)phys->data; 148 149 PetscFunctionBeginUser; 150 xG[0] = advect->inflowState; 151 PetscFunctionReturn(PETSC_SUCCESS); 152 } 153 154 static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 155 { 156 PetscFunctionBeginUser; 157 xG[0] = xI[0]; 158 PetscFunctionReturn(PETSC_SUCCESS); 159 } 160 161 static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 162 { 163 Physics_Advect *advect = (Physics_Advect *)phys->data; 164 PetscReal wind[DIM], wn; 165 166 switch (advect->soltype) { 167 case ADVECT_SOL_TILTED: { 168 Physics_Advect_Tilted *tilted = &advect->sol.tilted; 169 wind[0] = tilted->wind[0]; 170 wind[1] = tilted->wind[1]; 171 } break; 172 case ADVECT_SOL_BUMP: 173 wind[0] = -qp[1]; 174 wind[1] = qp[0]; 175 break; 176 case ADVECT_SOL_BUMP_CAVITY: { 177 PetscInt i; 178 PetscReal comp2[3] = {0., 0., 0.}, rad2; 179 180 rad2 = 0.; 181 for (i = 0; i < dim; i++) { 182 comp2[i] = qp[i] * qp[i]; 183 rad2 += comp2[i]; 184 } 185 186 wind[0] = -qp[1]; 187 wind[1] = qp[0]; 188 if (rad2 > 1.) { 189 PetscInt maxI = 0; 190 PetscReal maxComp2 = comp2[0]; 191 192 for (i = 1; i < dim; i++) { 193 if (comp2[i] > maxComp2) { 194 maxI = i; 195 maxComp2 = comp2[i]; 196 } 197 } 198 wind[maxI] = 0.; 199 } 200 } break; 201 default: { 202 PetscInt i; 203 for (i = 0; i < DIM; ++i) wind[i] = 0.0; 204 } 205 /* default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */ 206 } 207 wn = Dot2Real(wind, n); 208 flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn; 209 } 210 211 static PetscErrorCode PhysicsSolution_Advect(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 212 { 213 Physics phys = (Physics)ctx; 214 Physics_Advect *advect = (Physics_Advect *)phys->data; 215 216 PetscFunctionBeginUser; 217 switch (advect->soltype) { 218 case ADVECT_SOL_TILTED: { 219 PetscReal x0[DIM]; 220 Physics_Advect_Tilted *tilted = &advect->sol.tilted; 221 Waxpy2Real(-time, tilted->wind, x, x0); 222 if (x0[1] > 0) u[0] = 1. * x[0] + 3. * x[1]; 223 else u[0] = advect->inflowState; 224 } break; 225 case ADVECT_SOL_BUMP_CAVITY: 226 case ADVECT_SOL_BUMP: { 227 Physics_Advect_Bump *bump = &advect->sol.bump; 228 PetscReal x0[DIM], v[DIM], r, cost, sint; 229 cost = PetscCosReal(time); 230 sint = PetscSinReal(time); 231 x0[0] = cost * x[0] + sint * x[1]; 232 x0[1] = -sint * x[0] + cost * x[1]; 233 Waxpy2Real(-1, bump->center, x0, v); 234 r = Norm2Real(v); 235 switch (bump->type) { 236 case ADVECT_SOL_BUMP_CONE: 237 u[0] = PetscMax(1 - r / bump->radius, 0); 238 break; 239 case ADVECT_SOL_BUMP_COS: 240 u[0] = 0.5 + 0.5 * PetscCosReal(PetscMin(r / bump->radius, 1) * PETSC_PI); 241 break; 242 } 243 } break; 244 default: 245 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 246 } 247 PetscFunctionReturn(PETSC_SUCCESS); 248 } 249 250 static PetscErrorCode PhysicsFunctional_Advect(Model mod, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 251 { 252 Physics phys = (Physics)ctx; 253 Physics_Advect *advect = (Physics_Advect *)phys->data; 254 PetscScalar yexact[1] = {0.0}; 255 256 PetscFunctionBeginUser; 257 PetscCall(PhysicsSolution_Advect(mod, time, x, yexact, phys)); 258 f[advect->functional.Solution] = PetscRealPart(y[0]); 259 f[advect->functional.Error] = PetscAbsScalar(y[0] - yexact[0]); 260 PetscFunctionReturn(PETSC_SUCCESS); 261 } 262 263 static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys) 264 { 265 const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101}; 266 DMLabel label; 267 268 PetscFunctionBeginUser; 269 /* Register "canned" boundary conditions and defaults for where to apply. */ 270 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 271 PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow", label, PETSC_STATIC_ARRAY_LENGTH(inflowids), inflowids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Advect_Inflow, NULL, phys, NULL)); 272 PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", label, PETSC_STATIC_ARRAY_LENGTH(outflowids), outflowids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Advect_Outflow, NULL, phys, NULL)); 273 PetscFunctionReturn(PETSC_SUCCESS); 274 } 275 276 static PetscErrorCode PhysicsCreate_Advect(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject) 277 { 278 Physics_Advect *advect; 279 280 PetscFunctionBeginUser; 281 phys->field_desc = PhysicsFields_Advect; 282 phys->riemann = (PetscRiemannFunc)PhysicsRiemann_Advect; 283 PetscCall(PetscNew(&advect)); 284 phys->data = advect; 285 mod->setupbc = SetUpBC_Advect; 286 287 PetscOptionsHeadBegin(PetscOptionsObject, "Advect options"); 288 { 289 PetscInt two = 2, dof = 1; 290 advect->soltype = ADVECT_SOL_TILTED; 291 PetscCall(PetscOptionsEnum("-advect_sol_type", "solution type", "", AdvectSolTypes, (PetscEnum)advect->soltype, (PetscEnum *)&advect->soltype, NULL)); 292 switch (advect->soltype) { 293 case ADVECT_SOL_TILTED: { 294 Physics_Advect_Tilted *tilted = &advect->sol.tilted; 295 two = 2; 296 tilted->wind[0] = 0.0; 297 tilted->wind[1] = 1.0; 298 PetscCall(PetscOptionsRealArray("-advect_tilted_wind", "background wind vx,vy", "", tilted->wind, &two, NULL)); 299 advect->inflowState = -2.0; 300 PetscCall(PetscOptionsRealArray("-advect_tilted_inflow", "Inflow state", "", &advect->inflowState, &dof, NULL)); 301 phys->maxspeed = Norm2Real(tilted->wind); 302 } break; 303 case ADVECT_SOL_BUMP_CAVITY: 304 case ADVECT_SOL_BUMP: { 305 Physics_Advect_Bump *bump = &advect->sol.bump; 306 two = 2; 307 bump->center[0] = 2.; 308 bump->center[1] = 0.; 309 PetscCall(PetscOptionsRealArray("-advect_bump_center", "location of center of bump x,y", "", bump->center, &two, NULL)); 310 bump->radius = 0.9; 311 PetscCall(PetscOptionsReal("-advect_bump_radius", "radius of bump", "", bump->radius, &bump->radius, NULL)); 312 bump->type = ADVECT_SOL_BUMP_CONE; 313 PetscCall(PetscOptionsEnum("-advect_bump_type", "type of bump", "", AdvectSolBumpTypes, (PetscEnum)bump->type, (PetscEnum *)&bump->type, NULL)); 314 phys->maxspeed = 3.; /* radius of mesh, kludge */ 315 } break; 316 } 317 } 318 PetscOptionsHeadEnd(); 319 /* Initial/transient solution with default boundary conditions */ 320 PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Advect, phys)); 321 /* Register "canned" functionals */ 322 PetscCall(ModelFunctionalRegister(mod, "Solution", &advect->functional.Solution, PhysicsFunctional_Advect, phys)); 323 PetscCall(ModelFunctionalRegister(mod, "Error", &advect->functional.Error, PhysicsFunctional_Advect, phys)); 324 PetscFunctionReturn(PETSC_SUCCESS); 325 } 326 327 /******************* Shallow Water ********************/ 328 static const struct FieldDescription PhysicsFields_SW[] = { 329 {"Height", 1 }, 330 {"Momentum", DIM}, 331 {NULL, 0 } 332 }; 333 334 static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 335 { 336 PetscFunctionBeginUser; 337 xG[0] = xI[0]; 338 xG[1] = -xI[1]; 339 xG[2] = -xI[2]; 340 PetscFunctionReturn(PETSC_SUCCESS); 341 } 342 343 static PetscErrorCode PhysicsSolution_SW(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 344 { 345 PetscReal dx[2], r, sigma; 346 347 PetscFunctionBeginUser; 348 PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time); 349 dx[0] = x[0] - 1.5; 350 dx[1] = x[1] - 1.0; 351 r = Norm2Real(dx); 352 sigma = 0.5; 353 u[0] = 1 + 2 * PetscExpReal(-PetscSqr(r) / (2 * PetscSqr(sigma))); 354 u[1] = 0.0; 355 u[2] = 0.0; 356 PetscFunctionReturn(PETSC_SUCCESS); 357 } 358 359 static PetscErrorCode PhysicsFunctional_SW(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx) 360 { 361 Physics phys = (Physics)ctx; 362 Physics_SW *sw = (Physics_SW *)phys->data; 363 const SWNode *x = (const SWNode *)xx; 364 PetscReal u[2]; 365 PetscReal h; 366 367 PetscFunctionBeginUser; 368 h = x->h; 369 Scale2Real(1. / x->h, x->uh, u); 370 f[sw->functional.Height] = h; 371 f[sw->functional.Speed] = Norm2Real(u) + PetscSqrtReal(sw->gravity * h); 372 f[sw->functional.Energy] = 0.5 * (Dot2Real(x->uh, u) + sw->gravity * PetscSqr(h)); 373 PetscFunctionReturn(PETSC_SUCCESS); 374 } 375 376 static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob, Physics phys) 377 { 378 const PetscInt wallids[] = {100, 101, 200, 300}; 379 DMLabel label; 380 381 PetscFunctionBeginUser; 382 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 383 PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_SW_Wall, NULL, phys, NULL)); 384 PetscFunctionReturn(PETSC_SUCCESS); 385 } 386 387 #ifdef PETSC_HAVE_LIBCEED 388 static PetscErrorCode CreateQFunctionContext_SW(Physics phys, Ceed ceed, CeedQFunctionContext *qfCtx) 389 { 390 Physics_SW *in = (Physics_SW *)phys->data; 391 Physics_SW *sw; 392 393 PetscFunctionBeginUser; 394 PetscCall(PetscCalloc1(1, &sw)); 395 396 sw->gravity = in->gravity; 397 398 PetscCallCEED(CeedQFunctionContextCreate(ceed, qfCtx)); 399 PetscCallCEED(CeedQFunctionContextSetData(*qfCtx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*sw), sw)); 400 PetscCallCEED(CeedQFunctionContextSetDataDestroy(*qfCtx, CEED_MEM_HOST, FreeContextPetsc)); 401 PetscCallCEED(CeedQFunctionContextRegisterDouble(*qfCtx, "gravity", offsetof(Physics_SW, gravity), 1, "Accelaration due to gravity")); 402 PetscFunctionReturn(0); 403 } 404 #endif 405 406 static PetscErrorCode SetupCEED_SW(DM dm, Physics physics) 407 { 408 #ifdef PETSC_HAVE_LIBCEED 409 Ceed ceed; 410 CeedQFunctionContext qfCtx; 411 #endif 412 413 PetscFunctionBegin; 414 #ifdef PETSC_HAVE_LIBCEED 415 PetscCall(DMGetCeed(dm, &ceed)); 416 PetscCall(CreateQFunctionContext_SW(physics, ceed, &qfCtx)); 417 PetscCall(DMCeedCreateFVM(dm, PETSC_TRUE, PhysicsRiemann_SW_Rusanov_CEED, PhysicsRiemann_SW_Rusanov_CEED_loc, qfCtx)); 418 PetscCallCEED(CeedQFunctionContextDestroy(&qfCtx)); 419 #endif 420 PetscFunctionReturn(PETSC_SUCCESS); 421 } 422 423 static PetscErrorCode PhysicsCreate_SW(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject) 424 { 425 Physics_SW *sw; 426 char sw_riemann[64] = "rusanov"; 427 428 PetscFunctionBeginUser; 429 phys->field_desc = PhysicsFields_SW; 430 PetscCall(PetscNew(&sw)); 431 phys->data = sw; 432 mod->setupbc = SetUpBC_SW; 433 mod->setupCEED = SetupCEED_SW; 434 435 PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov)); 436 PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL)); 437 #ifdef PETSC_HAVE_LIBCEED 438 PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov_ceed", PhysicsRiemann_SW_Rusanov_CEED)); 439 #endif 440 441 PetscOptionsHeadBegin(PetscOptionsObject, "SW options"); 442 { 443 void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics); 444 sw->gravity = 1.0; 445 PetscCall(PetscOptionsReal("-sw_gravity", "Gravitational constant", "", sw->gravity, &sw->gravity, NULL)); 446 PetscCall(PetscOptionsFList("-sw_riemann", "Riemann solver", "", PhysicsRiemannList_SW, sw_riemann, sw_riemann, sizeof sw_riemann, NULL)); 447 PetscCall(PetscFunctionListFind(PhysicsRiemannList_SW, sw_riemann, &PhysicsRiemann_SW)); 448 phys->riemann = (PetscRiemannFunc)PhysicsRiemann_SW; 449 } 450 PetscOptionsHeadEnd(); 451 phys->maxspeed = PetscSqrtReal(2.0 * sw->gravity); /* Mach 1 for depth of 2 */ 452 453 PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_SW, phys)); 454 PetscCall(ModelFunctionalRegister(mod, "Height", &sw->functional.Height, PhysicsFunctional_SW, phys)); 455 PetscCall(ModelFunctionalRegister(mod, "Speed", &sw->functional.Speed, PhysicsFunctional_SW, phys)); 456 PetscCall(ModelFunctionalRegister(mod, "Energy", &sw->functional.Energy, PhysicsFunctional_SW, phys)); 457 458 PetscFunctionReturn(PETSC_SUCCESS); 459 } 460 461 /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/ 462 /* An initial-value and self-similar solutions of the compressible Euler equations */ 463 /* Ravi Samtaney and D. I. Pullin */ 464 /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */ 465 static const struct FieldDescription PhysicsFields_Euler[] = { 466 {"Density", 1 }, 467 {"Momentum", DIM}, 468 {"Energy", 1 }, 469 {NULL, 0 } 470 }; 471 472 /* initial condition */ 473 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx); 474 static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 475 { 476 PetscInt i; 477 Physics phys = (Physics)ctx; 478 Physics_Euler *eu = (Physics_Euler *)phys->data; 479 EulerNode *uu = (EulerNode *)u; 480 PetscReal p0, gamma, c = 0.; 481 PetscFunctionBeginUser; 482 PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time); 483 484 for (i = 0; i < DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */ 485 /* set E and rho */ 486 gamma = eu->gamma; 487 488 if (eu->type == EULER_IV_SHOCK || eu->type == EULER_SS_SHOCK) { 489 /******************* Euler Density Shock ********************/ 490 /* On initial-value and self-similar solutions of the compressible Euler equations */ 491 /* Ravi Samtaney and D. I. Pullin */ 492 /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */ 493 /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity, */ 494 p0 = 1.; 495 if (x[0] < 0.0 + x[1] * eu->itana) { 496 if (x[0] < mod->bounds[0] * 0.5) { /* left of shock (1) */ 497 PetscReal amach, rho, press, gas1, p1; 498 amach = eu->amach; 499 rho = 1.; 500 press = p0; 501 p1 = press * (1.0 + 2.0 * gamma / (gamma + 1.0) * (amach * amach - 1.0)); 502 gas1 = (gamma - 1.0) / (gamma + 1.0); 503 uu->r = rho * (p1 / press + gas1) / (gas1 * p1 / press + 1.0); 504 uu->ru[0] = ((uu->r - rho) * PetscSqrtReal(gamma * press / rho) * amach); 505 uu->E = p1 / (gamma - 1.0) + .5 / uu->r * uu->ru[0] * uu->ru[0]; 506 } else { /* left of discontinuity (0) */ 507 uu->r = 1.; /* rho = 1 */ 508 uu->E = p0 / (gamma - 1.0); 509 } 510 } else { /* right of discontinuity (2) */ 511 uu->r = eu->rhoR; 512 uu->E = p0 / (gamma - 1.0); 513 } 514 } else if (eu->type == EULER_SHOCK_TUBE) { 515 /* For (x<x0) set (rho,u,p)=(8,0,10) and for (x>x0) set (rho,u,p)=(1,0,1). Choose x0 to the midpoint of the domain in the x-direction. */ 516 if (x[0] < 0.0) { 517 uu->r = 8.; 518 uu->E = 10. / (gamma - 1.); 519 } else { 520 uu->r = 1.; 521 uu->E = 1. / (gamma - 1.); 522 } 523 } else if (eu->type == EULER_LINEAR_WAVE) { 524 initLinearWave(uu, gamma, x, mod->bounds[1] - mod->bounds[0]); 525 } else SETERRQ(mod->comm, PETSC_ERR_SUP, "Unknown type %d", eu->type); 526 527 /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */ 528 PetscCall(SpeedOfSound_PG(gamma, uu, &c)); 529 c = (uu->ru[0] / uu->r) + c; 530 if (c > phys->maxspeed) phys->maxspeed = c; 531 532 PetscFunctionReturn(PETSC_SUCCESS); 533 } 534 535 /* PetscReal* => EulerNode* conversion */ 536 static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx) 537 { 538 PetscInt i; 539 const EulerNode *xI = (const EulerNode *)a_xI; 540 EulerNode *xG = (EulerNode *)a_xG; 541 Physics phys = (Physics)ctx; 542 Physics_Euler *eu = (Physics_Euler *)phys->data; 543 PetscFunctionBeginUser; 544 xG->r = xI->r; /* ghost cell density - same */ 545 xG->E = xI->E; /* ghost cell energy - same */ 546 if (n[1] != 0.) { /* top and bottom */ 547 xG->ru[0] = xI->ru[0]; /* copy tang to wall */ 548 xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */ 549 } else { /* sides */ 550 for (i = 0; i < DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */ 551 } 552 if (eu->type == EULER_LINEAR_WAVE) { /* debug */ 553 #if 0 554 PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,(double)c[0],(double)c[1]); 555 #endif 556 } 557 PetscFunctionReturn(PETSC_SUCCESS); 558 } 559 560 static PetscErrorCode PhysicsFunctional_Euler(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, void *ctx) 561 { 562 Physics phys = (Physics)ctx; 563 Physics_Euler *eu = (Physics_Euler *)phys->data; 564 const EulerNode *x = (const EulerNode *)xx; 565 PetscReal p; 566 567 PetscFunctionBeginUser; 568 f[eu->monitor.Density] = x->r; 569 f[eu->monitor.Momentum] = NormDIM(x->ru); 570 f[eu->monitor.Energy] = x->E; 571 f[eu->monitor.Speed] = NormDIM(x->ru) / x->r; 572 PetscCall(Pressure_PG(eu->gamma, x, &p)); 573 f[eu->monitor.Pressure] = p; 574 PetscFunctionReturn(PETSC_SUCCESS); 575 } 576 577 static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob, Physics phys) 578 { 579 Physics_Euler *eu = (Physics_Euler *)phys->data; 580 DMLabel label; 581 582 PetscFunctionBeginUser; 583 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 584 if (eu->type == EULER_LINEAR_WAVE) { 585 const PetscInt wallids[] = {100, 101}; 586 PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Euler_Wall, NULL, phys, NULL)); 587 } else { 588 const PetscInt wallids[] = {100, 101, 200, 300}; 589 PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (void (*)(void))PhysicsBoundary_Euler_Wall, NULL, phys, NULL)); 590 } 591 PetscFunctionReturn(PETSC_SUCCESS); 592 } 593 594 #ifdef PETSC_HAVE_LIBCEED 595 static PetscErrorCode CreateQFunctionContext_Euler(Physics phys, Ceed ceed, CeedQFunctionContext *qfCtx) 596 { 597 Physics_Euler *in = (Physics_Euler *)phys->data; 598 Physics_Euler *eu; 599 600 PetscFunctionBeginUser; 601 PetscCall(PetscCalloc1(1, &eu)); 602 603 eu->gamma = in->gamma; 604 605 PetscCallCEED(CeedQFunctionContextCreate(ceed, qfCtx)); 606 PetscCallCEED(CeedQFunctionContextSetData(*qfCtx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*eu), eu)); 607 PetscCallCEED(CeedQFunctionContextSetDataDestroy(*qfCtx, CEED_MEM_HOST, FreeContextPetsc)); 608 PetscCallCEED(CeedQFunctionContextRegisterDouble(*qfCtx, "gamma", offsetof(Physics_Euler, gamma), 1, "Heat capacity ratio")); 609 PetscFunctionReturn(0); 610 } 611 #endif 612 613 static PetscErrorCode SetupCEED_Euler(DM dm, Physics physics) 614 { 615 #ifdef PETSC_HAVE_LIBCEED 616 Ceed ceed; 617 CeedQFunctionContext qfCtx; 618 #endif 619 620 PetscFunctionBegin; 621 #ifdef PETSC_HAVE_LIBCEED 622 PetscCall(DMGetCeed(dm, &ceed)); 623 PetscCall(CreateQFunctionContext_Euler(physics, ceed, &qfCtx)); 624 PetscCall(DMCeedCreateFVM(dm, PETSC_TRUE, PhysicsRiemann_Euler_Godunov_CEED, PhysicsRiemann_Euler_Godunov_CEED_loc, qfCtx)); 625 PetscCallCEED(CeedQFunctionContextDestroy(&qfCtx)); 626 #endif 627 PetscFunctionReturn(PETSC_SUCCESS); 628 } 629 630 static PetscErrorCode PhysicsCreate_Euler(Model mod, Physics phys, PetscOptionItems *PetscOptionsObject) 631 { 632 Physics_Euler *eu; 633 634 PetscFunctionBeginUser; 635 phys->field_desc = PhysicsFields_Euler; 636 phys->riemann = (PetscRiemannFunc)PhysicsRiemann_Euler_Godunov; 637 PetscCall(PetscNew(&eu)); 638 phys->data = eu; 639 mod->setupbc = SetUpBC_Euler; 640 mod->setupCEED = SetupCEED_Euler; 641 642 PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_Euler, "godunov", PhysicsRiemann_Euler_Godunov)); 643 #ifdef PETSC_HAVE_LIBCEED 644 PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_Euler, "godunov_ceed", PhysicsRiemann_Euler_Godunov_CEED)); 645 #endif 646 647 PetscOptionsHeadBegin(PetscOptionsObject, "Euler options"); 648 { 649 void (*PhysicsRiemann_Euler)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics); 650 PetscReal alpha; 651 char type[64] = "linear_wave"; 652 char eu_riemann[64] = "godunov"; 653 PetscBool is; 654 eu->gamma = 1.4; 655 eu->amach = 2.02; 656 eu->rhoR = 3.0; 657 eu->itana = 0.57735026918963; /* angle of Euler self similar (SS) shock */ 658 PetscCall(PetscOptionsFList("-eu_riemann", "Riemann solver", "", PhysicsRiemannList_Euler, eu_riemann, eu_riemann, sizeof eu_riemann, NULL)); 659 PetscCall(PetscFunctionListFind(PhysicsRiemannList_Euler, eu_riemann, &PhysicsRiemann_Euler)); 660 phys->riemann = (PetscRiemannFunc)PhysicsRiemann_Euler; 661 PetscCall(PetscOptionsReal("-eu_gamma", "Heat capacity ratio", "", eu->gamma, &eu->gamma, NULL)); 662 PetscCall(PetscOptionsReal("-eu_amach", "Shock speed (Mach)", "", eu->amach, &eu->amach, NULL)); 663 PetscCall(PetscOptionsReal("-eu_rho2", "Density right of discontinuity", "", eu->rhoR, &eu->rhoR, NULL)); 664 alpha = 60.; 665 PetscCall(PetscOptionsReal("-eu_alpha", "Angle of discontinuity", "", alpha, &alpha, NULL)); 666 PetscCheck(alpha > 0. && alpha <= 90., PETSC_COMM_WORLD, PETSC_ERR_SUP, "Alpha bust be > 0 and <= 90 (%g)", (double)alpha); 667 eu->itana = 1. / PetscTanReal(alpha * PETSC_PI / 180.0); 668 PetscCall(PetscOptionsString("-eu_type", "Type of Euler test", "", type, type, sizeof(type), NULL)); 669 PetscCall(PetscStrcmp(type, "linear_wave", &is)); 670 if (is) { 671 /* Remember this should be periodic */ 672 eu->type = EULER_LINEAR_WAVE; 673 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "linear_wave")); 674 } else { 675 PetscCheck(DIM == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "DIM must be 2 unless linear wave test %s", type); 676 PetscCall(PetscStrcmp(type, "iv_shock", &is)); 677 if (is) { 678 eu->type = EULER_IV_SHOCK; 679 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "iv_shock")); 680 } else { 681 PetscCall(PetscStrcmp(type, "ss_shock", &is)); 682 if (is) { 683 eu->type = EULER_SS_SHOCK; 684 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "ss_shock")); 685 } else { 686 PetscCall(PetscStrcmp(type, "shock_tube", &is)); 687 if (is) eu->type = EULER_SHOCK_TUBE; 688 else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unknown Euler type %s", type); 689 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "shock_tube")); 690 } 691 } 692 } 693 } 694 PetscOptionsHeadEnd(); 695 phys->maxspeed = 0.; /* will get set in solution */ 696 PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Euler, phys)); 697 PetscCall(ModelFunctionalRegister(mod, "Speed", &eu->monitor.Speed, PhysicsFunctional_Euler, phys)); 698 PetscCall(ModelFunctionalRegister(mod, "Energy", &eu->monitor.Energy, PhysicsFunctional_Euler, phys)); 699 PetscCall(ModelFunctionalRegister(mod, "Density", &eu->monitor.Density, PhysicsFunctional_Euler, phys)); 700 PetscCall(ModelFunctionalRegister(mod, "Momentum", &eu->monitor.Momentum, PhysicsFunctional_Euler, phys)); 701 PetscCall(ModelFunctionalRegister(mod, "Pressure", &eu->monitor.Pressure, PhysicsFunctional_Euler, phys)); 702 703 PetscFunctionReturn(PETSC_SUCCESS); 704 } 705 706 static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx) 707 { 708 PetscReal err = 0.; 709 PetscInt i, j; 710 711 PetscFunctionBeginUser; 712 for (i = 0; i < numComps; i++) { 713 for (j = 0; j < dim; j++) err += PetscSqr(PetscRealPart(grad[i * dim + j])); 714 } 715 *error = volume * err; 716 PetscFunctionReturn(PETSC_SUCCESS); 717 } 718 719 PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition) 720 { 721 PetscSF sfPoint; 722 PetscSection coordSection; 723 Vec coordinates; 724 PetscSection sectionCell; 725 PetscScalar *part; 726 PetscInt cStart, cEnd, c; 727 PetscMPIInt rank; 728 729 PetscFunctionBeginUser; 730 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 731 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 732 PetscCall(DMClone(dm, dmCell)); 733 PetscCall(DMGetPointSF(dm, &sfPoint)); 734 PetscCall(DMSetPointSF(*dmCell, sfPoint)); 735 PetscCall(DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection)); 736 PetscCall(DMSetCoordinatesLocal(*dmCell, coordinates)); 737 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 738 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionCell)); 739 PetscCall(DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd)); 740 PetscCall(PetscSectionSetChart(sectionCell, cStart, cEnd)); 741 for (c = cStart; c < cEnd; ++c) PetscCall(PetscSectionSetDof(sectionCell, c, 1)); 742 PetscCall(PetscSectionSetUp(sectionCell)); 743 PetscCall(DMSetLocalSection(*dmCell, sectionCell)); 744 PetscCall(PetscSectionDestroy(§ionCell)); 745 PetscCall(DMCreateLocalVector(*dmCell, partition)); 746 PetscCall(PetscObjectSetName((PetscObject)*partition, "partition")); 747 PetscCall(VecGetArray(*partition, &part)); 748 for (c = cStart; c < cEnd; ++c) { 749 PetscScalar *p; 750 751 PetscCall(DMPlexPointLocalRef(*dmCell, c, part, &p)); 752 p[0] = rank; 753 } 754 PetscCall(VecRestoreArray(*partition, &part)); 755 PetscFunctionReturn(PETSC_SUCCESS); 756 } 757 758 PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user) 759 { 760 DM plex, dmMass, dmFace, dmCell, dmCoord; 761 PetscSection coordSection; 762 Vec coordinates, facegeom, cellgeom; 763 PetscSection sectionMass; 764 PetscScalar *m; 765 const PetscScalar *fgeom, *cgeom, *coords; 766 PetscInt vStart, vEnd, v; 767 768 PetscFunctionBeginUser; 769 PetscCall(DMConvert(dm, DMPLEX, &plex)); 770 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 771 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 772 PetscCall(DMClone(dm, &dmMass)); 773 PetscCall(DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection)); 774 PetscCall(DMSetCoordinatesLocal(dmMass, coordinates)); 775 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionMass)); 776 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 777 PetscCall(PetscSectionSetChart(sectionMass, vStart, vEnd)); 778 for (v = vStart; v < vEnd; ++v) { 779 PetscInt numFaces; 780 781 PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces)); 782 PetscCall(PetscSectionSetDof(sectionMass, v, numFaces * numFaces)); 783 } 784 PetscCall(PetscSectionSetUp(sectionMass)); 785 PetscCall(DMSetLocalSection(dmMass, sectionMass)); 786 PetscCall(PetscSectionDestroy(§ionMass)); 787 PetscCall(DMGetLocalVector(dmMass, massMatrix)); 788 PetscCall(VecGetArray(*massMatrix, &m)); 789 PetscCall(DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL)); 790 PetscCall(VecGetDM(facegeom, &dmFace)); 791 PetscCall(VecGetArrayRead(facegeom, &fgeom)); 792 PetscCall(VecGetDM(cellgeom, &dmCell)); 793 PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 794 PetscCall(DMGetCoordinateDM(dm, &dmCoord)); 795 PetscCall(VecGetArrayRead(coordinates, &coords)); 796 for (v = vStart; v < vEnd; ++v) { 797 const PetscInt *faces; 798 PetscFVFaceGeom *fgA, *fgB, *cg; 799 PetscScalar *vertex; 800 PetscInt numFaces, sides[2], f, g; 801 802 PetscCall(DMPlexPointLocalRead(dmCoord, v, coords, &vertex)); 803 PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces)); 804 PetscCall(DMPlexGetSupport(dmMass, v, &faces)); 805 for (f = 0; f < numFaces; ++f) { 806 sides[0] = faces[f]; 807 PetscCall(DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA)); 808 for (g = 0; g < numFaces; ++g) { 809 const PetscInt *cells = NULL; 810 PetscReal area = 0.0; 811 PetscInt numCells; 812 813 sides[1] = faces[g]; 814 PetscCall(DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB)); 815 PetscCall(DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells)); 816 PetscCheck(numCells == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces"); 817 PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg)); 818 area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgA->centroid[0] - cg->centroid[0])); 819 area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgB->centroid[0] - cg->centroid[0])); 820 m[f * numFaces + g] = Dot2Real(fgA->normal, fgB->normal) * area * 0.5; 821 PetscCall(DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells)); 822 } 823 } 824 } 825 PetscCall(VecRestoreArrayRead(facegeom, &fgeom)); 826 PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 827 PetscCall(VecRestoreArrayRead(coordinates, &coords)); 828 PetscCall(VecRestoreArray(*massMatrix, &m)); 829 PetscCall(DMDestroy(&dmMass)); 830 PetscCall(DMDestroy(&plex)); 831 PetscFunctionReturn(PETSC_SUCCESS); 832 } 833 834 /* Behavior will be different for multi-physics or when using non-default boundary conditions */ 835 static PetscErrorCode ModelSolutionSetDefault(Model mod, SolutionFunction func, void *ctx) 836 { 837 PetscFunctionBeginUser; 838 mod->solution = func; 839 mod->solutionctx = ctx; 840 PetscFunctionReturn(PETSC_SUCCESS); 841 } 842 843 static PetscErrorCode ModelFunctionalRegister(Model mod, const char *name, PetscInt *offset, FunctionalFunction func, void *ctx) 844 { 845 FunctionalLink link, *ptr; 846 PetscInt lastoffset = -1; 847 848 PetscFunctionBeginUser; 849 for (ptr = &mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 850 PetscCall(PetscNew(&link)); 851 PetscCall(PetscStrallocpy(name, &link->name)); 852 link->offset = lastoffset + 1; 853 link->func = func; 854 link->ctx = ctx; 855 link->next = NULL; 856 *ptr = link; 857 *offset = link->offset; 858 PetscFunctionReturn(PETSC_SUCCESS); 859 } 860 861 static PetscErrorCode ModelFunctionalSetFromOptions(Model mod, PetscOptionItems *PetscOptionsObject) 862 { 863 PetscInt i, j; 864 FunctionalLink link; 865 char *names[256]; 866 867 PetscFunctionBeginUser; 868 mod->numMonitored = PETSC_STATIC_ARRAY_LENGTH(names); 869 PetscCall(PetscOptionsStringArray("-monitor", "list of functionals to monitor", "", names, &mod->numMonitored, NULL)); 870 /* Create list of functionals that will be computed somehow */ 871 PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalMonitored)); 872 /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */ 873 PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalCall)); 874 mod->numCall = 0; 875 for (i = 0; i < mod->numMonitored; i++) { 876 for (link = mod->functionalRegistry; link; link = link->next) { 877 PetscBool match; 878 PetscCall(PetscStrcasecmp(names[i], link->name, &match)); 879 if (match) break; 880 } 881 PetscCheck(link, mod->comm, PETSC_ERR_USER, "No known functional '%s'", names[i]); 882 mod->functionalMonitored[i] = link; 883 for (j = 0; j < i; j++) { 884 if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name; 885 } 886 mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */ 887 next_name: 888 PetscCall(PetscFree(names[i])); 889 } 890 891 /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 892 mod->maxComputed = -1; 893 for (link = mod->functionalRegistry; link; link = link->next) { 894 for (i = 0; i < mod->numCall; i++) { 895 FunctionalLink call = mod->functionalCall[i]; 896 if (link->func == call->func && link->ctx == call->ctx) mod->maxComputed = PetscMax(mod->maxComputed, link->offset); 897 } 898 } 899 PetscFunctionReturn(PETSC_SUCCESS); 900 } 901 902 static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link) 903 { 904 FunctionalLink l, next; 905 906 PetscFunctionBeginUser; 907 if (!link) PetscFunctionReturn(PETSC_SUCCESS); 908 l = *link; 909 *link = NULL; 910 for (; l; l = next) { 911 next = l->next; 912 PetscCall(PetscFree(l->name)); 913 PetscCall(PetscFree(l)); 914 } 915 PetscFunctionReturn(PETSC_SUCCESS); 916 } 917 918 /* put the solution callback into a functional callback */ 919 static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx) 920 { 921 Model mod; 922 PetscFunctionBeginUser; 923 mod = (Model)modctx; 924 PetscCall((*mod->solution)(mod, time, x, u, mod->solutionctx)); 925 PetscFunctionReturn(PETSC_SUCCESS); 926 } 927 928 PetscErrorCode SetInitialCondition(DM dm, Vec X, User user) 929 { 930 PetscErrorCode (*func[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 931 void *ctx[1]; 932 Model mod = user->model; 933 934 PetscFunctionBeginUser; 935 func[0] = SolutionFunctional; 936 ctx[0] = (void *)mod; 937 PetscCall(DMProjectFunction(dm, 0.0, func, ctx, INSERT_ALL_VALUES, X)); 938 PetscFunctionReturn(PETSC_SUCCESS); 939 } 940 941 static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer) 942 { 943 PetscFunctionBeginUser; 944 PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer)); 945 PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERVTK)); 946 PetscCall(PetscViewerFileSetName(*viewer, filename)); 947 PetscFunctionReturn(PETSC_SUCCESS); 948 } 949 950 static PetscErrorCode MonitorVTK(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 951 { 952 User user = (User)ctx; 953 DM dm, plex; 954 PetscViewer viewer; 955 char filename[PETSC_MAX_PATH_LEN], *ftable = NULL; 956 PetscReal xnorm; 957 958 PetscFunctionBeginUser; 959 PetscCall(PetscObjectSetName((PetscObject)X, "u")); 960 PetscCall(VecGetDM(X, &dm)); 961 PetscCall(VecNorm(X, NORM_INFINITY, &xnorm)); 962 963 if (stepnum >= 0) stepnum += user->monitorStepOffset; 964 if (stepnum >= 0) { /* No summary for final time */ 965 Model mod = user->model; 966 Vec cellgeom; 967 PetscInt c, cStart, cEnd, fcount, i; 968 size_t ftableused, ftablealloc; 969 const PetscScalar *cgeom, *x; 970 DM dmCell; 971 DMLabel vtkLabel; 972 PetscReal *fmin, *fmax, *fintegral, *ftmp; 973 974 PetscCall(DMConvert(dm, DMPLEX, &plex)); 975 PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL)); 976 fcount = mod->maxComputed + 1; 977 PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fintegral, fcount, &ftmp)); 978 for (i = 0; i < fcount; i++) { 979 fmin[i] = PETSC_MAX_REAL; 980 fmax[i] = PETSC_MIN_REAL; 981 fintegral[i] = 0; 982 } 983 PetscCall(VecGetDM(cellgeom, &dmCell)); 984 PetscCall(DMPlexGetSimplexOrBoxCells(dmCell, 0, &cStart, &cEnd)); 985 PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 986 PetscCall(VecGetArrayRead(X, &x)); 987 PetscCall(DMGetLabel(dm, "vtk", &vtkLabel)); 988 for (c = cStart; c < cEnd; ++c) { 989 PetscFVCellGeom *cg; 990 const PetscScalar *cx = NULL; 991 PetscInt vtkVal = 0; 992 993 /* not that these two routines as currently implemented work for any dm with a 994 * localSection/globalSection */ 995 PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 996 PetscCall(DMPlexPointGlobalRead(dm, c, x, &cx)); 997 if (vtkLabel) PetscCall(DMLabelGetValue(vtkLabel, c, &vtkVal)); 998 if (!vtkVal || !cx) continue; /* ghost, or not a global cell */ 999 for (i = 0; i < mod->numCall; i++) { 1000 FunctionalLink flink = mod->functionalCall[i]; 1001 PetscCall((*flink->func)(mod, time, cg->centroid, cx, ftmp, flink->ctx)); 1002 } 1003 for (i = 0; i < fcount; i++) { 1004 fmin[i] = PetscMin(fmin[i], ftmp[i]); 1005 fmax[i] = PetscMax(fmax[i], ftmp[i]); 1006 fintegral[i] += cg->volume * ftmp[i]; 1007 } 1008 } 1009 PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 1010 PetscCall(VecRestoreArrayRead(X, &x)); 1011 PetscCall(DMDestroy(&plex)); 1012 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 1013 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 1014 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fintegral, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 1015 1016 ftablealloc = fcount * 100; 1017 ftableused = 0; 1018 PetscCall(PetscMalloc1(ftablealloc, &ftable)); 1019 for (i = 0; i < mod->numMonitored; i++) { 1020 size_t countused; 1021 char buffer[256], *p; 1022 FunctionalLink flink = mod->functionalMonitored[i]; 1023 PetscInt id = flink->offset; 1024 if (i % 3) { 1025 PetscCall(PetscArraycpy(buffer, " ", 2)); 1026 p = buffer + 2; 1027 } else if (i) { 1028 char newline[] = "\n"; 1029 PetscCall(PetscMemcpy(buffer, newline, sizeof(newline) - 1)); 1030 p = buffer + sizeof(newline) - 1; 1031 } else { 1032 p = buffer; 1033 } 1034 PetscCall(PetscSNPrintfCount(p, sizeof buffer - (p - buffer), "%12s [%10.7g,%10.7g] int %10.7g", &countused, flink->name, (double)fmin[id], (double)fmax[id], (double)fintegral[id])); 1035 countused--; 1036 countused += p - buffer; 1037 if (countused > ftablealloc - ftableused - 1) { /* reallocate */ 1038 char *ftablenew; 1039 ftablealloc = 2 * ftablealloc + countused; 1040 PetscCall(PetscMalloc(ftablealloc, &ftablenew)); 1041 PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); 1042 PetscCall(PetscFree(ftable)); 1043 ftable = ftablenew; 1044 } 1045 PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused)); 1046 ftableused += countused; 1047 ftable[ftableused] = 0; 1048 } 1049 PetscCall(PetscFree4(fmin, fmax, fintegral, ftmp)); 1050 1051 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT " time %8.4g |x| %8.4g %s\n", stepnum, (double)time, (double)xnorm, ftable ? ftable : "")); 1052 PetscCall(PetscFree(ftable)); 1053 } 1054 if (user->vtkInterval < 1) PetscFunctionReturn(PETSC_SUCCESS); 1055 if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) { 1056 if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */ 1057 PetscCall(TSGetStepNumber(ts, &stepnum)); 1058 } 1059 PetscCall(PetscSNPrintf(filename, sizeof filename, "%s-%03" PetscInt_FMT ".vtu", user->outputBasename, stepnum)); 1060 PetscCall(OutputVTK(dm, filename, &viewer)); 1061 PetscCall(VecView(X, viewer)); 1062 PetscCall(PetscViewerDestroy(&viewer)); 1063 } 1064 PetscFunctionReturn(PETSC_SUCCESS); 1065 } 1066 1067 static PetscErrorCode initializeTS(DM dm, User user, TS *ts) 1068 { 1069 PetscBool useCeed; 1070 1071 PetscFunctionBeginUser; 1072 PetscCall(TSCreate(PetscObjectComm((PetscObject)dm), ts)); 1073 PetscCall(TSSetType(*ts, TSSSP)); 1074 PetscCall(TSSetDM(*ts, dm)); 1075 if (user->vtkmon) PetscCall(TSMonitorSet(*ts, MonitorVTK, user, NULL)); 1076 PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 1077 PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user)); 1078 if (useCeed) PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVMCEED, user)); 1079 else PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user)); 1080 PetscCall(TSSetMaxTime(*ts, 2.0)); 1081 PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 1082 PetscFunctionReturn(PETSC_SUCCESS); 1083 } 1084 1085 typedef struct { 1086 PetscFV fvm; 1087 VecTagger refineTag; 1088 VecTagger coarsenTag; 1089 DM adaptedDM; 1090 User user; 1091 PetscReal cfl; 1092 PetscLimiter limiter; 1093 PetscLimiter noneLimiter; 1094 } TransferCtx; 1095 1096 static PetscErrorCode adaptToleranceFVMSetUp(TS ts, PetscInt nstep, PetscReal time, Vec sol, PetscBool *resize, void *ctx) 1097 { 1098 TransferCtx *tctx = (TransferCtx *)ctx; 1099 PetscFV fvm = tctx->fvm; 1100 VecTagger refineTag = tctx->refineTag; 1101 VecTagger coarsenTag = tctx->coarsenTag; 1102 User user = tctx->user; 1103 DM dm, gradDM, plex, cellDM, adaptedDM = NULL; 1104 Vec cellGeom, faceGeom; 1105 PetscBool isForest, computeGradient; 1106 Vec grad, locGrad, locX, errVec; 1107 PetscInt cStart, cEnd, c, dim, nRefine, nCoarsen; 1108 PetscReal minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2]; 1109 PetscScalar *errArray; 1110 const PetscScalar *pointVals; 1111 const PetscScalar *pointGrads; 1112 const PetscScalar *pointGeom; 1113 DMLabel adaptLabel = NULL; 1114 IS refineIS, coarsenIS; 1115 1116 PetscFunctionBeginUser; 1117 *resize = PETSC_FALSE; 1118 PetscCall(VecGetDM(sol, &dm)); 1119 PetscCall(DMGetDimension(dm, &dim)); 1120 PetscCall(PetscFVSetLimiter(fvm, tctx->noneLimiter)); 1121 PetscCall(PetscFVGetComputeGradients(fvm, &computeGradient)); 1122 PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE)); 1123 PetscCall(DMIsForest(dm, &isForest)); 1124 PetscCall(DMConvert(dm, DMPLEX, &plex)); 1125 PetscCall(DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM)); 1126 PetscCall(DMCreateLocalVector(plex, &locX)); 1127 PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL)); 1128 PetscCall(DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX)); 1129 PetscCall(DMGlobalToLocalEnd(plex, sol, INSERT_VALUES, locX)); 1130 PetscCall(DMCreateGlobalVector(gradDM, &grad)); 1131 PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad)); 1132 PetscCall(DMCreateLocalVector(gradDM, &locGrad)); 1133 PetscCall(DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad)); 1134 PetscCall(DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad)); 1135 PetscCall(VecDestroy(&grad)); 1136 PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd)); 1137 PetscCall(VecGetArrayRead(locGrad, &pointGrads)); 1138 PetscCall(VecGetArrayRead(cellGeom, &pointGeom)); 1139 PetscCall(VecGetArrayRead(locX, &pointVals)); 1140 PetscCall(VecGetDM(cellGeom, &cellDM)); 1141 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 1142 PetscCall(VecCreateFromOptions(PetscObjectComm((PetscObject)plex), NULL, 1, cEnd - cStart, PETSC_DETERMINE, &errVec)); 1143 PetscCall(VecSetUp(errVec)); 1144 PetscCall(VecGetArray(errVec, &errArray)); 1145 for (c = cStart; c < cEnd; c++) { 1146 PetscReal errInd = 0.; 1147 PetscScalar *pointGrad; 1148 PetscScalar *pointVal; 1149 PetscFVCellGeom *cg; 1150 1151 PetscCall(DMPlexPointLocalRead(gradDM, c, pointGrads, &pointGrad)); 1152 PetscCall(DMPlexPointLocalRead(cellDM, c, pointGeom, &cg)); 1153 PetscCall(DMPlexPointLocalRead(plex, c, pointVals, &pointVal)); 1154 1155 PetscCall((user->model->errorIndicator)(dim, cg->volume, user->model->physics->dof, pointVal, pointGrad, &errInd, user->model->errorCtx)); 1156 errArray[c - cStart] = errInd; 1157 minMaxInd[0] = PetscMin(minMaxInd[0], errInd); 1158 minMaxInd[1] = PetscMax(minMaxInd[1], errInd); 1159 } 1160 PetscCall(VecRestoreArray(errVec, &errArray)); 1161 PetscCall(VecRestoreArrayRead(locX, &pointVals)); 1162 PetscCall(VecRestoreArrayRead(cellGeom, &pointGeom)); 1163 PetscCall(VecRestoreArrayRead(locGrad, &pointGrads)); 1164 PetscCall(VecDestroy(&locGrad)); 1165 PetscCall(VecDestroy(&locX)); 1166 PetscCall(DMDestroy(&plex)); 1167 1168 PetscCall(VecTaggerComputeIS(refineTag, errVec, &refineIS, NULL)); 1169 PetscCall(VecTaggerComputeIS(coarsenTag, errVec, &coarsenIS, NULL)); 1170 PetscCall(ISGetSize(refineIS, &nRefine)); 1171 PetscCall(ISGetSize(coarsenIS, &nCoarsen)); 1172 if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS)); 1173 if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS)); 1174 PetscCall(ISDestroy(&coarsenIS)); 1175 PetscCall(ISDestroy(&refineIS)); 1176 PetscCall(VecDestroy(&errVec)); 1177 1178 PetscCall(PetscFVSetComputeGradients(fvm, computeGradient)); 1179 PetscCall(PetscFVSetLimiter(fvm, tctx->limiter)); 1180 minMaxInd[1] = -minMaxInd[1]; 1181 PetscCall(MPIU_Allreduce(minMaxInd, minMaxIndGlobal, 2, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm))); 1182 PetscCall(PetscInfo(ts, "error indicator range (%E, %E)\n", (double)minMaxIndGlobal[0], (double)(-minMaxIndGlobal[1]))); 1183 if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */ 1184 PetscCall(DMAdaptLabel(dm, adaptLabel, &adaptedDM)); 1185 } 1186 PetscCall(DMLabelDestroy(&adaptLabel)); 1187 if (adaptedDM) { 1188 tctx->adaptedDM = adaptedDM; 1189 PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nstep, nRefine, nCoarsen)); 1190 *resize = PETSC_TRUE; 1191 } else { 1192 PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " no adaptation\n", nstep)); 1193 } 1194 PetscFunctionReturn(PETSC_SUCCESS); 1195 } 1196 1197 static PetscErrorCode Transfer(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], void *ctx) 1198 { 1199 TransferCtx *tctx = (TransferCtx *)ctx; 1200 DM dm; 1201 PetscReal time; 1202 1203 PetscFunctionBeginUser; 1204 PetscCall(TSGetDM(ts, &dm)); 1205 PetscCall(TSGetTime(ts, &time)); 1206 PetscCheck(tctx->adaptedDM, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Missing adaptedDM"); 1207 for (PetscInt i = 0; i < nv; i++) { 1208 const char *name; 1209 1210 PetscCall(DMCreateGlobalVector(tctx->adaptedDM, &vecsout[i])); 1211 PetscCall(DMForestTransferVec(dm, vecsin[i], tctx->adaptedDM, vecsout[i], PETSC_TRUE, time)); 1212 PetscCall(PetscObjectGetName((PetscObject)vecsin[i], &name)); 1213 PetscCall(PetscObjectSetName((PetscObject)vecsout[i], name)); 1214 } 1215 PetscCall(DMForestSetAdaptivityForest(tctx->adaptedDM, NULL)); /* clear internal references to the previous dm */ 1216 1217 Model mod = tctx->user->model; 1218 Physics phys = mod->physics; 1219 PetscReal minRadius; 1220 1221 PetscCall(DMPlexGetGeometryFVM(tctx->adaptedDM, NULL, NULL, &minRadius)); 1222 PetscCall(MPIU_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 1223 PetscCheck(mod->maxspeed > 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Physics did not set maxspeed"); 1224 1225 PetscReal dt = tctx->cfl * minRadius / mod->maxspeed; 1226 PetscCall(TSSetTimeStep(ts, dt)); 1227 1228 PetscCall(TSSetDM(ts, tctx->adaptedDM)); 1229 PetscCall(DMDestroy(&tctx->adaptedDM)); 1230 1231 PetscFunctionReturn(PETSC_SUCCESS); 1232 } 1233 1234 int main(int argc, char **argv) 1235 { 1236 MPI_Comm comm; 1237 PetscDS prob; 1238 PetscFV fvm; 1239 PetscLimiter limiter = NULL, noneLimiter = NULL; 1240 User user; 1241 Model mod; 1242 Physics phys; 1243 DM dm, plex; 1244 PetscReal ftime, cfl, dt, minRadius; 1245 PetscInt dim, nsteps; 1246 TS ts; 1247 TSConvergedReason reason; 1248 Vec X; 1249 PetscViewer viewer; 1250 PetscBool vtkCellGeom, useAMR; 1251 PetscInt adaptInterval; 1252 char physname[256] = "advect"; 1253 VecTagger refineTag = NULL, coarsenTag = NULL; 1254 TransferCtx tctx; 1255 1256 PetscFunctionBeginUser; 1257 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1258 comm = PETSC_COMM_WORLD; 1259 1260 PetscCall(PetscNew(&user)); 1261 PetscCall(PetscNew(&user->model)); 1262 PetscCall(PetscNew(&user->model->physics)); 1263 mod = user->model; 1264 phys = mod->physics; 1265 mod->comm = comm; 1266 useAMR = PETSC_FALSE; 1267 adaptInterval = 1; 1268 1269 /* Register physical models to be available on the command line */ 1270 PetscCall(PetscFunctionListAdd(&PhysicsList, "advect", PhysicsCreate_Advect)); 1271 PetscCall(PetscFunctionListAdd(&PhysicsList, "sw", PhysicsCreate_SW)); 1272 PetscCall(PetscFunctionListAdd(&PhysicsList, "euler", PhysicsCreate_Euler)); 1273 1274 PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Mesh Options", ""); 1275 { 1276 cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */ 1277 PetscCall(PetscOptionsReal("-ufv_cfl", "CFL number per step", "", cfl, &cfl, NULL)); 1278 user->vtkInterval = 1; 1279 PetscCall(PetscOptionsInt("-ufv_vtk_interval", "VTK output interval (0 to disable)", "", user->vtkInterval, &user->vtkInterval, NULL)); 1280 user->vtkmon = PETSC_TRUE; 1281 PetscCall(PetscOptionsBool("-ufv_vtk_monitor", "Use VTKMonitor routine", "", user->vtkmon, &user->vtkmon, NULL)); 1282 vtkCellGeom = PETSC_FALSE; 1283 PetscCall(PetscStrncpy(user->outputBasename, "ex11", sizeof(user->outputBasename))); 1284 PetscCall(PetscOptionsString("-ufv_vtk_basename", "VTK output basename", "", user->outputBasename, user->outputBasename, sizeof(user->outputBasename), NULL)); 1285 PetscCall(PetscOptionsBool("-ufv_vtk_cellgeom", "Write cell geometry (for debugging)", "", vtkCellGeom, &vtkCellGeom, NULL)); 1286 PetscCall(PetscOptionsBool("-ufv_use_amr", "use local adaptive mesh refinement", "", useAMR, &useAMR, NULL)); 1287 PetscCall(PetscOptionsInt("-ufv_adapt_interval", "time steps between AMR", "", adaptInterval, &adaptInterval, NULL)); 1288 } 1289 PetscOptionsEnd(); 1290 1291 if (useAMR) { 1292 VecTaggerBox refineBox, coarsenBox; 1293 1294 refineBox.min = refineBox.max = PETSC_MAX_REAL; 1295 coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL; 1296 1297 PetscCall(VecTaggerCreate(comm, &refineTag)); 1298 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)refineTag, "refine_")); 1299 PetscCall(VecTaggerSetType(refineTag, VECTAGGERABSOLUTE)); 1300 PetscCall(VecTaggerAbsoluteSetBox(refineTag, &refineBox)); 1301 PetscCall(VecTaggerSetFromOptions(refineTag)); 1302 PetscCall(VecTaggerSetUp(refineTag)); 1303 PetscCall(PetscObjectViewFromOptions((PetscObject)refineTag, NULL, "-tag_view")); 1304 1305 PetscCall(VecTaggerCreate(comm, &coarsenTag)); 1306 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)coarsenTag, "coarsen_")); 1307 PetscCall(VecTaggerSetType(coarsenTag, VECTAGGERABSOLUTE)); 1308 PetscCall(VecTaggerAbsoluteSetBox(coarsenTag, &coarsenBox)); 1309 PetscCall(VecTaggerSetFromOptions(coarsenTag)); 1310 PetscCall(VecTaggerSetUp(coarsenTag)); 1311 PetscCall(PetscObjectViewFromOptions((PetscObject)coarsenTag, NULL, "-tag_view")); 1312 } 1313 1314 PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Physics Options", ""); 1315 { 1316 PetscErrorCode (*physcreate)(Model, Physics, PetscOptionItems *); 1317 PetscCall(PetscOptionsFList("-physics", "Physics module to solve", "", PhysicsList, physname, physname, sizeof physname, NULL)); 1318 PetscCall(PetscFunctionListFind(PhysicsList, physname, &physcreate)); 1319 PetscCall(PetscMemzero(phys, sizeof(struct _n_Physics))); 1320 PetscCall((*physcreate)(mod, phys, PetscOptionsObject)); 1321 /* Count number of fields and dofs */ 1322 for (phys->nfields = 0, phys->dof = 0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof; 1323 PetscCheck(phys->dof > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set dof", physname); 1324 PetscCall(ModelFunctionalSetFromOptions(mod, PetscOptionsObject)); 1325 } 1326 PetscOptionsEnd(); 1327 1328 /* Create mesh */ 1329 { 1330 PetscInt i; 1331 1332 PetscCall(DMCreate(comm, &dm)); 1333 PetscCall(DMSetType(dm, DMPLEX)); 1334 PetscCall(DMSetFromOptions(dm)); 1335 for (i = 0; i < DIM; i++) { 1336 mod->bounds[2 * i] = 0.; 1337 mod->bounds[2 * i + 1] = 1.; 1338 }; 1339 dim = DIM; 1340 { /* a null name means just do a hex box */ 1341 PetscInt cells[3] = {1, 1, 1}, n = 3; 1342 PetscBool flg2, skew = PETSC_FALSE; 1343 PetscInt nret2 = 2 * DIM; 1344 PetscOptionsBegin(comm, NULL, "Rectangular mesh options", ""); 1345 PetscCall(PetscOptionsRealArray("-grid_bounds", "bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max", "", mod->bounds, &nret2, &flg2)); 1346 PetscCall(PetscOptionsBool("-grid_skew_60", "Skew grid for 60 degree shock mesh", "", skew, &skew, NULL)); 1347 PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", cells, &n, NULL)); 1348 PetscOptionsEnd(); 1349 /* TODO Rewrite this with Mark, and remove grid_bounds at that time */ 1350 if (flg2) { 1351 PetscInt dimEmbed, i; 1352 PetscInt nCoords; 1353 PetscScalar *coords; 1354 Vec coordinates; 1355 1356 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1357 PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 1358 PetscCall(VecGetLocalSize(coordinates, &nCoords)); 1359 PetscCheck(!(nCoords % dimEmbed), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Coordinate vector the wrong size"); 1360 PetscCall(VecGetArray(coordinates, &coords)); 1361 for (i = 0; i < nCoords; i += dimEmbed) { 1362 PetscInt j; 1363 1364 PetscScalar *coord = &coords[i]; 1365 for (j = 0; j < dimEmbed; j++) { 1366 coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]); 1367 if (dim == 2 && cells[1] == 1 && j == 0 && skew) { 1368 if (cells[0] == 2 && i == 8) { 1369 coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */ 1370 } else if (cells[0] == 3) { 1371 if (i == 2 || i == 10) coord[j] = mod->bounds[1] / 4.; 1372 else if (i == 4) coord[j] = mod->bounds[1] / 2.; 1373 else if (i == 12) coord[j] = 1.57735026918963 * mod->bounds[1] / 2.; 1374 } 1375 } 1376 } 1377 } 1378 PetscCall(VecRestoreArray(coordinates, &coords)); 1379 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 1380 } 1381 } 1382 } 1383 PetscCall(DMViewFromOptions(dm, NULL, "-orig_dm_view")); 1384 PetscCall(DMGetDimension(dm, &dim)); 1385 1386 /* set up BCs, functions, tags */ 1387 PetscCall(DMCreateLabel(dm, "Face Sets")); 1388 mod->errorIndicator = ErrorIndicator_Simple; 1389 1390 { 1391 DM gdm; 1392 1393 PetscCall(DMPlexConstructGhostCells(dm, NULL, NULL, &gdm)); 1394 PetscCall(DMDestroy(&dm)); 1395 dm = gdm; 1396 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 1397 } 1398 1399 PetscCall(PetscFVCreate(comm, &fvm)); 1400 PetscCall(PetscFVSetFromOptions(fvm)); 1401 PetscCall(PetscFVSetNumComponents(fvm, phys->dof)); 1402 PetscCall(PetscFVSetSpatialDimension(fvm, dim)); 1403 PetscCall(PetscObjectSetName((PetscObject)fvm, "")); 1404 { 1405 PetscInt f, dof; 1406 for (f = 0, dof = 0; f < phys->nfields; f++) { 1407 PetscInt newDof = phys->field_desc[f].dof; 1408 1409 if (newDof == 1) { 1410 PetscCall(PetscFVSetComponentName(fvm, dof, phys->field_desc[f].name)); 1411 } else { 1412 PetscInt j; 1413 1414 for (j = 0; j < newDof; j++) { 1415 char compName[256] = "Unknown"; 1416 1417 PetscCall(PetscSNPrintf(compName, sizeof(compName), "%s_%" PetscInt_FMT, phys->field_desc[f].name, j)); 1418 PetscCall(PetscFVSetComponentName(fvm, dof + j, compName)); 1419 } 1420 } 1421 dof += newDof; 1422 } 1423 } 1424 /* FV is now structured with one field having all physics as components */ 1425 PetscCall(DMAddField(dm, NULL, (PetscObject)fvm)); 1426 PetscCall(DMCreateDS(dm)); 1427 PetscCall(DMGetDS(dm, &prob)); 1428 PetscCall(PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann)); 1429 PetscCall(PetscDSSetContext(prob, 0, user->model->physics)); 1430 PetscCall((*mod->setupbc)(dm, prob, phys)); 1431 PetscCall(PetscDSSetFromOptions(prob)); 1432 { 1433 char convType[256]; 1434 PetscBool flg; 1435 1436 PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX"); 1437 PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex12", DMList, DMPLEX, convType, 256, &flg)); 1438 PetscOptionsEnd(); 1439 if (flg) { 1440 DM dmConv; 1441 1442 PetscCall(DMConvert(dm, convType, &dmConv)); 1443 if (dmConv) { 1444 PetscCall(DMViewFromOptions(dmConv, NULL, "-dm_conv_view")); 1445 PetscCall(DMDestroy(&dm)); 1446 dm = dmConv; 1447 PetscCall(DMSetFromOptions(dm)); 1448 } 1449 } 1450 } 1451 #ifdef PETSC_HAVE_LIBCEED 1452 { 1453 PetscBool useCeed; 1454 PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 1455 if (useCeed) PetscCall((*user->model->setupCEED)(dm, user->model->physics)); 1456 } 1457 #endif 1458 1459 PetscCall(initializeTS(dm, user, &ts)); 1460 1461 PetscCall(DMCreateGlobalVector(dm, &X)); 1462 PetscCall(PetscObjectSetName((PetscObject)X, "solution")); 1463 PetscCall(SetInitialCondition(dm, X, user)); 1464 if (useAMR) { 1465 PetscInt adaptIter; 1466 1467 /* use no limiting when reconstructing gradients for adaptivity */ 1468 PetscCall(PetscFVGetLimiter(fvm, &limiter)); 1469 PetscCall(PetscObjectReference((PetscObject)limiter)); 1470 PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter)); 1471 PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE)); 1472 1473 /* Refinement context */ 1474 tctx.fvm = fvm; 1475 tctx.refineTag = refineTag; 1476 tctx.coarsenTag = coarsenTag; 1477 tctx.adaptedDM = NULL; 1478 tctx.user = user; 1479 tctx.noneLimiter = noneLimiter; 1480 tctx.limiter = limiter; 1481 tctx.cfl = cfl; 1482 1483 /* Do some initial refinement steps */ 1484 for (adaptIter = 0;; ++adaptIter) { 1485 PetscLogDouble bytes; 1486 PetscBool resize; 1487 1488 PetscCall(PetscMemoryGetCurrentUsage(&bytes)); 1489 PetscCall(PetscInfo(ts, "refinement loop %" PetscInt_FMT ": memory used %g\n", adaptIter, (double)bytes)); 1490 PetscCall(DMViewFromOptions(dm, NULL, "-initial_dm_view")); 1491 PetscCall(VecViewFromOptions(X, NULL, "-initial_vec_view")); 1492 1493 PetscCall(adaptToleranceFVMSetUp(ts, -1, 0.0, X, &resize, &tctx)); 1494 if (!resize) break; 1495 PetscCall(DMDestroy(&dm)); 1496 PetscCall(VecDestroy(&X)); 1497 dm = tctx.adaptedDM; 1498 tctx.adaptedDM = NULL; 1499 PetscCall(TSSetDM(ts, dm)); 1500 PetscCall(DMCreateGlobalVector(dm, &X)); 1501 PetscCall(PetscObjectSetName((PetscObject)X, "solution")); 1502 PetscCall(SetInitialCondition(dm, X, user)); 1503 } 1504 } 1505 1506 PetscCall(DMConvert(dm, DMPLEX, &plex)); 1507 if (vtkCellGeom) { 1508 DM dmCell; 1509 Vec cellgeom, partition; 1510 1511 PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL)); 1512 PetscCall(OutputVTK(dm, "ex11-cellgeom.vtk", &viewer)); 1513 PetscCall(VecView(cellgeom, viewer)); 1514 PetscCall(PetscViewerDestroy(&viewer)); 1515 PetscCall(CreatePartitionVec(dm, &dmCell, &partition)); 1516 PetscCall(OutputVTK(dmCell, "ex11-partition.vtk", &viewer)); 1517 PetscCall(VecView(partition, viewer)); 1518 PetscCall(PetscViewerDestroy(&viewer)); 1519 PetscCall(VecDestroy(&partition)); 1520 PetscCall(DMDestroy(&dmCell)); 1521 } 1522 /* collect max maxspeed from all processes -- todo */ 1523 PetscCall(DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius)); 1524 PetscCall(DMDestroy(&plex)); 1525 PetscCall(MPIU_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 1526 PetscCheck(mod->maxspeed > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set maxspeed", physname); 1527 dt = cfl * minRadius / mod->maxspeed; 1528 PetscCall(TSSetTimeStep(ts, dt)); 1529 PetscCall(TSSetFromOptions(ts)); 1530 1531 /* When using adaptive mesh refinement 1532 specify callbacks to refine the solution 1533 and interpolate data from old to new mesh */ 1534 if (useAMR) { PetscCall(TSSetResize(ts, adaptToleranceFVMSetUp, Transfer, &tctx)); } 1535 PetscCall(TSSetSolution(ts, X)); 1536 PetscCall(VecDestroy(&X)); 1537 PetscCall(TSSolve(ts, NULL)); 1538 PetscCall(TSGetSolveTime(ts, &ftime)); 1539 PetscCall(TSGetStepNumber(ts, &nsteps)); 1540 PetscCall(TSGetConvergedReason(ts, &reason)); 1541 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps)); 1542 PetscCall(TSDestroy(&ts)); 1543 1544 PetscCall(VecTaggerDestroy(&refineTag)); 1545 PetscCall(VecTaggerDestroy(&coarsenTag)); 1546 PetscCall(PetscFunctionListDestroy(&PhysicsList)); 1547 PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_SW)); 1548 PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_Euler)); 1549 PetscCall(FunctionalLinkDestroy(&user->model->functionalRegistry)); 1550 PetscCall(PetscFree(user->model->functionalMonitored)); 1551 PetscCall(PetscFree(user->model->functionalCall)); 1552 PetscCall(PetscFree(user->model->physics->data)); 1553 PetscCall(PetscFree(user->model->physics)); 1554 PetscCall(PetscFree(user->model)); 1555 PetscCall(PetscFree(user)); 1556 PetscCall(PetscLimiterDestroy(&limiter)); 1557 PetscCall(PetscLimiterDestroy(&noneLimiter)); 1558 PetscCall(PetscFVDestroy(&fvm)); 1559 PetscCall(DMDestroy(&dm)); 1560 PetscCall(PetscFinalize()); 1561 return 0; 1562 } 1563 1564 /* Subroutine to set up the initial conditions for the */ 1565 /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */ 1566 /* ----------------------------------------------------------------------- */ 1567 int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3]) 1568 { 1569 int j, k; 1570 /* Wc=matmul(lv,Ueq) 3 vars */ 1571 for (k = 0; k < 3; ++k) { 1572 wc[k] = 0.; 1573 for (j = 0; j < 3; ++j) wc[k] += lv[k][j] * ueq[j]; 1574 } 1575 return 0; 1576 } 1577 /* ----------------------------------------------------------------------- */ 1578 int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3]) 1579 { 1580 int k, j; 1581 /* V=matmul(rv,WC) 3 vars */ 1582 for (k = 0; k < 3; ++k) { 1583 v[k] = 0.; 1584 for (j = 0; j < 3; ++j) v[k] += rv[k][j] * wc[j]; 1585 } 1586 return 0; 1587 } 1588 /* ---------------------------------------------------------------------- */ 1589 int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma) 1590 { 1591 int j, k; 1592 PetscReal rho, csnd, p0; 1593 /* PetscScalar u; */ 1594 1595 for (k = 0; k < 3; ++k) 1596 for (j = 0; j < 3; ++j) { 1597 lv[k][j] = 0.; 1598 rv[k][j] = 0.; 1599 } 1600 rho = ueq[0]; 1601 /* u = ueq[1]; */ 1602 p0 = ueq[2]; 1603 csnd = PetscSqrtReal(gamma * p0 / rho); 1604 lv[0][1] = rho * .5; 1605 lv[0][2] = -.5 / csnd; 1606 lv[1][0] = csnd; 1607 lv[1][2] = -1. / csnd; 1608 lv[2][1] = rho * .5; 1609 lv[2][2] = .5 / csnd; 1610 rv[0][0] = -1. / csnd; 1611 rv[1][0] = 1. / rho; 1612 rv[2][0] = -csnd; 1613 rv[0][1] = 1. / csnd; 1614 rv[0][2] = 1. / csnd; 1615 rv[1][2] = 1. / rho; 1616 rv[2][2] = csnd; 1617 return 0; 1618 } 1619 1620 int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx) 1621 { 1622 PetscReal p0, u0, wcp[3], wc[3]; 1623 PetscReal lv[3][3]; 1624 PetscReal vp[3]; 1625 PetscReal rv[3][3]; 1626 PetscReal eps, ueq[3], rho0, twopi; 1627 1628 /* Function Body */ 1629 twopi = 2. * PETSC_PI; 1630 eps = 1e-4; /* perturbation */ 1631 rho0 = 1e3; /* density of water */ 1632 p0 = 101325.; /* init pressure of 1 atm (?) */ 1633 u0 = 0.; 1634 ueq[0] = rho0; 1635 ueq[1] = u0; 1636 ueq[2] = p0; 1637 /* Project initial state to characteristic variables */ 1638 eigenvectors(rv, lv, ueq, gamma); 1639 projecteqstate(wc, ueq, lv); 1640 wcp[0] = wc[0]; 1641 wcp[1] = wc[1]; 1642 wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx); 1643 projecttoprim(vp, wcp, rv); 1644 ux->r = vp[0]; /* density */ 1645 ux->ru[0] = vp[0] * vp[1]; /* x momentum */ 1646 ux->ru[1] = 0.; 1647 #if defined DIM > 2 1648 if (dim > 2) ux->ru[2] = 0.; 1649 #endif 1650 /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */ 1651 ux->E = vp[2] / (gamma - 1.) + 0.5 * vp[0] * vp[1] * vp[1]; 1652 return 0; 1653 } 1654 1655 /*TEST 1656 1657 testset: 1658 args: -dm_plex_adj_cone -dm_plex_adj_closure 0 1659 1660 test: 1661 suffix: adv_2d_tri_0 1662 requires: triangle 1663 TODO: how did this ever get in main when there is no support for this 1664 args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_faces 1,1 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 1665 1666 test: 1667 suffix: adv_2d_tri_1 1668 requires: triangle 1669 TODO: how did this ever get in main when there is no support for this 1670 args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -dm_plex_faces 1,1 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 1671 1672 test: 1673 suffix: tut_1 1674 requires: exodusii 1675 nsize: 1 1676 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo - 1677 1678 test: 1679 suffix: tut_2 1680 requires: exodusii 1681 nsize: 1 1682 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw 1683 1684 test: 1685 suffix: tut_3 1686 requires: exodusii 1687 nsize: 4 1688 args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin 1689 1690 test: 1691 suffix: tut_4 1692 requires: exodusii 1693 nsize: 4 1694 args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod 1695 1696 testset: 1697 args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 1698 1699 # 2D Advection 0-10 1700 test: 1701 suffix: 0 1702 requires: exodusii 1703 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo 1704 1705 test: 1706 suffix: 1 1707 requires: exodusii 1708 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo 1709 1710 test: 1711 suffix: 2 1712 requires: exodusii 1713 nsize: 2 1714 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo 1715 1716 test: 1717 suffix: 3 1718 requires: exodusii 1719 nsize: 2 1720 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo 1721 1722 test: 1723 suffix: 4 1724 requires: exodusii 1725 nsize: 4 1726 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -petscpartitioner_type simple 1727 1728 test: 1729 suffix: 5 1730 requires: exodusii 1731 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1 1732 1733 test: 1734 suffix: 7 1735 requires: exodusii 1736 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 1737 1738 test: 1739 suffix: 8 1740 requires: exodusii 1741 nsize: 2 1742 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 1743 1744 test: 1745 suffix: 9 1746 requires: exodusii 1747 nsize: 8 1748 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 1749 1750 test: 1751 suffix: 10 1752 requires: exodusii 1753 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo 1754 1755 # 2D Shallow water 1756 testset: 1757 args: -physics sw -ufv_vtk_interval 0 -dm_plex_adj_cone -dm_plex_adj_closure 0 1758 1759 test: 1760 suffix: sw_0 1761 requires: exodusii 1762 args: -bc_wall 100,101 -ufv_cfl 5 -petscfv_type leastsquares -petsclimiter_type sin \ 1763 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo \ 1764 -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ 1765 -monitor height,energy 1766 1767 test: 1768 suffix: sw_ceed 1769 requires: exodusii libceed 1770 args: -sw_riemann rusanov_ceed -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \ 1771 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -dm_plex_use_ceed \ 1772 -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ 1773 -monitor height,energy 1774 1775 test: 1776 suffix: sw_1 1777 nsize: 2 1778 args: -bc_wall 1,3 -ufv_cfl 5 -petsclimiter_type sin \ 1779 -dm_plex_shape annulus -dm_plex_simplex 0 -dm_plex_box_faces 24,12 -dm_plex_box_lower 0,1 -dm_plex_box_upper 1,3 -dm_distribute_overlap 1 \ 1780 -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ 1781 -monitor height,energy 1782 1783 test: 1784 suffix: sw_hll 1785 args: -sw_riemann hll -bc_wall 1,2,3,4 -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin \ 1786 -grid_bounds 0,5,0,5 -dm_plex_simplex 0 -dm_plex_box_faces 25,25 \ 1787 -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ 1788 -monitor height,energy 1789 1790 # 2D Euler 1791 testset: 1792 args: -physics euler -eu_type linear_wave -eu_gamma 1.4 -dm_plex_adj_cone -dm_plex_adj_closure 0 \ 1793 -ufv_vtk_interval 0 -ufv_vtk_basename ${wPETSC_DIR}/ex11 -monitor density,energy 1794 1795 test: 1796 suffix: euler_0 1797 requires: exodusii !complex 1798 args: -eu_riemann godunov -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \ 1799 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo \ 1800 -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 1801 1802 test: 1803 suffix: euler_ceed 1804 requires: exodusii libceed 1805 args: -eu_riemann godunov_ceed -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \ 1806 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -dm_plex_use_ceed \ 1807 -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 1808 1809 testset: 1810 args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 1811 1812 # 2D Advection: p4est 1813 test: 1814 suffix: p4est_advec_2d 1815 requires: p4est 1816 args: -ufv_vtk_interval 0 -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash -dm_forest_maximum_refinement 5 1817 1818 # Advection in a box 1819 test: 1820 suffix: adv_2d_quad_0 1821 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 1822 1823 test: 1824 suffix: adv_2d_quad_1 1825 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 1826 timeoutfactor: 3 1827 1828 test: 1829 suffix: adv_2d_quad_p4est_0 1830 requires: p4est 1831 args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 1832 1833 test: 1834 suffix: adv_2d_quad_p4est_1 1835 requires: p4est 1836 args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 1837 timeoutfactor: 3 1838 1839 test: 1840 suffix: adv_2d_quad_p4est_adapt_0 1841 requires: p4est !__float128 #broken for quad precision 1842 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 -ufv_use_amr -refine_vec_tagger_box 0.005,inf -coarsen_vec_tagger_box 0,1.e-5 -petscfv_type leastsquares -ts_max_time 0.01 1843 timeoutfactor: 3 1844 1845 test: 1846 suffix: adv_0 1847 requires: exodusii 1848 args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201 1849 1850 test: 1851 suffix: shock_0 1852 requires: p4est !single !complex 1853 args: -dm_plex_box_faces 2,1 -grid_bounds -1,1.,0.,1 -grid_skew_60 \ 1854 -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 \ 1855 -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view \ 1856 -bc_wall 1,2,3,4 -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. \ 1857 -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 \ 1858 -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ 1859 -ufv_vtk_basename ${wPETSC_DIR}/ex11 -ufv_vtk_interval 0 -monitor density,energy 1860 timeoutfactor: 3 1861 1862 # Test GLVis visualization of PetscFV fields 1863 test: 1864 suffix: glvis_adv_2d_tet 1865 args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \ 1866 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \ 1867 -ts_monitor_solution glvis: -ts_max_steps 0 1868 1869 test: 1870 suffix: glvis_adv_2d_quad 1871 args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \ 1872 -dm_refine 5 -dm_plex_separate_marker \ 1873 -ts_monitor_solution glvis: -ts_max_steps 0 1874 1875 TEST*/ 1876