1 static char help[] = "Hybrid Finite Element-Finite Volume Example.\n"; 2 /*F 3 Here we are advecting a passive tracer in a harmonic velocity field, defined by 4 a forcing function $f$: 5 \begin{align} 6 -\Delta \mathbf{u} + f &= 0 \\ 7 \frac{\partial\phi}{\partial t} + \nabla\cdot \phi \mathbf{u} &= 0 8 \end{align} 9 F*/ 10 11 #include <petscdmplex.h> 12 #include <petscds.h> 13 #include <petscts.h> 14 15 #include <petsc/private/dmpleximpl.h> /* For DotD */ 16 17 typedef enum {VEL_ZERO, VEL_CONSTANT, VEL_HARMONIC, VEL_SHEAR} VelocityDistribution; 18 19 typedef enum {ZERO, CONSTANT, GAUSSIAN, TILTED, DELTA} PorosityDistribution; 20 21 static PetscErrorCode constant_u_2d(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 22 23 /* 24 FunctionalFunc - Calculates the value of a functional of the solution at a point 25 26 Input Parameters: 27 + dm - The DM 28 . time - The TS time 29 . x - The coordinates of the evaluation point 30 . u - The field values at point x 31 - ctx - A user context, or NULL 32 33 Output Parameter: 34 . f - The value of the functional at point x 35 36 */ 37 typedef PetscErrorCode (*FunctionalFunc)(DM, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); 38 39 typedef struct _n_Functional *Functional; 40 struct _n_Functional { 41 char *name; 42 FunctionalFunc func; 43 void *ctx; 44 PetscInt offset; 45 Functional next; 46 }; 47 48 typedef struct { 49 /* Problem definition */ 50 PetscBool useFV; /* Use a finite volume scheme for advection */ 51 PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 52 VelocityDistribution velocityDist; 53 PorosityDistribution porosityDist; 54 PetscReal inflowState; 55 PetscReal source[3]; 56 /* Monitoring */ 57 PetscInt numMonitorFuncs, maxMonitorFunc; 58 Functional *monitorFuncs; 59 PetscInt errorFunctional; 60 Functional functionalRegistry; 61 } AppCtx; 62 63 static AppCtx *globalUser; 64 65 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 66 { 67 const char *velocityDist[4] = {"zero", "constant", "harmonic", "shear"}; 68 const char *porosityDist[5] = {"zero", "constant", "gaussian", "tilted", "delta"}; 69 PetscInt vd, pd, d; 70 PetscBool flg; 71 72 PetscFunctionBeginUser; 73 options->useFV = PETSC_FALSE; 74 options->velocityDist = VEL_HARMONIC; 75 options->porosityDist = ZERO; 76 options->inflowState = -2.0; 77 options->numMonitorFuncs = 0; 78 options->source[0] = 0.5; 79 options->source[1] = 0.5; 80 options->source[2] = 0.5; 81 82 PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX"); 83 PetscCall(PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL)); 84 vd = options->velocityDist; 85 PetscCall(PetscOptionsEList("-velocity_dist","Velocity distribution type","ex18.c",velocityDist,4,velocityDist[options->velocityDist],&vd,NULL)); 86 options->velocityDist = (VelocityDistribution) vd; 87 pd = options->porosityDist; 88 PetscCall(PetscOptionsEList("-porosity_dist","Initial porosity distribution type","ex18.c",porosityDist,5,porosityDist[options->porosityDist],&pd,NULL)); 89 options->porosityDist = (PorosityDistribution) pd; 90 PetscCall(PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL)); 91 d = 2; 92 PetscCall(PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg)); 93 PetscCheck(!flg || d == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %" PetscInt_FMT, d); 94 PetscOptionsEnd(); 95 96 PetscFunctionReturn(0); 97 } 98 99 static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) 100 { 101 Functional func; 102 char *names[256]; 103 PetscInt f; 104 105 PetscFunctionBeginUser; 106 PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX"); 107 options->numMonitorFuncs = PETSC_STATIC_ARRAY_LENGTH(names); 108 PetscCall(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL)); 109 PetscCall(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs)); 110 for (f = 0; f < options->numMonitorFuncs; ++f) { 111 for (func = options->functionalRegistry; func; func = func->next) { 112 PetscBool match; 113 114 PetscCall(PetscStrcasecmp(names[f], func->name, &match)); 115 if (match) break; 116 } 117 PetscCheck(func,comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 118 options->monitorFuncs[f] = func; 119 /* Jed inserts a de-duplication of functionals here */ 120 PetscCall(PetscFree(names[f])); 121 } 122 /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 123 options->maxMonitorFunc = -1; 124 for (func = options->functionalRegistry; func; func = func->next) { 125 for (f = 0; f < options->numMonitorFuncs; ++f) { 126 Functional call = options->monitorFuncs[f]; 127 128 if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 129 } 130 } 131 PetscOptionsEnd(); 132 PetscFunctionReturn(0); 133 } 134 135 static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) 136 { 137 Functional *ptr, f; 138 PetscInt lastoffset = -1; 139 140 PetscFunctionBeginUser; 141 for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 142 PetscCall(PetscNew(&f)); 143 PetscCall(PetscStrallocpy(name, &f->name)); 144 f->offset = lastoffset + 1; 145 f->func = func; 146 f->ctx = ctx; 147 f->next = NULL; 148 *ptr = f; 149 *offset = f->offset; 150 PetscFunctionReturn(0); 151 } 152 153 static PetscErrorCode FunctionalDestroy(Functional *link) 154 { 155 Functional next, l; 156 157 PetscFunctionBeginUser; 158 if (!link) PetscFunctionReturn(0); 159 l = *link; 160 *link = NULL; 161 for (; l; l=next) { 162 next = l->next; 163 PetscCall(PetscFree(l->name)); 164 PetscCall(PetscFree(l)); 165 } 166 PetscFunctionReturn(0); 167 } 168 169 static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 170 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 171 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 172 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 173 { 174 PetscInt comp; 175 for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 176 } 177 178 static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 179 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 180 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 181 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 182 { 183 PetscScalar wind[3] = {0.0, 0.0, 0.0}; 184 PetscInt comp; 185 186 constant_u_2d(dim, t, x, Nf, wind, NULL); 187 for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp]; 188 } 189 190 static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 191 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 192 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 193 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 194 { 195 PetscInt comp; 196 for (comp = 0; comp < dim*dim; ++comp) f1[comp] = 0.0; 197 } 198 199 static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 200 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 201 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 202 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 203 { 204 PetscInt d; 205 for (d = 0; d < dim; ++d) g0[d*dim+d] = 1.0; 206 } 207 208 static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 209 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 210 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 211 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 212 { 213 g0[0] = 1.0; 214 } 215 216 static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 217 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 218 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 219 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 220 { 221 PetscInt comp; 222 for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 223 } 224 225 static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 226 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 227 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 228 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 229 { 230 PetscInt comp, d; 231 for (comp = 0; comp < dim; ++comp) { 232 for (d = 0; d < dim; ++d) { 233 f1[comp*dim+d] = u_x[comp*dim+d]; 234 } 235 } 236 } 237 238 static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 239 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 240 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 241 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 242 { 243 f0[0] = -PetscSinReal(2.0*PETSC_PI*x[0]); 244 f0[1] = 2.0*PETSC_PI*x[1]*PetscCosReal(2.0*PETSC_PI*x[0]); 245 } 246 247 static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 248 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 249 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 250 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 251 { 252 f0[0] = -2.0*PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1]); 253 f0[1] = 2.0*PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0]); 254 } 255 256 void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 257 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 258 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 259 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 260 { 261 const PetscInt Ncomp = dim; 262 PetscInt compI, d; 263 264 for (compI = 0; compI < Ncomp; ++compI) { 265 for (d = 0; d < dim; ++d) { 266 g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0; 267 } 268 } 269 } 270 271 /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 272 static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 273 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 274 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 275 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 276 { 277 PetscInt d; 278 f0[0] = u_t[dim]; 279 for (d = 0; d < dim; ++d) f0[0] += u[dim]*u_x[d*dim+d] + u_x[dim*dim+d]*u[d]; 280 } 281 282 static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 283 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 284 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 285 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 286 { 287 PetscInt d; 288 for (d = 0; d < dim; ++d) f1[0] = 0.0; 289 } 290 291 void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 292 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 293 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 294 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 295 { 296 PetscInt d; 297 g0[0] = u_tShift; 298 for (d = 0; d < dim; ++d) g0[0] += u_x[d*dim+d]; 299 } 300 301 void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 302 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 303 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 304 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 305 { 306 PetscInt d; 307 for (d = 0; d < dim; ++d) g1[d] = u[d]; 308 } 309 310 void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 311 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 312 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 313 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 314 { 315 PetscInt d; 316 for (d = 0; d < dim; ++d) g0[0] += u_x[dim*dim+d]; 317 } 318 319 void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 320 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 321 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 322 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 323 { 324 PetscInt d; 325 for (d = 0; d < dim; ++d) g1[d*dim+d] = u[dim]; 326 } 327 328 static void riemann_advection(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *uL, const PetscScalar *uR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, void *ctx) 329 { 330 PetscReal wind[3] = {0.0, 1.0, 0.0}; 331 PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim,3), wind, n); 332 333 flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 334 } 335 336 static void riemann_coupled_advection(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *uL, const PetscScalar *uR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, void *ctx) 337 { 338 PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 339 340 #if 1 341 flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 342 #else 343 /* if (fabs(uL[0] - wind[0]) > 1.0e-7 || fabs(uL[1] - wind[1]) > 1.0e-7) PetscPrintf(PETSC_COMM_SELF, "wind (%g, %g) uL (%g, %g) uR (%g, %g)\n", wind[0], wind[1], uL[0], uL[1], uR[0], uR[1]); */ 344 /* Smear it out */ 345 flux[0] = 0.5*((uL[dim] + uR[dim]) + (uL[dim] - uR[dim])*tanh(1.0e5*wn)) * wn; 346 #endif 347 } 348 349 static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 350 { 351 u[0] = 0.0; 352 u[1] = 0.0; 353 return 0; 354 } 355 356 static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 357 { 358 u[0] = 0.0; 359 u[1] = 1.0; 360 return 0; 361 } 362 363 /* Coordinates of the point which was at x at t = 0 */ 364 static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 365 { 366 const PetscReal t = *((PetscReal *) ctx); 367 u[0] = x[0]; 368 u[1] = x[1] + t; 369 #if 0 370 PetscCall(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u)); 371 #else 372 u[1] = u[1] - (int) PetscRealPart(u[1]); 373 #endif 374 return 0; 375 } 376 377 /* 378 In 2D we use the exact solution: 379 380 u = x^2 + y^2 381 v = 2 x^2 - 2xy 382 phi = h(x + y + (u + v) t) 383 f_x = f_y = 4 384 385 so that 386 387 -\Delta u + f = <-4, -4> + <4, 4> = 0 388 {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 389 h_t(x + y + (u + v) t) - u . grad phi - phi div u 390 = u h' + v h' - u h_x - v h_y 391 = 0 392 393 We will conserve phi since 394 395 \nabla \cdot u = 2x - 2x = 0 396 397 Also try h((x + ut)^2 + (y + vt)^2), so that 398 399 h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 400 = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 401 = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 402 = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 403 = 0 404 405 */ 406 static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 407 { 408 u[0] = x[0]*x[0] + x[1]*x[1]; 409 u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 410 return 0; 411 } 412 413 /* 414 In 2D we use the exact, periodic solution: 415 416 u = sin(2 pi x)/4 pi^2 417 v = -y cos(2 pi x)/2 pi 418 phi = h(x + y + (u + v) t) 419 f_x = -sin(2 pi x) 420 f_y = 2 pi y cos(2 pi x) 421 422 so that 423 424 -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 425 426 We will conserve phi since 427 428 \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 429 */ 430 static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 431 { 432 u[0] = PetscSinReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 433 u[1] = -x[1]*PetscCosReal(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI); 434 return 0; 435 } 436 437 /* 438 In 2D we use the exact, doubly periodic solution: 439 440 u = sin(2 pi x) cos(2 pi y)/4 pi^2 441 v = -sin(2 pi y) cos(2 pi x)/4 pi^2 442 phi = h(x + y + (u + v) t) 443 f_x = -2sin(2 pi x) cos(2 pi y) 444 f_y = 2sin(2 pi y) cos(2 pi x) 445 446 so that 447 448 -\Delta u + f = <2 sin(2pi x) cos(2pi y), -2 sin(2pi y) cos(2pi x)> + <-2 sin(2pi x) cos(2pi y), 2 sin(2pi y) cos(2pi x)> = 0 449 450 We will conserve phi since 451 452 \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 453 */ 454 static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 455 { 456 u[0] = PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1])/PetscSqr(2.0*PETSC_PI); 457 u[1] = -PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 458 return 0; 459 } 460 461 static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 462 { 463 u[0] = x[1] - 0.5; 464 u[1] = 0.0; 465 return 0; 466 } 467 468 static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 469 { 470 PetscInt d; 471 for (d = 0; d < dim; ++d) u[d] = 0.0; 472 return 0; 473 } 474 475 static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 476 { 477 u[0] = 0.0; 478 return 0; 479 } 480 481 static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 482 { 483 u[0] = 1.0; 484 return 0; 485 } 486 487 static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 488 { 489 PetscReal x0[2]; 490 PetscScalar xn[2]; 491 492 x0[0] = globalUser->source[0]; 493 x0[1] = globalUser->source[1]; 494 constant_x_2d(dim, time, x0, Nf, xn, ctx); 495 { 496 const PetscReal xi = x[0] - PetscRealPart(xn[0]); 497 const PetscReal eta = x[1] - PetscRealPart(xn[1]); 498 const PetscReal r2 = xi*xi + eta*eta; 499 500 u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 501 } 502 return 0; 503 } 504 505 /* 506 Gaussian blob, initially centered on (0.5, 0.5) 507 508 xi = x(t) - x0, eta = y(t) - y0 509 510 where x(t), y(t) are the integral curves of v(t), 511 512 dx/dt . grad f = v . f 513 514 Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 515 516 v0 f_x + w0 f_y = v . f 517 */ 518 static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 519 { 520 const PetscReal x0[2] = {0.5, 0.5}; 521 const PetscReal sigma = 1.0/6.0; 522 PetscScalar xn[2]; 523 524 constant_x_2d(dim, time, x0, Nf, xn, ctx); 525 { 526 /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 527 /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 528 const PetscReal xi = x[0] - PetscRealPart(xn[0]); 529 const PetscReal eta = x[1] - PetscRealPart(xn[1]); 530 const PetscReal r2 = xi*xi + eta*eta; 531 532 u[0] = PetscExpReal(-r2/(2.0*sigma*sigma))/(sigma*PetscSqrtReal(2.0*PETSC_PI)); 533 } 534 return 0; 535 } 536 537 static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 538 { 539 PetscReal x0[3]; 540 const PetscReal wind[3] = {0.0, 1.0, 0.0}; 541 const PetscReal t = *((PetscReal *) ctx); 542 543 DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 544 if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 545 else u[0] = -2.0; /* Inflow state */ 546 return 0; 547 } 548 549 static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 550 { 551 PetscReal ur[3]; 552 PetscReal x0[3]; 553 const PetscReal t = *((PetscReal *) ctx); 554 555 ur[0] = PetscRealPart(u[0]); ur[1] = PetscRealPart(u[1]); ur[2] = PetscRealPart(u[2]); 556 DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 557 if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 558 else u[0] = -2.0; /* Inflow state */ 559 return 0; 560 } 561 562 static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 563 { 564 AppCtx *user = (AppCtx *) ctx; 565 566 PetscFunctionBeginUser; 567 xG[0] = user->inflowState; 568 PetscFunctionReturn(0); 569 } 570 571 static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 572 { 573 PetscFunctionBeginUser; 574 //xG[0] = xI[dim]; 575 xG[0] = xI[2]; 576 PetscFunctionReturn(0); 577 } 578 579 static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 580 { 581 AppCtx *user = (AppCtx *) ctx; 582 PetscInt dim; 583 584 PetscFunctionBeginUser; 585 PetscCall(DMGetDimension(dm, &dim)); 586 switch (user->porosityDist) { 587 case TILTED: 588 if (user->velocityDist == VEL_ZERO) tilted_phi_2d(dim, time, x, 2, u, (void *) &time); 589 else tilted_phi_coupled_2d(dim, time, x, 2, u, (void *) &time); 590 break; 591 case GAUSSIAN: 592 gaussian_phi_2d(dim, time, x, 2, u, (void *) &time); 593 break; 594 case DELTA: 595 delta_phi_2d(dim, time, x, 2, u, (void *) &time); 596 break; 597 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 598 } 599 PetscFunctionReturn(0); 600 } 601 602 static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 603 { 604 AppCtx *user = (AppCtx *) ctx; 605 PetscScalar yexact[3]={0,0,0}; 606 607 PetscFunctionBeginUser; 608 PetscCall(ExactSolution(dm, time, x, yexact, ctx)); 609 f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 610 PetscFunctionReturn(0); 611 } 612 613 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 614 { 615 PetscFunctionBeginUser; 616 PetscCall(DMCreate(comm, dm)); 617 PetscCall(DMSetType(*dm, DMPLEX)); 618 PetscCall(DMSetFromOptions(*dm)); 619 PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 620 PetscFunctionReturn(0); 621 } 622 623 static PetscErrorCode SetupBC(DM dm, AppCtx *user) 624 { 625 PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 626 DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 627 PetscDS prob; 628 DMLabel label; 629 PetscBool check; 630 PetscInt dim, n = 3; 631 const char *prefix; 632 633 PetscFunctionBeginUser; 634 PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 635 PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 636 PetscCall(DMGetDimension(dm, &dim)); 637 /* Set initial guesses and exact solutions */ 638 switch (dim) { 639 case 2: 640 user->initialGuess[0] = initialVelocity; 641 switch(user->porosityDist) { 642 case ZERO: user->initialGuess[1] = zero_phi;break; 643 case CONSTANT: user->initialGuess[1] = constant_phi;break; 644 case GAUSSIAN: user->initialGuess[1] = gaussian_phi_2d;break; 645 case DELTA: user->initialGuess[1] = delta_phi_2d;break; 646 case TILTED: 647 if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 648 else user->initialGuess[1] = tilted_phi_coupled_2d; 649 break; 650 } 651 break; 652 default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 653 } 654 exactFuncs[0] = user->initialGuess[0]; 655 exactFuncs[1] = user->initialGuess[1]; 656 switch (dim) { 657 case 2: 658 switch (user->velocityDist) { 659 case VEL_ZERO: 660 exactFuncs[0] = zero_u_2d; break; 661 case VEL_CONSTANT: 662 exactFuncs[0] = constant_u_2d; break; 663 case VEL_HARMONIC: 664 switch (bdt[0]) { 665 case DM_BOUNDARY_PERIODIC: 666 switch (bdt[1]) { 667 case DM_BOUNDARY_PERIODIC: 668 exactFuncs[0] = doubly_periodic_u_2d; break; 669 default: 670 exactFuncs[0] = periodic_u_2d; break; 671 } 672 break; 673 default: 674 exactFuncs[0] = quadratic_u_2d; break; 675 } 676 break; 677 case VEL_SHEAR: 678 exactFuncs[0] = shear_bc; break; 679 default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 680 } 681 break; 682 default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 683 } 684 { 685 PetscBool isImplicit = PETSC_FALSE; 686 687 PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 688 if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 689 } 690 PetscCall(PetscOptionsHasName(NULL,NULL, "-dmts_check", &check)); 691 if (check) { 692 user->initialGuess[0] = exactFuncs[0]; 693 user->initialGuess[1] = exactFuncs[1]; 694 } 695 /* Set BC */ 696 PetscCall(DMGetDS(dm, &prob)); 697 PetscCall(DMGetLabel(dm, "marker", &label)); 698 PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user)); 699 PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user)); 700 if (label) { 701 const PetscInt id = 1; 702 703 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL)); 704 } 705 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 706 if (label && user->useFV) { 707 const PetscInt inflowids[] = {100,200,300}, outflowids[] = {101}; 708 709 PetscCall(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "inflow", label, PETSC_STATIC_ARRAY_LENGTH(inflowids), inflowids, 1, 0, NULL, (void (*)(void)) advect_inflow, NULL, user, NULL)); 710 PetscCall(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "outflow", label, PETSC_STATIC_ARRAY_LENGTH(outflowids), outflowids, 1, 0, NULL, (void (*)(void)) advect_outflow, NULL, user, NULL)); 711 } 712 PetscFunctionReturn(0); 713 } 714 715 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 716 { 717 DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 718 PetscDS prob; 719 PetscInt n = 3; 720 const char *prefix; 721 722 PetscFunctionBeginUser; 723 PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 724 PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 725 PetscCall(DMGetDS(dm, &prob)); 726 switch (user->velocityDist) { 727 case VEL_ZERO: 728 PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u)); 729 break; 730 case VEL_CONSTANT: 731 PetscCall(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u)); 732 PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL)); 733 PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL)); 734 break; 735 case VEL_HARMONIC: 736 switch (bdt[0]) { 737 case DM_BOUNDARY_PERIODIC: 738 switch (bdt[1]) { 739 case DM_BOUNDARY_PERIODIC: 740 PetscCall(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u)); 741 break; 742 default: 743 PetscCall(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u)); 744 break; 745 } 746 break; 747 default: 748 PetscCall(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u)); 749 break; 750 } 751 PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 752 break; 753 case VEL_SHEAR: 754 PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u)); 755 PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 756 break; 757 } 758 PetscCall(PetscDSSetResidual(prob, 1, f0_advection, f1_advection)); 759 PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL)); 760 PetscCall(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL)); 761 if (user->velocityDist == VEL_ZERO) PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_advection)); 762 else PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection)); 763 764 PetscCall(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user)); 765 PetscFunctionReturn(0); 766 } 767 768 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 769 { 770 DM cdm = dm; 771 PetscQuadrature q; 772 PetscFE fe[2]; 773 PetscFV fv; 774 MPI_Comm comm; 775 PetscInt dim; 776 777 PetscFunctionBeginUser; 778 /* Create finite element */ 779 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 780 PetscCall(DMGetDimension(dm, &dim)); 781 PetscCall(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0])); 782 PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity")); 783 PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1])); 784 PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 785 PetscCall(PetscObjectSetName((PetscObject) fe[1], "porosity")); 786 787 PetscCall(PetscFVCreate(PetscObjectComm((PetscObject) dm), &fv)); 788 PetscCall(PetscObjectSetName((PetscObject) fv, "porosity")); 789 PetscCall(PetscFVSetFromOptions(fv)); 790 PetscCall(PetscFVSetNumComponents(fv, 1)); 791 PetscCall(PetscFVSetSpatialDimension(fv, dim)); 792 PetscCall(PetscFEGetQuadrature(fe[0], &q)); 793 PetscCall(PetscFVSetQuadrature(fv, q)); 794 795 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 796 if (user->useFV) PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fv)); 797 else PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 798 PetscCall(DMCreateDS(dm)); 799 PetscCall(SetupProblem(dm, user)); 800 801 /* Set discretization and boundary conditions for each mesh */ 802 while (cdm) { 803 PetscCall(DMCopyDisc(dm, cdm)); 804 PetscCall(DMGetCoarseDM(cdm, &cdm)); 805 /* Coordinates were never localized for coarse meshes */ 806 if (cdm) PetscCall(DMLocalizeCoordinates(cdm)); 807 } 808 PetscCall(PetscFEDestroy(&fe[0])); 809 PetscCall(PetscFEDestroy(&fe[1])); 810 PetscCall(PetscFVDestroy(&fv)); 811 PetscFunctionReturn(0); 812 } 813 814 static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 815 { 816 PetscFunctionBeginUser; 817 PetscCall(CreateMesh(comm, user, dm)); 818 /* Handle refinement, etc. */ 819 PetscCall(DMSetFromOptions(*dm)); 820 /* Construct ghost cells */ 821 if (user->useFV) { 822 DM gdm; 823 824 PetscCall(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm)); 825 PetscCall(DMDestroy(dm)); 826 *dm = gdm; 827 } 828 /* Localize coordinates */ 829 PetscCall(DMLocalizeCoordinates(*dm)); 830 PetscCall(PetscObjectSetName((PetscObject)(*dm),"Mesh")); 831 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 832 /* Setup problem */ 833 PetscCall(SetupDiscretization(*dm, user)); 834 /* Setup BC */ 835 PetscCall(SetupBC(*dm, user)); 836 PetscFunctionReturn(0); 837 } 838 839 static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx) 840 { 841 PetscDS prob; 842 DM dmCell; 843 Vec cellgeom; 844 const PetscScalar *cgeom; 845 PetscScalar *x; 846 PetscInt dim, Nf, cStart, cEnd, c; 847 848 PetscFunctionBeginUser; 849 PetscCall(DMGetDS(dm, &prob)); 850 PetscCall(DMGetDimension(dm, &dim)); 851 PetscCall(PetscDSGetNumFields(prob, &Nf)); 852 PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 853 PetscCall(VecGetDM(cellgeom, &dmCell)); 854 PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 855 PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 856 PetscCall(VecGetArray(X, &x)); 857 for (c = cStart; c < cEnd; ++c) { 858 PetscFVCellGeom *cg; 859 PetscScalar *xc; 860 861 PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 862 PetscCall(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc)); 863 if (xc) PetscCall((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx)); 864 } 865 PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 866 PetscCall(VecRestoreArray(X, &x)); 867 PetscFunctionReturn(0); 868 } 869 870 static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 871 { 872 AppCtx *user = (AppCtx *) ctx; 873 char *ftable = NULL; 874 DM dm; 875 PetscSection s; 876 Vec cellgeom; 877 const PetscScalar *x; 878 PetscScalar *a; 879 PetscReal *xnorms; 880 PetscInt pStart, pEnd, p, Nf, f; 881 882 PetscFunctionBeginUser; 883 PetscCall(VecViewFromOptions(X, (PetscObject) ts, "-view_solution")); 884 PetscCall(VecGetDM(X, &dm)); 885 PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 886 PetscCall(DMGetLocalSection(dm, &s)); 887 PetscCall(PetscSectionGetNumFields(s, &Nf)); 888 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 889 PetscCall(PetscCalloc1(Nf*2, &xnorms)); 890 PetscCall(VecGetArrayRead(X, &x)); 891 for (p = pStart; p < pEnd; ++p) { 892 for (f = 0; f < Nf; ++f) { 893 PetscInt dof, cdof, d; 894 895 PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 896 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 897 PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, x, &a)); 898 /* TODO Use constrained indices here */ 899 for (d = 0; d < dof-cdof; ++d) xnorms[f*2+0] = PetscMax(xnorms[f*2+0], PetscAbsScalar(a[d])); 900 for (d = 0; d < dof-cdof; ++d) xnorms[f*2+1] += PetscAbsScalar(a[d]); 901 } 902 } 903 PetscCall(VecRestoreArrayRead(X, &x)); 904 if (stepnum >= 0) { /* No summary for final time */ 905 DM dmCell, *fdm; 906 Vec *fv; 907 const PetscScalar *cgeom; 908 PetscScalar **fx; 909 PetscReal *fmin, *fmax, *fint, *ftmp, t; 910 PetscInt cStart, cEnd, c, fcount, f, num; 911 912 size_t ftableused,ftablealloc; 913 914 /* Functionals have indices after registering, this is an upper bound */ 915 fcount = user->numMonitorFuncs; 916 PetscCall(PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fint,fcount,&ftmp)); 917 PetscCall(PetscMalloc3(fcount,&fdm,fcount,&fv,fcount,&fx)); 918 for (f = 0; f < fcount; ++f) { 919 PetscSection fs; 920 const char *name = user->monitorFuncs[f]->name; 921 922 fmin[f] = PETSC_MAX_REAL; 923 fmax[f] = PETSC_MIN_REAL; 924 fint[f] = 0; 925 /* Make monitor vecs */ 926 PetscCall(DMClone(dm, &fdm[f])); 927 PetscCall(DMGetOutputSequenceNumber(dm, &num, &t)); 928 PetscCall(DMSetOutputSequenceNumber(fdm[f], num, t)); 929 PetscCall(PetscSectionClone(s, &fs)); 930 PetscCall(PetscSectionSetFieldName(fs, 0, NULL)); 931 PetscCall(PetscSectionSetFieldName(fs, 1, name)); 932 PetscCall(DMSetLocalSection(fdm[f], fs)); 933 PetscCall(PetscSectionDestroy(&fs)); 934 PetscCall(DMGetGlobalVector(fdm[f], &fv[f])); 935 PetscCall(PetscObjectSetName((PetscObject) fv[f], name)); 936 PetscCall(VecGetArray(fv[f], &fx[f])); 937 } 938 PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 939 PetscCall(VecGetDM(cellgeom, &dmCell)); 940 PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 941 PetscCall(VecGetArrayRead(X, &x)); 942 for (c = cStart; c < cEnd; ++c) { 943 PetscFVCellGeom *cg; 944 PetscScalar *cx; 945 946 PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 947 PetscCall(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx)); 948 if (!cx) continue; /* not a global cell */ 949 for (f = 0; f < user->numMonitorFuncs; ++f) { 950 Functional func = user->monitorFuncs[f]; 951 PetscScalar *fxc; 952 953 PetscCall(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc)); 954 /* I need to make it easier to get interpolated values here */ 955 PetscCall((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx)); 956 fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 957 } 958 for (f = 0; f < fcount; ++f) { 959 fmin[f] = PetscMin(fmin[f], ftmp[f]); 960 fmax[f] = PetscMax(fmax[f], ftmp[f]); 961 fint[f] += cg->volume * ftmp[f]; 962 } 963 } 964 PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 965 PetscCall(VecRestoreArrayRead(X, &x)); 966 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 967 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 968 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 969 /* Output functional data */ 970 ftablealloc = fcount * 100; 971 ftableused = 0; 972 PetscCall(PetscCalloc1(ftablealloc, &ftable)); 973 for (f = 0; f < user->numMonitorFuncs; ++f) { 974 Functional func = user->monitorFuncs[f]; 975 PetscInt id = func->offset; 976 char newline[] = "\n"; 977 char buffer[256], *p, *prefix; 978 size_t countused, len; 979 980 /* Create string with functional outputs */ 981 if (f % 3) { 982 PetscCall(PetscArraycpy(buffer, " ", 2)); 983 p = buffer + 2; 984 } else if (f) { 985 PetscCall(PetscArraycpy(buffer, newline, sizeof(newline)-1)); 986 p = buffer + sizeof(newline) - 1; 987 } else { 988 p = buffer; 989 } 990 PetscCall(PetscSNPrintfCount(p, sizeof buffer-(p-buffer), "%12s [%12.6g,%12.6g] int %12.6g", &countused, func->name, (double) fmin[id], (double) fmax[id], (double) fint[id])); 991 countused += p - buffer; 992 /* reallocate */ 993 if (countused > ftablealloc-ftableused-1) { 994 char *ftablenew; 995 996 ftablealloc = 2*ftablealloc + countused; 997 PetscCall(PetscMalloc1(ftablealloc, &ftablenew)); 998 PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); 999 PetscCall(PetscFree(ftable)); 1000 ftable = ftablenew; 1001 } 1002 PetscCall(PetscArraycpy(ftable+ftableused, buffer, countused)); 1003 ftableused += countused; 1004 ftable[ftableused] = 0; 1005 /* Output vecs */ 1006 PetscCall(VecRestoreArray(fv[f], &fx[f])); 1007 PetscCall(PetscStrlen(func->name, &len)); 1008 PetscCall(PetscMalloc1(len+2,&prefix)); 1009 PetscCall(PetscStrcpy(prefix, func->name)); 1010 PetscCall(PetscStrcat(prefix, "_")); 1011 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix)); 1012 PetscCall(VecViewFromOptions(fv[f], NULL, "-vec_view")); 1013 PetscCall(PetscFree(prefix)); 1014 PetscCall(DMRestoreGlobalVector(fdm[f], &fv[f])); 1015 PetscCall(DMDestroy(&fdm[f])); 1016 } 1017 PetscCall(PetscFree4(fmin, fmax, fint, ftmp)); 1018 PetscCall(PetscFree3(fdm, fv, fx)); 1019 PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "% 3" PetscInt_FMT " time %8.4g |x| (", stepnum, (double) time)); 1020 for (f = 0; f < Nf; ++f) { 1021 if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 1022 PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+0])); 1023 } 1024 PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ") |x|_1 (")); 1025 for (f = 0; f < Nf; ++f) { 1026 if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 1027 PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+1])); 1028 } 1029 PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ") %s\n", ftable ? ftable : "")); 1030 PetscCall(PetscFree(ftable)); 1031 } 1032 PetscCall(PetscFree(xnorms)); 1033 PetscFunctionReturn(0); 1034 } 1035 1036 int main(int argc, char **argv) 1037 { 1038 MPI_Comm comm; 1039 TS ts; 1040 DM dm; 1041 Vec u; 1042 AppCtx user; 1043 PetscReal t0, t = 0.0; 1044 void *ctxs[2]; 1045 1046 ctxs[0] = &t; 1047 ctxs[1] = &t; 1048 PetscFunctionBeginUser; 1049 PetscCall(PetscInitialize(&argc, &argv, (char*) 0, help)); 1050 comm = PETSC_COMM_WORLD; 1051 user.functionalRegistry = NULL; 1052 globalUser = &user; 1053 PetscCall(ProcessOptions(comm, &user)); 1054 PetscCall(TSCreate(comm, &ts)); 1055 PetscCall(TSSetType(ts, TSBEULER)); 1056 PetscCall(CreateDM(comm, &user, &dm)); 1057 PetscCall(TSSetDM(ts, dm)); 1058 PetscCall(ProcessMonitorOptions(comm, &user)); 1059 1060 PetscCall(DMCreateGlobalVector(dm, &u)); 1061 PetscCall(PetscObjectSetName((PetscObject) u, "solution")); 1062 if (user.useFV) { 1063 PetscBool isImplicit = PETSC_FALSE; 1064 1065 PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 1066 if (isImplicit) { 1067 PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 1068 PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1069 } 1070 PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 1071 PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user)); 1072 } else { 1073 PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 1074 PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 1075 PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1076 } 1077 if (user.useFV) PetscCall(TSMonitorSet(ts, MonitorFunctionals, &user, NULL)); 1078 PetscCall(TSSetMaxSteps(ts, 1)); 1079 PetscCall(TSSetMaxTime(ts, 2.0)); 1080 PetscCall(TSSetTimeStep(ts,0.01)); 1081 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 1082 PetscCall(TSSetFromOptions(ts)); 1083 1084 PetscCall(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u)); 1085 if (user.useFV) PetscCall(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1])); 1086 PetscCall(VecViewFromOptions(u, NULL, "-init_vec_view")); 1087 PetscCall(TSGetTime(ts, &t)); 1088 t0 = t; 1089 PetscCall(DMTSCheckFromOptions(ts, u)); 1090 PetscCall(TSSolve(ts, u)); 1091 PetscCall(TSGetTime(ts, &t)); 1092 if (t > t0) PetscCall(DMTSCheckFromOptions(ts, u)); 1093 PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 1094 { 1095 PetscReal ftime; 1096 PetscInt nsteps; 1097 TSConvergedReason reason; 1098 1099 PetscCall(TSGetSolveTime(ts, &ftime)); 1100 PetscCall(TSGetStepNumber(ts, &nsteps)); 1101 PetscCall(TSGetConvergedReason(ts, &reason)); 1102 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double) ftime, nsteps)); 1103 } 1104 1105 PetscCall(VecDestroy(&u)); 1106 PetscCall(DMDestroy(&dm)); 1107 PetscCall(TSDestroy(&ts)); 1108 PetscCall(PetscFree(user.monitorFuncs)); 1109 PetscCall(FunctionalDestroy(&user.functionalRegistry)); 1110 PetscCall(PetscFinalize()); 1111 return 0; 1112 } 1113 1114 /*TEST 1115 1116 testset: 1117 args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 1118 1119 # 2D harmonic velocity, no porosity 1120 test: 1121 suffix: p1p1 1122 requires: !complex !single 1123 args: -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_factor_shift_type nonzero -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1124 test: 1125 suffix: p1p1_xper 1126 requires: !complex !single 1127 args: -dm_refine 1 -dm_plex_box_bd periodic,none -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1128 test: 1129 suffix: p1p1_xper_ref 1130 requires: !complex !single 1131 args: -dm_refine 2 -dm_plex_box_bd periodic,none -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1132 test: 1133 suffix: p1p1_xyper 1134 requires: !complex !single 1135 args: -dm_refine 1 -dm_plex_box_bd periodic,periodic -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1136 test: 1137 suffix: p1p1_xyper_ref 1138 requires: !complex !single 1139 args: -dm_refine 2 -dm_plex_box_bd periodic,periodic -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1140 test: 1141 suffix: p2p1 1142 requires: !complex !single 1143 args: -velocity_petscspace_degree 2 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1144 test: 1145 suffix: p2p1_xyper 1146 requires: !complex !single 1147 args: -dm_refine 1 -dm_plex_box_bd periodic,periodic -velocity_petscspace_degree 2 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -pc_factor_shift_type nonzero -ksp_rtol 1.0e-8 -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1148 1149 test: 1150 suffix: adv_1 1151 requires: !complex !single 1152 args: -use_fv -velocity_dist zero -porosity_dist tilted -ts_type ssp -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 1,2,4 -bc_outflow 3 -ts_view -dm_view 1153 1154 test: 1155 suffix: adv_2 1156 requires: !complex 1157 TODO: broken memory corruption 1158 args: -use_fv -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 3,4 -bc_outflow 1,2 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it 100 -ts_view -dm_view -snes_converged_reason -ksp_converged_reason 1159 1160 test: 1161 suffix: adv_3 1162 requires: !complex 1163 TODO: broken memory corruption 1164 args: -dm_plex_box_bd periodic,none -use_fv -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 3 -bc_outflow 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it 100 -ts_view -dm_view -snes_converged_reason 1165 1166 test: 1167 suffix: adv_3_ex 1168 requires: !complex 1169 args: -dm_plex_box_bd periodic,none -use_fv -velocity_dist zero -porosity_dist tilted -ts_type ssp -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.1 -bc_inflow 3 -bc_outflow 1 -snes_fd_color -ksp_max_it 100 -ts_view -dm_view 1170 1171 test: 1172 suffix: adv_4 1173 requires: !complex 1174 TODO: broken memory corruption 1175 args: -use_fv -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 3 -bc_outflow 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it 100 -ts_view -dm_view -snes_converged_reason 1176 1177 # 2D Advection, box, delta 1178 test: 1179 suffix: adv_delta_yper_0 1180 requires: !complex 1181 TODO: broken 1182 args: -dm_plex_box_bd none,periodic -use_fv -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type euler -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -dm_view -monitor Error 1183 1184 test: 1185 suffix: adv_delta_yper_1 1186 requires: !complex 1187 TODO: broken 1188 args: -dm_plex_box_bd none,periodic -use_fv -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type euler -ts_max_time 5.0 -ts_max_steps 40 -ts_dt 0.166666 -bc_inflow 2 -bc_outflow 4 -ts_view -dm_view -monitor Error -dm_refine 1 -source_loc 0.416666,0.416666 1189 1190 test: 1191 suffix: adv_delta_yper_2 1192 requires: !complex 1193 TODO: broken 1194 args: -dm_plex_box_bd none,periodic -use_fv -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type euler -ts_max_time 5.0 -ts_max_steps 80 -ts_dt 0.083333 -bc_inflow 2 -bc_outflow 4 -ts_view -dm_view -monitor Error -dm_refine 2 -source_loc 0.458333,0.458333 1195 1196 test: 1197 suffix: adv_delta_yper_fim_0 1198 requires: !complex 1199 TODO: broken 1200 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_greedy_symmetric 0 -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason 1201 1202 test: 1203 suffix: adv_delta_yper_fim_1 1204 requires: !complex 1205 TODO: broken 1206 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_greedy_symmetric 0 -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -snes_linesearch_type basic 1207 1208 test: 1209 suffix: adv_delta_yper_fim_2 1210 requires: !complex 1211 TODO: broken 1212 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -mat_coloring_greedy_symmetric 0 -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -snes_linesearch_type basic 1213 1214 test: 1215 suffix: adv_delta_yper_im_0 1216 requires: !complex 1217 TODO: broken 1218 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason 1219 1220 test: 1221 suffix: adv_delta_yper_im_1 1222 requires: !complex 1223 TODO: broken 1224 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 40 -ts_dt 0.166666 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 1 -source_loc 0.416666,0.416666 1225 1226 test: 1227 suffix: adv_delta_yper_im_2 1228 requires: !complex 1229 TODO: broken 1230 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 80 -ts_dt 0.083333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 2 -source_loc 0.458333,0.458333 1231 1232 test: 1233 suffix: adv_delta_yper_im_3 1234 requires: !complex 1235 TODO: broken 1236 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason 1237 1238 # I believe the nullspace is sin(pi y) 1239 test: 1240 suffix: adv_delta_yper_im_4 1241 requires: !complex 1242 TODO: broken 1243 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 40 -ts_dt 0.166666 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 1 -source_loc 0.416666,0.416666 1244 1245 test: 1246 suffix: adv_delta_yper_im_5 1247 requires: !complex 1248 TODO: broken 1249 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_mimex_version 0 -ts_max_time 5.0 -ts_max_steps 80 -ts_dt 0.083333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -dm_refine 2 -source_loc 0.458333,0.458333 1250 1251 test: 1252 suffix: adv_delta_yper_im_6 1253 requires: !complex 1254 TODO: broken 1255 args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 2 -bc_outflow 4 -ts_view -monitor Error -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type svd -snes_converged_reason 1256 # 2D Advection, magma benchmark 1 1257 1258 test: 1259 suffix: adv_delta_shear_im_0 1260 requires: !complex 1261 TODO: broken 1262 args: -dm_plex_box_bd periodic,none -dm_refine 2 -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist shear -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_max_time 5.0 -ts_max_steps 20 -ts_dt 0.333333 -bc_inflow 1,3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 -pc_type lu -snes_converged_reason -source_loc 0.458333,0.708333 1263 # 2D Advection, box, gaussian 1264 1265 test: 1266 suffix: adv_gauss 1267 requires: !complex 1268 TODO: broken 1269 args: -use_fv -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type ssp -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view 1270 1271 test: 1272 suffix: adv_gauss_im 1273 requires: !complex 1274 TODO: broken 1275 args: -use_fv -use_implicit -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type beuler -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 1276 1277 test: 1278 suffix: adv_gauss_im_1 1279 requires: !complex 1280 TODO: broken 1281 args: -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type beuler -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 1282 1283 test: 1284 suffix: adv_gauss_im_2 1285 requires: !complex 1286 TODO: broken 1287 args: -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type beuler -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -snes_rtol 1.0e-7 1288 1289 test: 1290 suffix: adv_gauss_corner 1291 requires: !complex 1292 TODO: broken 1293 args: -use_fv -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type ssp -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 2 -ts_view -dm_view 1294 1295 # 2D Advection+Harmonic 12- 1296 test: 1297 suffix: adv_harm_0 1298 requires: !complex 1299 TODO: broken memory corruption 1300 args: -velocity_petscspace_degree 2 -use_fv -velocity_dist harmonic -porosity_dist gaussian -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -bc_inflow 1,2,4 -bc_outflow 3 -use_implicit -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -ksp_max_it 100 -ts_view -dm_view -snes_converged_reason -ksp_converged_reason -snes_monitor -dmts_check 1301 1302 # Must check that FV BCs propagate to coarse meshes 1303 # Must check that FV BC ids propagate to coarse meshes 1304 # Must check that FE+FV BCs work at the same time 1305 # 2D Advection, matching wind in ex11 8-11 1306 # NOTE implicit solves are limited by accuracy of FD Jacobian 1307 test: 1308 suffix: adv_0 1309 requires: !complex !single exodusii 1310 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -velocity_dist zero -porosity_dist tilted -ts_type ssp -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view 1311 1312 test: 1313 suffix: adv_0_im 1314 requires: !complex exodusii 1315 TODO: broken memory corruption 1316 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_dist zero -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu 1317 1318 test: 1319 suffix: adv_0_im_2 1320 requires: !complex exodusii 1321 TODO: broken 1322 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_dist constant -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type lu -snes_rtol 1.0e-7 1323 1324 test: 1325 suffix: adv_0_im_3 1326 requires: !complex exodusii 1327 TODO: broken 1328 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type svd -snes_rtol 1.0e-7 1329 1330 test: 1331 suffix: adv_0_im_4 1332 requires: !complex exodusii 1333 TODO: broken 1334 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -use_fv -use_implicit -velocity_petscspace_degree 2 -velocity_dist constant -porosity_dist tilted -ts_type beuler -ts_max_time 2.0 -ts_max_steps 1000 -ts_dt 0.993392 -ts_view -dm_view -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type svd -snes_rtol 1.0e-7 1335 # 2D Advection, misc 1336 1337 TEST*/ 1338