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 break; 700 default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim); 701 } 702 break; 703 default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %D not supported", dim); 704 } 705 { 706 PetscBool isImplicit = PETSC_FALSE; 707 708 ierr = PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit);CHKERRQ(ierr); 709 if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 710 } 711 ierr = PetscOptionsHasName(NULL,NULL, "-dmts_check", &check);CHKERRQ(ierr); 712 if (check) { 713 user->initialGuess[0] = exactFuncs[0]; 714 user->initialGuess[1] = exactFuncs[1]; 715 } 716 /* Set BC */ 717 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 718 ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 719 ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], user);CHKERRQ(ierr); 720 ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], user);CHKERRQ(ierr); 721 if (label) { 722 const PetscInt id = 1; 723 724 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL);CHKERRQ(ierr); 725 } 726 ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr); 727 if (label && user->useFV) { 728 const PetscInt inflowids[] = {100,200,300}, outflowids[] = {101}; 729 730 ierr = DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "inflow", label, ALEN(inflowids), inflowids, 1, 0, NULL, (void (*)(void)) advect_inflow, NULL, user, NULL);CHKERRQ(ierr); 731 ierr = DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "outflow", label, ALEN(outflowids), outflowids, 1, 0, NULL, (void (*)(void)) advect_outflow, NULL, user, NULL);CHKERRQ(ierr); 732 } 733 PetscFunctionReturn(0); 734 } 735 736 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 737 { 738 DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 739 PetscDS prob; 740 PetscInt n = 3; 741 const char *prefix; 742 PetscErrorCode ierr; 743 744 PetscFunctionBeginUser; 745 ierr = PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix);CHKERRQ(ierr); 746 ierr = PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL);CHKERRQ(ierr); 747 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 748 switch (user->velocityDist) { 749 case VEL_ZERO: 750 ierr = PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u);CHKERRQ(ierr); 751 break; 752 case VEL_CONSTANT: 753 ierr = PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u);CHKERRQ(ierr); 754 ierr = PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL);CHKERRQ(ierr); 755 ierr = PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL);CHKERRQ(ierr); 756 break; 757 case VEL_HARMONIC: 758 switch (bdt[0]) { 759 case DM_BOUNDARY_PERIODIC: 760 switch (bdt[1]) { 761 case DM_BOUNDARY_PERIODIC: 762 ierr = PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u);CHKERRQ(ierr); 763 break; 764 default: 765 ierr = PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u);CHKERRQ(ierr); 766 break; 767 } 768 break; 769 default: 770 ierr = PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u);CHKERRQ(ierr); 771 break; 772 } 773 ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 774 break; 775 case VEL_SHEAR: 776 ierr = PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u);CHKERRQ(ierr); 777 ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 778 break; 779 } 780 ierr = PetscDSSetResidual(prob, 1, f0_advection, f1_advection);CHKERRQ(ierr); 781 ierr = PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL);CHKERRQ(ierr); 782 ierr = PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL);CHKERRQ(ierr); 783 if (user->velocityDist == VEL_ZERO) {ierr = PetscDSSetRiemannSolver(prob, 1, riemann_advection);CHKERRQ(ierr);} 784 else {ierr = PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection);CHKERRQ(ierr);} 785 786 ierr = FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user);CHKERRQ(ierr); 787 PetscFunctionReturn(0); 788 } 789 790 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 791 { 792 DM cdm = dm; 793 PetscQuadrature q; 794 PetscFE fe[2]; 795 PetscFV fv; 796 MPI_Comm comm; 797 PetscInt dim; 798 PetscErrorCode ierr; 799 800 PetscFunctionBeginUser; 801 /* Create finite element */ 802 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 803 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 804 ierr = PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 805 ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 806 ierr = PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 807 ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 808 ierr = PetscObjectSetName((PetscObject) fe[1], "porosity");CHKERRQ(ierr); 809 810 ierr = PetscFVCreate(PetscObjectComm((PetscObject) dm), &fv);CHKERRQ(ierr); 811 ierr = PetscObjectSetName((PetscObject) fv, "porosity");CHKERRQ(ierr); 812 ierr = PetscFVSetFromOptions(fv);CHKERRQ(ierr); 813 ierr = PetscFVSetNumComponents(fv, 1);CHKERRQ(ierr); 814 ierr = PetscFVSetSpatialDimension(fv, dim);CHKERRQ(ierr); 815 ierr = PetscFEGetQuadrature(fe[0], &q);CHKERRQ(ierr); 816 ierr = PetscFVSetQuadrature(fv, q);CHKERRQ(ierr); 817 818 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 819 if (user->useFV) {ierr = DMSetField(dm, 1, NULL, (PetscObject) fv);CHKERRQ(ierr);} 820 else {ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);} 821 ierr = DMCreateDS(dm);CHKERRQ(ierr); 822 ierr = SetupProblem(dm, user);CHKERRQ(ierr); 823 824 /* Set discretization and boundary conditions for each mesh */ 825 while (cdm) { 826 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 827 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 828 /* Coordinates were never localized for coarse meshes */ 829 if (cdm) {ierr = DMLocalizeCoordinates(cdm);CHKERRQ(ierr);} 830 } 831 ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 832 ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 833 ierr = PetscFVDestroy(&fv);CHKERRQ(ierr); 834 PetscFunctionReturn(0); 835 } 836 837 static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 838 { 839 PetscErrorCode ierr; 840 841 PetscFunctionBeginUser; 842 ierr = CreateMesh(comm, user, dm);CHKERRQ(ierr); 843 /* Handle refinement, etc. */ 844 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 845 /* Construct ghost cells */ 846 if (user->useFV) { 847 DM gdm; 848 849 ierr = DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm);CHKERRQ(ierr); 850 ierr = DMDestroy(dm);CHKERRQ(ierr); 851 *dm = gdm; 852 } 853 /* Localize coordinates */ 854 ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 855 ierr = PetscObjectSetName((PetscObject)(*dm),"Mesh");CHKERRQ(ierr); 856 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 857 /* Setup problem */ 858 ierr = SetupDiscretization(*dm, user);CHKERRQ(ierr); 859 /* Setup BC */ 860 ierr = SetupBC(*dm, user);CHKERRQ(ierr); 861 PetscFunctionReturn(0); 862 } 863 864 static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx) 865 { 866 PetscDS prob; 867 DM dmCell; 868 Vec cellgeom; 869 const PetscScalar *cgeom; 870 PetscScalar *x; 871 PetscInt dim, Nf, cStart, cEnd, c; 872 PetscErrorCode ierr; 873 874 PetscFunctionBeginUser; 875 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 876 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 877 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 878 ierr = DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL);CHKERRQ(ierr); 879 ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr); 880 ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 881 ierr = VecGetArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 882 ierr = VecGetArray(X, &x);CHKERRQ(ierr); 883 for (c = cStart; c < cEnd; ++c) { 884 PetscFVCellGeom *cg; 885 PetscScalar *xc; 886 887 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 888 ierr = DMPlexPointGlobalFieldRef(dm, c, field, x, &xc);CHKERRQ(ierr); 889 if (xc) {ierr = (*func)(dim, 0.0, cg->centroid, Nf, xc, ctx);CHKERRQ(ierr);} 890 } 891 ierr = VecRestoreArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 892 ierr = VecRestoreArray(X, &x);CHKERRQ(ierr); 893 PetscFunctionReturn(0); 894 } 895 896 static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 897 { 898 AppCtx *user = (AppCtx *) ctx; 899 char *ftable = NULL; 900 DM dm; 901 PetscSection s; 902 Vec cellgeom; 903 const PetscScalar *x; 904 PetscScalar *a; 905 PetscReal *xnorms; 906 PetscInt pStart, pEnd, p, Nf, f; 907 PetscErrorCode ierr; 908 909 PetscFunctionBeginUser; 910 ierr = VecViewFromOptions(X, (PetscObject) ts, "-view_solution");CHKERRQ(ierr); 911 ierr = VecGetDM(X, &dm);CHKERRQ(ierr); 912 ierr = DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL);CHKERRQ(ierr); 913 ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 914 ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr); 915 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 916 ierr = PetscCalloc1(Nf*2, &xnorms);CHKERRQ(ierr); 917 ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 918 for (p = pStart; p < pEnd; ++p) { 919 for (f = 0; f < Nf; ++f) { 920 PetscInt dof, cdof, d; 921 922 ierr = PetscSectionGetFieldDof(s, p, f, &dof);CHKERRQ(ierr); 923 ierr = PetscSectionGetFieldConstraintDof(s, p, f, &cdof);CHKERRQ(ierr); 924 ierr = DMPlexPointGlobalFieldRead(dm, p, f, x, &a);CHKERRQ(ierr); 925 /* TODO Use constrained indices here */ 926 for (d = 0; d < dof-cdof; ++d) xnorms[f*2+0] = PetscMax(xnorms[f*2+0], PetscAbsScalar(a[d])); 927 for (d = 0; d < dof-cdof; ++d) xnorms[f*2+1] += PetscAbsScalar(a[d]); 928 } 929 } 930 ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 931 if (stepnum >= 0) { /* No summary for final time */ 932 DM dmCell, *fdm; 933 Vec *fv; 934 const PetscScalar *cgeom; 935 PetscScalar **fx; 936 PetscReal *fmin, *fmax, *fint, *ftmp, t; 937 PetscInt cStart, cEnd, c, fcount, f, num; 938 939 size_t ftableused,ftablealloc; 940 941 /* Functionals have indices after registering, this is an upper bound */ 942 fcount = user->numMonitorFuncs; 943 ierr = PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fint,fcount,&ftmp);CHKERRQ(ierr); 944 ierr = PetscMalloc3(fcount,&fdm,fcount,&fv,fcount,&fx);CHKERRQ(ierr); 945 for (f = 0; f < fcount; ++f) { 946 PetscSection fs; 947 const char *name = user->monitorFuncs[f]->name; 948 949 fmin[f] = PETSC_MAX_REAL; 950 fmax[f] = PETSC_MIN_REAL; 951 fint[f] = 0; 952 /* Make monitor vecs */ 953 ierr = DMClone(dm, &fdm[f]);CHKERRQ(ierr); 954 ierr = DMGetOutputSequenceNumber(dm, &num, &t);CHKERRQ(ierr); 955 ierr = DMSetOutputSequenceNumber(fdm[f], num, t);CHKERRQ(ierr); 956 ierr = PetscSectionClone(s, &fs);CHKERRQ(ierr); 957 ierr = PetscSectionSetFieldName(fs, 0, NULL);CHKERRQ(ierr); 958 ierr = PetscSectionSetFieldName(fs, 1, name);CHKERRQ(ierr); 959 ierr = DMSetLocalSection(fdm[f], fs);CHKERRQ(ierr); 960 ierr = PetscSectionDestroy(&fs);CHKERRQ(ierr); 961 ierr = DMGetGlobalVector(fdm[f], &fv[f]);CHKERRQ(ierr); 962 ierr = PetscObjectSetName((PetscObject) fv[f], name);CHKERRQ(ierr); 963 ierr = VecGetArray(fv[f], &fx[f]);CHKERRQ(ierr); 964 } 965 ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 966 ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr); 967 ierr = VecGetArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 968 ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 969 for (c = cStart; c < cEnd; ++c) { 970 PetscFVCellGeom *cg; 971 PetscScalar *cx; 972 973 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 974 ierr = DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx);CHKERRQ(ierr); 975 if (!cx) continue; /* not a global cell */ 976 for (f = 0; f < user->numMonitorFuncs; ++f) { 977 Functional func = user->monitorFuncs[f]; 978 PetscScalar *fxc; 979 980 ierr = DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc);CHKERRQ(ierr); 981 /* I need to make it easier to get interpolated values here */ 982 ierr = (*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx);CHKERRQ(ierr); 983 fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 984 } 985 for (f = 0; f < fcount; ++f) { 986 fmin[f] = PetscMin(fmin[f], ftmp[f]); 987 fmax[f] = PetscMax(fmax[f], ftmp[f]); 988 fint[f] += cg->volume * ftmp[f]; 989 } 990 } 991 ierr = VecRestoreArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 992 ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 993 ierr = MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 994 ierr = MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 995 ierr = MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts));CHKERRMPI(ierr); 996 /* Output functional data */ 997 ftablealloc = fcount * 100; 998 ftableused = 0; 999 ierr = PetscCalloc1(ftablealloc, &ftable);CHKERRQ(ierr); 1000 for (f = 0; f < user->numMonitorFuncs; ++f) { 1001 Functional func = user->monitorFuncs[f]; 1002 PetscInt id = func->offset; 1003 char newline[] = "\n"; 1004 char buffer[256], *p, *prefix; 1005 size_t countused, len; 1006 1007 /* Create string with functional outputs */ 1008 if (f % 3) { 1009 ierr = PetscArraycpy(buffer, " ", 2);CHKERRQ(ierr); 1010 p = buffer + 2; 1011 } else if (f) { 1012 ierr = PetscArraycpy(buffer, newline, sizeof(newline)-1);CHKERRQ(ierr); 1013 p = buffer + sizeof(newline) - 1; 1014 } else { 1015 p = buffer; 1016 } 1017 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); 1018 countused += p - buffer; 1019 /* reallocate */ 1020 if (countused > ftablealloc-ftableused-1) { 1021 char *ftablenew; 1022 1023 ftablealloc = 2*ftablealloc + countused; 1024 ierr = PetscMalloc1(ftablealloc, &ftablenew);CHKERRQ(ierr); 1025 ierr = PetscArraycpy(ftablenew, ftable, ftableused);CHKERRQ(ierr); 1026 ierr = PetscFree(ftable);CHKERRQ(ierr); 1027 ftable = ftablenew; 1028 } 1029 ierr = PetscArraycpy(ftable+ftableused, buffer, countused);CHKERRQ(ierr); 1030 ftableused += countused; 1031 ftable[ftableused] = 0; 1032 /* Output vecs */ 1033 ierr = VecRestoreArray(fv[f], &fx[f]);CHKERRQ(ierr); 1034 ierr = PetscStrlen(func->name, &len);CHKERRQ(ierr); 1035 ierr = PetscMalloc1(len+2,&prefix);CHKERRQ(ierr); 1036 ierr = PetscStrcpy(prefix, func->name);CHKERRQ(ierr); 1037 ierr = PetscStrcat(prefix, "_");CHKERRQ(ierr); 1038 ierr = PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix);CHKERRQ(ierr); 1039 ierr = VecViewFromOptions(fv[f], NULL, "-vec_view");CHKERRQ(ierr); 1040 ierr = PetscFree(prefix);CHKERRQ(ierr); 1041 ierr = DMRestoreGlobalVector(fdm[f], &fv[f]);CHKERRQ(ierr); 1042 ierr = DMDestroy(&fdm[f]);CHKERRQ(ierr); 1043 } 1044 ierr = PetscFree4(fmin, fmax, fint, ftmp);CHKERRQ(ierr); 1045 ierr = PetscFree3(fdm, fv, fx);CHKERRQ(ierr); 1046 ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), "% 3D time %8.4g |x| (", stepnum, (double) time);CHKERRQ(ierr); 1047 for (f = 0; f < Nf; ++f) { 1048 if (f > 0) {ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ", ");CHKERRQ(ierr);} 1049 ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+0]);CHKERRQ(ierr); 1050 } 1051 ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ") |x|_1 (");CHKERRQ(ierr); 1052 for (f = 0; f < Nf; ++f) { 1053 if (f > 0) {ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ", ");CHKERRQ(ierr);} 1054 ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+1]);CHKERRQ(ierr); 1055 } 1056 ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ") %s\n", ftable ? ftable : "");CHKERRQ(ierr); 1057 ierr = PetscFree(ftable);CHKERRQ(ierr); 1058 } 1059 ierr = PetscFree(xnorms);CHKERRQ(ierr); 1060 PetscFunctionReturn(0); 1061 } 1062 1063 int main(int argc, char **argv) 1064 { 1065 MPI_Comm comm; 1066 TS ts; 1067 DM dm; 1068 Vec u; 1069 AppCtx user; 1070 PetscReal t0, t = 0.0; 1071 void *ctxs[2]; 1072 PetscErrorCode ierr; 1073 1074 ctxs[0] = &t; 1075 ctxs[1] = &t; 1076 ierr = PetscInitialize(&argc, &argv, (char*) 0, help);if (ierr) return ierr; 1077 comm = PETSC_COMM_WORLD; 1078 user.functionalRegistry = NULL; 1079 globalUser = &user; 1080 ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 1081 ierr = TSCreate(comm, &ts);CHKERRQ(ierr); 1082 ierr = TSSetType(ts, TSBEULER);CHKERRQ(ierr); 1083 ierr = CreateDM(comm, &user, &dm);CHKERRQ(ierr); 1084 ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 1085 ierr = ProcessMonitorOptions(comm, &user);CHKERRQ(ierr); 1086 1087 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 1088 ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr); 1089 if (user.useFV) { 1090 PetscBool isImplicit = PETSC_FALSE; 1091 1092 ierr = PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit);CHKERRQ(ierr); 1093 if (isImplicit) { 1094 ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 1095 ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 1096 } 1097 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 1098 ierr = DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user);CHKERRQ(ierr); 1099 } else { 1100 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 1101 ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 1102 ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 1103 } 1104 if (user.useFV) {ierr = TSMonitorSet(ts, MonitorFunctionals, &user, NULL);CHKERRQ(ierr);} 1105 ierr = TSSetMaxSteps(ts, 1);CHKERRQ(ierr); 1106 ierr = TSSetMaxTime(ts, 2.0);CHKERRQ(ierr); 1107 ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr); 1108 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 1109 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1110 1111 ierr = DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u);CHKERRQ(ierr); 1112 if (user.useFV) {ierr = SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1]);CHKERRQ(ierr);} 1113 ierr = VecViewFromOptions(u, NULL, "-init_vec_view");CHKERRQ(ierr); 1114 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 1115 t0 = t; 1116 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 1117 ierr = TSSolve(ts, u);CHKERRQ(ierr); 1118 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 1119 if (t > t0) {ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);} 1120 ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 1121 { 1122 PetscReal ftime; 1123 PetscInt nsteps; 1124 TSConvergedReason reason; 1125 1126 ierr = TSGetSolveTime(ts, &ftime);CHKERRQ(ierr); 1127 ierr = TSGetStepNumber(ts, &nsteps);CHKERRQ(ierr); 1128 ierr = TSGetConvergedReason(ts, &reason);CHKERRQ(ierr); 1129 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %D steps\n", TSConvergedReasons[reason], (double) ftime, nsteps);CHKERRQ(ierr); 1130 } 1131 1132 ierr = VecDestroy(&u);CHKERRQ(ierr); 1133 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1134 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1135 ierr = PetscFree(user.monitorFuncs);CHKERRQ(ierr); 1136 ierr = FunctionalDestroy(&user.functionalRegistry);CHKERRQ(ierr); 1137 ierr = PetscFinalize(); 1138 return ierr; 1139 } 1140 1141 /*TEST 1142 1143 testset: 1144 args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 1145 1146 # 2D harmonic velocity, no porosity 1147 test: 1148 suffix: p1p1 1149 requires: !complex !single 1150 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 1151 test: 1152 suffix: p1p1_xper 1153 requires: !complex !single 1154 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 1155 test: 1156 suffix: p1p1_xper_ref 1157 requires: !complex !single 1158 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 1159 test: 1160 suffix: p1p1_xyper 1161 requires: !complex !single 1162 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 1163 test: 1164 suffix: p1p1_xyper_ref 1165 requires: !complex !single 1166 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 1167 test: 1168 suffix: p2p1 1169 requires: !complex !single 1170 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 1171 test: 1172 suffix: p2p1_xyper 1173 requires: !complex !single 1174 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 1175 1176 test: 1177 suffix: adv_1 1178 requires: !complex !single 1179 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 1180 1181 test: 1182 suffix: adv_2 1183 requires: !complex 1184 TODO: broken memory corruption 1185 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 1186 1187 test: 1188 suffix: adv_3 1189 requires: !complex 1190 TODO: broken memory corruption 1191 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 1192 1193 test: 1194 suffix: adv_3_ex 1195 requires: !complex 1196 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 1197 1198 test: 1199 suffix: adv_4 1200 requires: !complex 1201 TODO: broken memory corruption 1202 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 1203 1204 # 2D Advection, box, delta 1205 test: 1206 suffix: adv_delta_yper_0 1207 requires: !complex 1208 TODO: broken 1209 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 1210 1211 test: 1212 suffix: adv_delta_yper_1 1213 requires: !complex 1214 TODO: broken 1215 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 1216 1217 test: 1218 suffix: adv_delta_yper_2 1219 requires: !complex 1220 TODO: broken 1221 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 1222 1223 test: 1224 suffix: adv_delta_yper_fim_0 1225 requires: !complex 1226 TODO: broken 1227 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 1228 1229 test: 1230 suffix: adv_delta_yper_fim_1 1231 requires: !complex 1232 TODO: broken 1233 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 1234 1235 test: 1236 suffix: adv_delta_yper_fim_2 1237 requires: !complex 1238 TODO: broken 1239 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 1240 1241 test: 1242 suffix: adv_delta_yper_im_0 1243 requires: !complex 1244 TODO: broken 1245 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 1246 1247 test: 1248 suffix: adv_delta_yper_im_1 1249 requires: !complex 1250 TODO: broken 1251 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 1252 1253 test: 1254 suffix: adv_delta_yper_im_2 1255 requires: !complex 1256 TODO: broken 1257 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 1258 1259 test: 1260 suffix: adv_delta_yper_im_3 1261 requires: !complex 1262 TODO: broken 1263 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 1264 1265 # I believe the nullspace is sin(pi y) 1266 test: 1267 suffix: adv_delta_yper_im_4 1268 requires: !complex 1269 TODO: broken 1270 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 1271 1272 test: 1273 suffix: adv_delta_yper_im_5 1274 requires: !complex 1275 TODO: broken 1276 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 1277 1278 test: 1279 suffix: adv_delta_yper_im_6 1280 requires: !complex 1281 TODO: broken 1282 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 1283 # 2D Advection, magma benchmark 1 1284 1285 test: 1286 suffix: adv_delta_shear_im_0 1287 requires: !complex 1288 TODO: broken 1289 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 1290 # 2D Advection, box, gaussian 1291 1292 test: 1293 suffix: adv_gauss 1294 requires: !complex 1295 TODO: broken 1296 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 1297 1298 test: 1299 suffix: adv_gauss_im 1300 requires: !complex 1301 TODO: broken 1302 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 1303 1304 test: 1305 suffix: adv_gauss_im_1 1306 requires: !complex 1307 TODO: broken 1308 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 1309 1310 test: 1311 suffix: adv_gauss_im_2 1312 requires: !complex 1313 TODO: broken 1314 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 1315 1316 test: 1317 suffix: adv_gauss_corner 1318 requires: !complex 1319 TODO: broken 1320 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 1321 1322 # 2D Advection+Harmonic 12- 1323 test: 1324 suffix: adv_harm_0 1325 requires: !complex 1326 TODO: broken memory corruption 1327 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 1328 1329 # Must check that FV BCs propagate to coarse meshes 1330 # Must check that FV BC ids propagate to coarse meshes 1331 # Must check that FE+FV BCs work at the same time 1332 # 2D Advection, matching wind in ex11 8-11 1333 # NOTE implicit solves are limited by accuracy of FD Jacobian 1334 test: 1335 suffix: adv_0 1336 requires: !complex !single exodusii 1337 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 1338 1339 test: 1340 suffix: adv_0_im 1341 requires: !complex exodusii 1342 TODO: broken memory corruption 1343 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 1344 1345 test: 1346 suffix: adv_0_im_2 1347 requires: !complex exodusii 1348 TODO: broken 1349 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 1350 1351 test: 1352 suffix: adv_0_im_3 1353 requires: !complex exodusii 1354 TODO: broken 1355 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 1356 1357 test: 1358 suffix: adv_0_im_4 1359 requires: !complex exodusii 1360 TODO: broken 1361 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 1362 # 2D Advection, misc 1363 1364 TEST*/ 1365