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