1c4762a1bSJed Brown static char help[] = "Hybrid Finite Element-Finite Volume Example.\n"; 2c4762a1bSJed Brown /*F 3c4762a1bSJed Brown Here we are advecting a passive tracer in a harmonic velocity field, defined by 4c4762a1bSJed Brown a forcing function $f$: 5c4762a1bSJed Brown \begin{align} 6c4762a1bSJed Brown -\Delta \mathbf{u} + f &= 0 \\ 7c4762a1bSJed Brown \frac{\partial\phi}{\partial t} + \nabla\cdot \phi \mathbf{u} &= 0 8c4762a1bSJed Brown \end{align} 9c4762a1bSJed Brown F*/ 10c4762a1bSJed Brown 11c4762a1bSJed Brown #include <petscdmplex.h> 12c4762a1bSJed Brown #include <petscds.h> 13c4762a1bSJed Brown #include <petscts.h> 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h> /* For DotD */ 16c4762a1bSJed Brown 179371c9d4SSatish Balay typedef enum { 189371c9d4SSatish Balay VEL_ZERO, 199371c9d4SSatish Balay VEL_CONSTANT, 209371c9d4SSatish Balay VEL_HARMONIC, 219371c9d4SSatish Balay VEL_SHEAR 229371c9d4SSatish Balay } VelocityDistribution; 23c4762a1bSJed Brown 249371c9d4SSatish Balay typedef enum { 259371c9d4SSatish Balay ZERO, 269371c9d4SSatish Balay CONSTANT, 279371c9d4SSatish Balay GAUSSIAN, 289371c9d4SSatish Balay TILTED, 299371c9d4SSatish Balay DELTA 309371c9d4SSatish Balay } PorosityDistribution; 31c4762a1bSJed Brown 32c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown FunctionalFunc - Calculates the value of a functional of the solution at a point 36c4762a1bSJed Brown 37c4762a1bSJed Brown Input Parameters: 38c4762a1bSJed Brown + dm - The DM 39c4762a1bSJed Brown . time - The TS time 40c4762a1bSJed Brown . x - The coordinates of the evaluation point 41c4762a1bSJed Brown . u - The field values at point x 42c4762a1bSJed Brown - ctx - A user context, or NULL 43c4762a1bSJed Brown 44c4762a1bSJed Brown Output Parameter: 45c4762a1bSJed Brown . f - The value of the functional at point x 46c4762a1bSJed Brown 47c4762a1bSJed Brown */ 48c4762a1bSJed Brown typedef PetscErrorCode (*FunctionalFunc)(DM, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); 49c4762a1bSJed Brown 50c4762a1bSJed Brown typedef struct _n_Functional *Functional; 51c4762a1bSJed Brown struct _n_Functional { 52c4762a1bSJed Brown char *name; 53c4762a1bSJed Brown FunctionalFunc func; 54c4762a1bSJed Brown void *ctx; 55c4762a1bSJed Brown PetscInt offset; 56c4762a1bSJed Brown Functional next; 57c4762a1bSJed Brown }; 58c4762a1bSJed Brown 59c4762a1bSJed Brown typedef struct { 60c4762a1bSJed Brown /* Problem definition */ 61c4762a1bSJed Brown PetscBool useFV; /* Use a finite volume scheme for advection */ 62c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 63c4762a1bSJed Brown VelocityDistribution velocityDist; 64c4762a1bSJed Brown PorosityDistribution porosityDist; 65c4762a1bSJed Brown PetscReal inflowState; 66c4762a1bSJed Brown PetscReal source[3]; 67c4762a1bSJed Brown /* Monitoring */ 68c4762a1bSJed Brown PetscInt numMonitorFuncs, maxMonitorFunc; 69c4762a1bSJed Brown Functional *monitorFuncs; 70c4762a1bSJed Brown PetscInt errorFunctional; 71c4762a1bSJed Brown Functional functionalRegistry; 72c4762a1bSJed Brown } AppCtx; 73c4762a1bSJed Brown 74c4762a1bSJed Brown static AppCtx *globalUser; 75c4762a1bSJed Brown 76d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 77d71ae5a4SJacob Faibussowitsch { 78c4762a1bSJed Brown const char *velocityDist[4] = {"zero", "constant", "harmonic", "shear"}; 79c4762a1bSJed Brown const char *porosityDist[5] = {"zero", "constant", "gaussian", "tilted", "delta"}; 8030602db0SMatthew G. Knepley PetscInt vd, pd, d; 81c4762a1bSJed Brown PetscBool flg; 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscFunctionBeginUser; 84c4762a1bSJed Brown options->useFV = PETSC_FALSE; 85c4762a1bSJed Brown options->velocityDist = VEL_HARMONIC; 86c4762a1bSJed Brown options->porosityDist = ZERO; 87c4762a1bSJed Brown options->inflowState = -2.0; 88c4762a1bSJed Brown options->numMonitorFuncs = 0; 89c4762a1bSJed Brown options->source[0] = 0.5; 90c4762a1bSJed Brown options->source[1] = 0.5; 91c4762a1bSJed Brown options->source[2] = 0.5; 92c4762a1bSJed Brown 93d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX"); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL)); 95c4762a1bSJed Brown vd = options->velocityDist; 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-velocity_dist", "Velocity distribution type", "ex18.c", velocityDist, 4, velocityDist[options->velocityDist], &vd, NULL)); 97c4762a1bSJed Brown options->velocityDist = (VelocityDistribution)vd; 98c4762a1bSJed Brown pd = options->porosityDist; 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-porosity_dist", "Initial porosity distribution type", "ex18.c", porosityDist, 5, porosityDist[options->porosityDist], &pd, NULL)); 100c4762a1bSJed Brown options->porosityDist = (PorosityDistribution)pd; 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL)); 10230602db0SMatthew G. Knepley d = 2; 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg)); 10463a3b9bcSJacob Faibussowitsch PetscCheck(!flg || d == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %" PetscInt_FMT, d); 105d0609cedSBarry Smith PetscOptionsEnd(); 106c4762a1bSJed Brown 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) 111d71ae5a4SJacob Faibussowitsch { 112c4762a1bSJed Brown Functional func; 113c4762a1bSJed Brown char *names[256]; 114c4762a1bSJed Brown PetscInt f; 115c4762a1bSJed Brown 116c4762a1bSJed Brown PetscFunctionBeginUser; 117d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX"); 118dd39110bSPierre Jolivet options->numMonitorFuncs = PETSC_STATIC_ARRAY_LENGTH(names); 1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL)); 1209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs)); 121c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 122c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 123c4762a1bSJed Brown PetscBool match; 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(names[f], func->name, &match)); 126c4762a1bSJed Brown if (match) break; 127c4762a1bSJed Brown } 1283c633725SBarry Smith PetscCheck(func, comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 129c4762a1bSJed Brown options->monitorFuncs[f] = func; 130c4762a1bSJed Brown /* Jed inserts a de-duplication of functionals here */ 1319566063dSJacob Faibussowitsch PetscCall(PetscFree(names[f])); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 134c4762a1bSJed Brown options->maxMonitorFunc = -1; 135c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 136c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 137c4762a1bSJed Brown Functional call = options->monitorFuncs[f]; 138c4762a1bSJed Brown 139c4762a1bSJed Brown if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } 142d0609cedSBarry Smith PetscOptionsEnd(); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146d71ae5a4SJacob Faibussowitsch static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) 147d71ae5a4SJacob Faibussowitsch { 148c4762a1bSJed Brown Functional *ptr, f; 149c4762a1bSJed Brown PetscInt lastoffset = -1; 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFunctionBeginUser; 152c4762a1bSJed Brown for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 1539566063dSJacob Faibussowitsch PetscCall(PetscNew(&f)); 1549566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &f->name)); 155c4762a1bSJed Brown f->offset = lastoffset + 1; 156c4762a1bSJed Brown f->func = func; 157c4762a1bSJed Brown f->ctx = ctx; 158c4762a1bSJed Brown f->next = NULL; 159c4762a1bSJed Brown *ptr = f; 160c4762a1bSJed Brown *offset = f->offset; 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164d71ae5a4SJacob Faibussowitsch static PetscErrorCode FunctionalDestroy(Functional *link) 165d71ae5a4SJacob Faibussowitsch { 166c4762a1bSJed Brown Functional next, l; 167c4762a1bSJed Brown 168c4762a1bSJed Brown PetscFunctionBeginUser; 1693ba16761SJacob Faibussowitsch if (!link) PetscFunctionReturn(PETSC_SUCCESS); 170c4762a1bSJed Brown l = *link; 171c4762a1bSJed Brown *link = NULL; 172c4762a1bSJed Brown for (; l; l = next) { 173c4762a1bSJed Brown next = l->next; 1749566063dSJacob Faibussowitsch PetscCall(PetscFree(l->name)); 1759566063dSJacob Faibussowitsch PetscCall(PetscFree(l)); 176c4762a1bSJed Brown } 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180d71ae5a4SJacob Faibussowitsch static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 181d71ae5a4SJacob Faibussowitsch { 182c4762a1bSJed Brown PetscInt comp; 183c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186d71ae5a4SJacob Faibussowitsch static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 187d71ae5a4SJacob Faibussowitsch { 188c4762a1bSJed Brown PetscScalar wind[3] = {0.0, 0.0, 0.0}; 189c4762a1bSJed Brown PetscInt comp; 190c4762a1bSJed Brown 1913ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, constant_u_2d(dim, t, x, Nf, wind, NULL)); 192c4762a1bSJed Brown for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp]; 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195d71ae5a4SJacob Faibussowitsch static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 196d71ae5a4SJacob Faibussowitsch { 197c4762a1bSJed Brown PetscInt comp; 198c4762a1bSJed Brown for (comp = 0; comp < dim * dim; ++comp) f1[comp] = 0.0; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201d71ae5a4SJacob Faibussowitsch static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 202d71ae5a4SJacob Faibussowitsch { 203c4762a1bSJed Brown PetscInt d; 204c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d * dim + d] = 1.0; 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207d71ae5a4SJacob Faibussowitsch static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 208d71ae5a4SJacob Faibussowitsch { 209c4762a1bSJed Brown g0[0] = 1.0; 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212d71ae5a4SJacob Faibussowitsch static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 213d71ae5a4SJacob Faibussowitsch { 214c4762a1bSJed Brown PetscInt comp; 215c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218d71ae5a4SJacob Faibussowitsch static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 219d71ae5a4SJacob Faibussowitsch { 220c4762a1bSJed Brown PetscInt comp, d; 221c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) { 222ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f1[comp * dim + d] = u_x[comp * dim + d]; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226d71ae5a4SJacob Faibussowitsch static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 227d71ae5a4SJacob Faibussowitsch { 228c4762a1bSJed Brown f0[0] = -PetscSinReal(2.0 * PETSC_PI * x[0]); 229c4762a1bSJed Brown f0[1] = 2.0 * PETSC_PI * x[1] * PetscCosReal(2.0 * PETSC_PI * x[0]); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown 232d71ae5a4SJacob Faibussowitsch static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 233d71ae5a4SJacob Faibussowitsch { 234c4762a1bSJed Brown f0[0] = -2.0 * PetscSinReal(2.0 * PETSC_PI * x[0]) * PetscCosReal(2.0 * PETSC_PI * x[1]); 235c4762a1bSJed Brown f0[1] = 2.0 * PetscSinReal(2.0 * PETSC_PI * x[1]) * PetscCosReal(2.0 * PETSC_PI * x[0]); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238d71ae5a4SJacob Faibussowitsch void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 239d71ae5a4SJacob Faibussowitsch { 240c4762a1bSJed Brown const PetscInt Ncomp = dim; 241c4762a1bSJed Brown PetscInt compI, d; 242c4762a1bSJed Brown 243c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 244ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g3[((compI * Ncomp + compI) * dim + d) * dim + d] = 1.0; 245c4762a1bSJed Brown } 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 249d71ae5a4SJacob Faibussowitsch static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 250d71ae5a4SJacob Faibussowitsch { 251c4762a1bSJed Brown PetscInt d; 252c4762a1bSJed Brown f0[0] = u_t[dim]; 253c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += u[dim] * u_x[d * dim + d] + u_x[dim * dim + d] * u[d]; 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256d71ae5a4SJacob Faibussowitsch static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 257d71ae5a4SJacob Faibussowitsch { 258c4762a1bSJed Brown PetscInt d; 259c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[0] = 0.0; 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262d71ae5a4SJacob Faibussowitsch void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 263d71ae5a4SJacob Faibussowitsch { 264c4762a1bSJed Brown PetscInt d; 265c4762a1bSJed Brown g0[0] = u_tShift; 266c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[d * dim + d]; 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 269d71ae5a4SJacob Faibussowitsch void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 270d71ae5a4SJacob Faibussowitsch { 271c4762a1bSJed Brown PetscInt d; 272c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = u[d]; 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275d71ae5a4SJacob Faibussowitsch void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 276d71ae5a4SJacob Faibussowitsch { 277c4762a1bSJed Brown PetscInt d; 278c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[dim * dim + d]; 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281d71ae5a4SJacob Faibussowitsch void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 282d71ae5a4SJacob Faibussowitsch { 283c4762a1bSJed Brown PetscInt d; 284c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = u[dim]; 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287d71ae5a4SJacob Faibussowitsch 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) 288d71ae5a4SJacob Faibussowitsch { 289c4762a1bSJed Brown PetscReal wind[3] = {0.0, 1.0, 0.0}; 290c4762a1bSJed Brown PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim, 3), wind, n); 291c4762a1bSJed Brown 292c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295d71ae5a4SJacob Faibussowitsch 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) 296d71ae5a4SJacob Faibussowitsch { 297c4762a1bSJed Brown PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 298c4762a1bSJed Brown 299c4762a1bSJed Brown #if 1 300c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 301c4762a1bSJed Brown #else 302c4762a1bSJed Brown /* 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]); */ 303c4762a1bSJed Brown /* Smear it out */ 304c4762a1bSJed Brown flux[0] = 0.5 * ((uL[dim] + uR[dim]) + (uL[dim] - uR[dim]) * tanh(1.0e5 * wn)) * wn; 305c4762a1bSJed Brown #endif 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 309d71ae5a4SJacob Faibussowitsch { 310c4762a1bSJed Brown u[0] = 0.0; 311c4762a1bSJed Brown u[1] = 0.0; 3123ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 313c4762a1bSJed Brown } 314c4762a1bSJed Brown 315d71ae5a4SJacob Faibussowitsch static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 316d71ae5a4SJacob Faibussowitsch { 317c4762a1bSJed Brown u[0] = 0.0; 318c4762a1bSJed Brown u[1] = 1.0; 3193ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 320c4762a1bSJed Brown } 321c4762a1bSJed Brown 322c4762a1bSJed Brown /* Coordinates of the point which was at x at t = 0 */ 323d71ae5a4SJacob Faibussowitsch static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 324d71ae5a4SJacob Faibussowitsch { 325c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 326c4762a1bSJed Brown u[0] = x[0]; 327c4762a1bSJed Brown u[1] = x[1] + t; 32830602db0SMatthew G. Knepley #if 0 3299566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u)); 33030602db0SMatthew G. Knepley #else 33130602db0SMatthew G. Knepley u[1] = u[1] - (int)PetscRealPart(u[1]); 33230602db0SMatthew G. Knepley #endif 3333ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* 337c4762a1bSJed Brown In 2D we use the exact solution: 338c4762a1bSJed Brown 339c4762a1bSJed Brown u = x^2 + y^2 340c4762a1bSJed Brown v = 2 x^2 - 2xy 341c4762a1bSJed Brown phi = h(x + y + (u + v) t) 342c4762a1bSJed Brown f_x = f_y = 4 343c4762a1bSJed Brown 344c4762a1bSJed Brown so that 345c4762a1bSJed Brown 346c4762a1bSJed Brown -\Delta u + f = <-4, -4> + <4, 4> = 0 347c4762a1bSJed Brown {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 348c4762a1bSJed Brown h_t(x + y + (u + v) t) - u . grad phi - phi div u 349c4762a1bSJed Brown = u h' + v h' - u h_x - v h_y 350c4762a1bSJed Brown = 0 351c4762a1bSJed Brown 352c4762a1bSJed Brown We will conserve phi since 353c4762a1bSJed Brown 354c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 355c4762a1bSJed Brown 356c4762a1bSJed Brown Also try h((x + ut)^2 + (y + vt)^2), so that 357c4762a1bSJed Brown 358c4762a1bSJed Brown h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 359c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 360c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 361c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 362c4762a1bSJed Brown = 0 363c4762a1bSJed Brown 364c4762a1bSJed Brown */ 365d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 366d71ae5a4SJacob Faibussowitsch { 367c4762a1bSJed Brown u[0] = x[0] * x[0] + x[1] * x[1]; 368c4762a1bSJed Brown u[1] = 2.0 * x[0] * x[0] - 2.0 * x[0] * x[1]; 3693ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 370c4762a1bSJed Brown } 371c4762a1bSJed Brown 372c4762a1bSJed Brown /* 373c4762a1bSJed Brown In 2D we use the exact, periodic solution: 374c4762a1bSJed Brown 375c4762a1bSJed Brown u = sin(2 pi x)/4 pi^2 376c4762a1bSJed Brown v = -y cos(2 pi x)/2 pi 377c4762a1bSJed Brown phi = h(x + y + (u + v) t) 378c4762a1bSJed Brown f_x = -sin(2 pi x) 379c4762a1bSJed Brown f_y = 2 pi y cos(2 pi x) 380c4762a1bSJed Brown 381c4762a1bSJed Brown so that 382c4762a1bSJed Brown 383c4762a1bSJed Brown -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 384c4762a1bSJed Brown 385c4762a1bSJed Brown We will conserve phi since 386c4762a1bSJed Brown 387c4762a1bSJed Brown \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 388c4762a1bSJed Brown */ 389d71ae5a4SJacob Faibussowitsch static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 390d71ae5a4SJacob Faibussowitsch { 391c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]) / PetscSqr(2.0 * PETSC_PI); 392c4762a1bSJed Brown u[1] = -x[1] * PetscCosReal(2.0 * PETSC_PI * x[0]) / (2.0 * PETSC_PI); 3933ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 396c4762a1bSJed Brown /* 397c4762a1bSJed Brown In 2D we use the exact, doubly periodic solution: 398c4762a1bSJed Brown 399c4762a1bSJed Brown u = sin(2 pi x) cos(2 pi y)/4 pi^2 400c4762a1bSJed Brown v = -sin(2 pi y) cos(2 pi x)/4 pi^2 401c4762a1bSJed Brown phi = h(x + y + (u + v) t) 402c4762a1bSJed Brown f_x = -2sin(2 pi x) cos(2 pi y) 403c4762a1bSJed Brown f_y = 2sin(2 pi y) cos(2 pi x) 404c4762a1bSJed Brown 405c4762a1bSJed Brown so that 406c4762a1bSJed Brown 407c4762a1bSJed Brown -\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 408c4762a1bSJed Brown 409c4762a1bSJed Brown We will conserve phi since 410c4762a1bSJed Brown 411c4762a1bSJed Brown \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 412c4762a1bSJed Brown */ 413d71ae5a4SJacob Faibussowitsch static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 414d71ae5a4SJacob Faibussowitsch { 415c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]) * PetscCosReal(2.0 * PETSC_PI * x[1]) / PetscSqr(2.0 * PETSC_PI); 416c4762a1bSJed Brown u[1] = -PetscSinReal(2.0 * PETSC_PI * x[1]) * PetscCosReal(2.0 * PETSC_PI * x[0]) / PetscSqr(2.0 * PETSC_PI); 4173ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 418c4762a1bSJed Brown } 419c4762a1bSJed Brown 420d71ae5a4SJacob Faibussowitsch static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 421d71ae5a4SJacob Faibussowitsch { 422c4762a1bSJed Brown u[0] = x[1] - 0.5; 423c4762a1bSJed Brown u[1] = 0.0; 4243ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427d71ae5a4SJacob Faibussowitsch static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 428d71ae5a4SJacob Faibussowitsch { 429c4762a1bSJed Brown PetscInt d; 430c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 4313ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 432c4762a1bSJed Brown } 433c4762a1bSJed Brown 434d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 435d71ae5a4SJacob Faibussowitsch { 436c4762a1bSJed Brown u[0] = 0.0; 4373ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 438c4762a1bSJed Brown } 439c4762a1bSJed Brown 440d71ae5a4SJacob Faibussowitsch static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 441d71ae5a4SJacob Faibussowitsch { 442c4762a1bSJed Brown u[0] = 1.0; 4433ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 444c4762a1bSJed Brown } 445c4762a1bSJed Brown 446d71ae5a4SJacob Faibussowitsch static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 447d71ae5a4SJacob Faibussowitsch { 448c4762a1bSJed Brown PetscReal x0[2]; 449c4762a1bSJed Brown PetscScalar xn[2]; 450c4762a1bSJed Brown 451c4762a1bSJed Brown x0[0] = globalUser->source[0]; 452c4762a1bSJed Brown x0[1] = globalUser->source[1]; 4533ba16761SJacob Faibussowitsch PetscCall(constant_x_2d(dim, time, x0, Nf, xn, ctx)); 454c4762a1bSJed Brown { 455c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 456c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 457c4762a1bSJed Brown const PetscReal r2 = xi * xi + eta * eta; 458c4762a1bSJed Brown 459c4762a1bSJed Brown u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 460c4762a1bSJed Brown } 4613ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 462c4762a1bSJed Brown } 463c4762a1bSJed Brown 464c4762a1bSJed Brown /* 465c4762a1bSJed Brown Gaussian blob, initially centered on (0.5, 0.5) 466c4762a1bSJed Brown 467c4762a1bSJed Brown xi = x(t) - x0, eta = y(t) - y0 468c4762a1bSJed Brown 469c4762a1bSJed Brown where x(t), y(t) are the integral curves of v(t), 470c4762a1bSJed Brown 471c4762a1bSJed Brown dx/dt . grad f = v . f 472c4762a1bSJed Brown 473c4762a1bSJed Brown Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 474c4762a1bSJed Brown 475c4762a1bSJed Brown v0 f_x + w0 f_y = v . f 476c4762a1bSJed Brown */ 477d71ae5a4SJacob Faibussowitsch static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 478d71ae5a4SJacob Faibussowitsch { 479c4762a1bSJed Brown const PetscReal x0[2] = {0.5, 0.5}; 480c4762a1bSJed Brown const PetscReal sigma = 1.0 / 6.0; 481c4762a1bSJed Brown PetscScalar xn[2]; 482c4762a1bSJed Brown 4833ba16761SJacob Faibussowitsch PetscCall(constant_x_2d(dim, time, x0, Nf, xn, ctx)); 484c4762a1bSJed Brown { 485c4762a1bSJed Brown /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 486c4762a1bSJed Brown /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 487c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 488c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 489c4762a1bSJed Brown const PetscReal r2 = xi * xi + eta * eta; 490c4762a1bSJed Brown 491c4762a1bSJed Brown u[0] = PetscExpReal(-r2 / (2.0 * sigma * sigma)) / (sigma * PetscSqrtReal(2.0 * PETSC_PI)); 492c4762a1bSJed Brown } 4933ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 494c4762a1bSJed Brown } 495c4762a1bSJed Brown 496d71ae5a4SJacob Faibussowitsch static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 497d71ae5a4SJacob Faibussowitsch { 498c4762a1bSJed Brown PetscReal x0[3]; 499c4762a1bSJed Brown const PetscReal wind[3] = {0.0, 1.0, 0.0}; 500c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 501c4762a1bSJed Brown 502c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 503c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0 * x[0] + 3.0 * x[1]; 504c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 5053ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 506c4762a1bSJed Brown } 507c4762a1bSJed Brown 508d71ae5a4SJacob Faibussowitsch static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 509d71ae5a4SJacob Faibussowitsch { 510c4762a1bSJed Brown PetscReal ur[3]; 511c4762a1bSJed Brown PetscReal x0[3]; 512c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 513c4762a1bSJed Brown 5149371c9d4SSatish Balay ur[0] = PetscRealPart(u[0]); 5159371c9d4SSatish Balay ur[1] = PetscRealPart(u[1]); 5169371c9d4SSatish Balay ur[2] = PetscRealPart(u[2]); 517c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 518c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0 * x[0] + 3.0 * x[1]; 519c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 5203ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 521c4762a1bSJed Brown } 522c4762a1bSJed Brown 523d71ae5a4SJacob Faibussowitsch static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 524d71ae5a4SJacob Faibussowitsch { 525c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 526c4762a1bSJed Brown 527c4762a1bSJed Brown PetscFunctionBeginUser; 528c4762a1bSJed Brown xG[0] = user->inflowState; 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 530c4762a1bSJed Brown } 531c4762a1bSJed Brown 532d71ae5a4SJacob Faibussowitsch static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 533d71ae5a4SJacob Faibussowitsch { 534c4762a1bSJed Brown PetscFunctionBeginUser; 53530602db0SMatthew G. Knepley //xG[0] = xI[dim]; 53630602db0SMatthew G. Knepley xG[0] = xI[2]; 5373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 538c4762a1bSJed Brown } 539c4762a1bSJed Brown 540d71ae5a4SJacob Faibussowitsch static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 541d71ae5a4SJacob Faibussowitsch { 542c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 543c4762a1bSJed Brown PetscInt dim; 544c4762a1bSJed Brown 545c4762a1bSJed Brown PetscFunctionBeginUser; 5469566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 547c4762a1bSJed Brown switch (user->porosityDist) { 548c4762a1bSJed Brown case TILTED: 5493ba16761SJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) PetscCall(tilted_phi_2d(dim, time, x, 2, u, (void *)&time)); 5503ba16761SJacob Faibussowitsch else PetscCall(tilted_phi_coupled_2d(dim, time, x, 2, u, (void *)&time)); 551c4762a1bSJed Brown break; 552d71ae5a4SJacob Faibussowitsch case GAUSSIAN: 5533ba16761SJacob Faibussowitsch PetscCall(gaussian_phi_2d(dim, time, x, 2, u, (void *)&time)); 554d71ae5a4SJacob Faibussowitsch break; 555d71ae5a4SJacob Faibussowitsch case DELTA: 5563ba16761SJacob Faibussowitsch PetscCall(delta_phi_2d(dim, time, x, 2, u, (void *)&time)); 557d71ae5a4SJacob Faibussowitsch break; 558d71ae5a4SJacob Faibussowitsch default: 559d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 560c4762a1bSJed Brown } 5613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 562c4762a1bSJed Brown } 563c4762a1bSJed Brown 564d71ae5a4SJacob Faibussowitsch static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 565d71ae5a4SJacob Faibussowitsch { 566c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 567c4762a1bSJed Brown PetscScalar yexact[3] = {0, 0, 0}; 568c4762a1bSJed Brown 569c4762a1bSJed Brown PetscFunctionBeginUser; 5709566063dSJacob Faibussowitsch PetscCall(ExactSolution(dm, time, x, yexact, ctx)); 571c4762a1bSJed Brown f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 573c4762a1bSJed Brown } 574c4762a1bSJed Brown 575d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 576d71ae5a4SJacob Faibussowitsch { 577c4762a1bSJed Brown PetscFunctionBeginUser; 5789566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 5799566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 5809566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 5819566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 5823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 583c4762a1bSJed Brown } 584c4762a1bSJed Brown 585d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupBC(DM dm, AppCtx *user) 586d71ae5a4SJacob Faibussowitsch { 587348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 58830602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 589c4762a1bSJed Brown PetscDS prob; 590c4762a1bSJed Brown DMLabel label; 591c4762a1bSJed Brown PetscBool check; 59230602db0SMatthew G. Knepley PetscInt dim, n = 3; 59330602db0SMatthew G. Knepley const char *prefix; 594c4762a1bSJed Brown 595c4762a1bSJed Brown PetscFunctionBeginUser; 5969566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 5979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *)bdt, &n, NULL)); 5989566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 599c4762a1bSJed Brown /* Set initial guesses and exact solutions */ 60030602db0SMatthew G. Knepley switch (dim) { 601c4762a1bSJed Brown case 2: 602c4762a1bSJed Brown user->initialGuess[0] = initialVelocity; 603c4762a1bSJed Brown switch (user->porosityDist) { 604d71ae5a4SJacob Faibussowitsch case ZERO: 605d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = zero_phi; 606d71ae5a4SJacob Faibussowitsch break; 607d71ae5a4SJacob Faibussowitsch case CONSTANT: 608d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = constant_phi; 609d71ae5a4SJacob Faibussowitsch break; 610d71ae5a4SJacob Faibussowitsch case GAUSSIAN: 611d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = gaussian_phi_2d; 612d71ae5a4SJacob Faibussowitsch break; 613d71ae5a4SJacob Faibussowitsch case DELTA: 614d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = delta_phi_2d; 615d71ae5a4SJacob Faibussowitsch break; 616c4762a1bSJed Brown case TILTED: 617c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 618c4762a1bSJed Brown else user->initialGuess[1] = tilted_phi_coupled_2d; 619c4762a1bSJed Brown break; 620c4762a1bSJed Brown } 621348a1646SMatthew G. Knepley break; 622d71ae5a4SJacob Faibussowitsch default: 623d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 624348a1646SMatthew G. Knepley } 625348a1646SMatthew G. Knepley exactFuncs[0] = user->initialGuess[0]; 626348a1646SMatthew G. Knepley exactFuncs[1] = user->initialGuess[1]; 62730602db0SMatthew G. Knepley switch (dim) { 628348a1646SMatthew G. Knepley case 2: 629c4762a1bSJed Brown switch (user->velocityDist) { 630d71ae5a4SJacob Faibussowitsch case VEL_ZERO: 631d71ae5a4SJacob Faibussowitsch exactFuncs[0] = zero_u_2d; 632d71ae5a4SJacob Faibussowitsch break; 633d71ae5a4SJacob Faibussowitsch case VEL_CONSTANT: 634d71ae5a4SJacob Faibussowitsch exactFuncs[0] = constant_u_2d; 635d71ae5a4SJacob Faibussowitsch break; 636c4762a1bSJed Brown case VEL_HARMONIC: 63730602db0SMatthew G. Knepley switch (bdt[0]) { 638c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 63930602db0SMatthew G. Knepley switch (bdt[1]) { 640d71ae5a4SJacob Faibussowitsch case DM_BOUNDARY_PERIODIC: 641d71ae5a4SJacob Faibussowitsch exactFuncs[0] = doubly_periodic_u_2d; 642d71ae5a4SJacob Faibussowitsch break; 643d71ae5a4SJacob Faibussowitsch default: 644d71ae5a4SJacob Faibussowitsch exactFuncs[0] = periodic_u_2d; 645d71ae5a4SJacob Faibussowitsch break; 646c4762a1bSJed Brown } 647c4762a1bSJed Brown break; 648d71ae5a4SJacob Faibussowitsch default: 649d71ae5a4SJacob Faibussowitsch exactFuncs[0] = quadratic_u_2d; 650d71ae5a4SJacob Faibussowitsch break; 651c4762a1bSJed Brown } 652c4762a1bSJed Brown break; 653d71ae5a4SJacob Faibussowitsch case VEL_SHEAR: 654d71ae5a4SJacob Faibussowitsch exactFuncs[0] = shear_bc; 655d71ae5a4SJacob Faibussowitsch break; 656d71ae5a4SJacob Faibussowitsch default: 657d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 658c4762a1bSJed Brown } 659348a1646SMatthew G. Knepley break; 660d71ae5a4SJacob Faibussowitsch default: 661d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 662c4762a1bSJed Brown } 663c4762a1bSJed Brown { 664c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 665c4762a1bSJed Brown 6669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, "", "-use_implicit", &isImplicit)); 667348a1646SMatthew G. Knepley if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 668c4762a1bSJed Brown } 6699566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-dmts_check", &check)); 670c4762a1bSJed Brown if (check) { 671348a1646SMatthew G. Knepley user->initialGuess[0] = exactFuncs[0]; 672348a1646SMatthew G. Knepley user->initialGuess[1] = exactFuncs[1]; 673c4762a1bSJed Brown } 674c4762a1bSJed Brown /* Set BC */ 6759566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 6769566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 6779566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user)); 6789566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user)); 679c4762a1bSJed Brown if (label) { 680c4762a1bSJed Brown const PetscInt id = 1; 681c4762a1bSJed Brown 6829566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, user, NULL)); 683c4762a1bSJed Brown } 6849566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &label)); 685c4762a1bSJed Brown if (label && user->useFV) { 686c4762a1bSJed Brown const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101}; 687c4762a1bSJed Brown 688dd39110bSPierre Jolivet PetscCall(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "inflow", label, PETSC_STATIC_ARRAY_LENGTH(inflowids), inflowids, 1, 0, NULL, (void (*)(void))advect_inflow, NULL, user, NULL)); 689dd39110bSPierre Jolivet PetscCall(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "outflow", label, PETSC_STATIC_ARRAY_LENGTH(outflowids), outflowids, 1, 0, NULL, (void (*)(void))advect_outflow, NULL, user, NULL)); 690c4762a1bSJed Brown } 6913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 692c4762a1bSJed Brown } 693c4762a1bSJed Brown 694d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 695d71ae5a4SJacob Faibussowitsch { 69630602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 697c4762a1bSJed Brown PetscDS prob; 69830602db0SMatthew G. Knepley PetscInt n = 3; 69930602db0SMatthew G. Knepley const char *prefix; 700c4762a1bSJed Brown 701c4762a1bSJed Brown PetscFunctionBeginUser; 7029566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 7039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *)bdt, &n, NULL)); 7049566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 705c4762a1bSJed Brown switch (user->velocityDist) { 706d71ae5a4SJacob Faibussowitsch case VEL_ZERO: 707d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u)); 708d71ae5a4SJacob Faibussowitsch break; 709c4762a1bSJed Brown case VEL_CONSTANT: 7109566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u)); 7119566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL)); 7129566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL)); 713c4762a1bSJed Brown break; 714c4762a1bSJed Brown case VEL_HARMONIC: 71530602db0SMatthew G. Knepley switch (bdt[0]) { 716c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 71730602db0SMatthew G. Knepley switch (bdt[1]) { 718d71ae5a4SJacob Faibussowitsch case DM_BOUNDARY_PERIODIC: 719d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u)); 720d71ae5a4SJacob Faibussowitsch break; 721d71ae5a4SJacob Faibussowitsch default: 722d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u)); 723d71ae5a4SJacob Faibussowitsch break; 724c4762a1bSJed Brown } 725c4762a1bSJed Brown break; 726d71ae5a4SJacob Faibussowitsch default: 727d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u)); 728d71ae5a4SJacob Faibussowitsch break; 729c4762a1bSJed Brown } 7309566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 731c4762a1bSJed Brown break; 732c4762a1bSJed Brown case VEL_SHEAR: 7339566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u)); 7349566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 735c4762a1bSJed Brown break; 736c4762a1bSJed Brown } 7379566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_advection, f1_advection)); 7389566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL)); 7399566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL)); 7409566063dSJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_advection)); 7419566063dSJacob Faibussowitsch else PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection)); 742c4762a1bSJed Brown 7439566063dSJacob Faibussowitsch PetscCall(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user)); 7443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 745c4762a1bSJed Brown } 746c4762a1bSJed Brown 747d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 748d71ae5a4SJacob Faibussowitsch { 749c4762a1bSJed Brown DM cdm = dm; 750c4762a1bSJed Brown PetscQuadrature q; 751c4762a1bSJed Brown PetscFE fe[2]; 752c4762a1bSJed Brown PetscFV fv; 753c4762a1bSJed Brown MPI_Comm comm; 75430602db0SMatthew G. Knepley PetscInt dim; 755c4762a1bSJed Brown 756c4762a1bSJed Brown PetscFunctionBeginUser; 757c4762a1bSJed Brown /* Create finite element */ 7589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 7599566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7609566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0])); 7619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 7629566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1])); 7639566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 7649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "porosity")); 765c4762a1bSJed Brown 7669566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)dm), &fv)); 7679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fv, "porosity")); 7689566063dSJacob Faibussowitsch PetscCall(PetscFVSetFromOptions(fv)); 7699566063dSJacob Faibussowitsch PetscCall(PetscFVSetNumComponents(fv, 1)); 7709566063dSJacob Faibussowitsch PetscCall(PetscFVSetSpatialDimension(fv, dim)); 7719566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe[0], &q)); 7729566063dSJacob Faibussowitsch PetscCall(PetscFVSetQuadrature(fv, q)); 773c4762a1bSJed Brown 7749566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 7759566063dSJacob Faibussowitsch if (user->useFV) PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fv)); 7769566063dSJacob Faibussowitsch else PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 7779566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 7789566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 779c4762a1bSJed Brown 780c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 781c4762a1bSJed Brown while (cdm) { 7829566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 7839566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 784c4762a1bSJed Brown /* Coordinates were never localized for coarse meshes */ 7859566063dSJacob Faibussowitsch if (cdm) PetscCall(DMLocalizeCoordinates(cdm)); 786c4762a1bSJed Brown } 7879566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 7889566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 7899566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fv)); 7903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 791c4762a1bSJed Brown } 792c4762a1bSJed Brown 793d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 794d71ae5a4SJacob Faibussowitsch { 795c4762a1bSJed Brown PetscFunctionBeginUser; 7969566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, user, dm)); 797c4762a1bSJed Brown /* Handle refinement, etc. */ 7989566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 799c4762a1bSJed Brown /* Construct ghost cells */ 800c4762a1bSJed Brown if (user->useFV) { 801c4762a1bSJed Brown DM gdm; 802c4762a1bSJed Brown 8039566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm)); 8049566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 805c4762a1bSJed Brown *dm = gdm; 806c4762a1bSJed Brown } 807c4762a1bSJed Brown /* Localize coordinates */ 8089566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 8099566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)(*dm), "Mesh")); 8109566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 811c4762a1bSJed Brown /* Setup problem */ 8129566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*dm, user)); 813c4762a1bSJed Brown /* Setup BC */ 8149566063dSJacob Faibussowitsch PetscCall(SetupBC(*dm, user)); 8153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 816c4762a1bSJed Brown } 817c4762a1bSJed Brown 818d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx) 819d71ae5a4SJacob Faibussowitsch { 820c4762a1bSJed Brown PetscDS prob; 821c4762a1bSJed Brown DM dmCell; 822c4762a1bSJed Brown Vec cellgeom; 823c4762a1bSJed Brown const PetscScalar *cgeom; 824c4762a1bSJed Brown PetscScalar *x; 825c4762a1bSJed Brown PetscInt dim, Nf, cStart, cEnd, c; 826c4762a1bSJed Brown 827c4762a1bSJed Brown PetscFunctionBeginUser; 8289566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 8299566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8309566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 8319566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8329566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 8339566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 8349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 8359566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 836c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 837c4762a1bSJed Brown PetscFVCellGeom *cg; 838c4762a1bSJed Brown PetscScalar *xc; 839c4762a1bSJed Brown 8409566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 8419566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc)); 8429566063dSJacob Faibussowitsch if (xc) PetscCall((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx)); 843c4762a1bSJed Brown } 8449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 8459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 8463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 847c4762a1bSJed Brown } 848c4762a1bSJed Brown 849d71ae5a4SJacob Faibussowitsch static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 850d71ae5a4SJacob Faibussowitsch { 851c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 852c4762a1bSJed Brown char *ftable = NULL; 853c4762a1bSJed Brown DM dm; 854c4762a1bSJed Brown PetscSection s; 855c4762a1bSJed Brown Vec cellgeom; 856c4762a1bSJed Brown const PetscScalar *x; 857c4762a1bSJed Brown PetscScalar *a; 858c4762a1bSJed Brown PetscReal *xnorms; 859c4762a1bSJed Brown PetscInt pStart, pEnd, p, Nf, f; 860c4762a1bSJed Brown 861c4762a1bSJed Brown PetscFunctionBeginUser; 8629566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X, (PetscObject)ts, "-view_solution")); 8639566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 8649566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8659566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 8669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &Nf)); 8679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 8689566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nf * 2, &xnorms)); 8699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 870c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 871c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 872c4762a1bSJed Brown PetscInt dof, cdof, d; 873c4762a1bSJed Brown 8749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 8759566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 8769566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, x, &a)); 877c4762a1bSJed Brown /* TODO Use constrained indices here */ 878c4762a1bSJed Brown for (d = 0; d < dof - cdof; ++d) xnorms[f * 2 + 0] = PetscMax(xnorms[f * 2 + 0], PetscAbsScalar(a[d])); 879c4762a1bSJed Brown for (d = 0; d < dof - cdof; ++d) xnorms[f * 2 + 1] += PetscAbsScalar(a[d]); 880c4762a1bSJed Brown } 881c4762a1bSJed Brown } 8829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 883c4762a1bSJed Brown if (stepnum >= 0) { /* No summary for final time */ 884c4762a1bSJed Brown DM dmCell, *fdm; 885c4762a1bSJed Brown Vec *fv; 886c4762a1bSJed Brown const PetscScalar *cgeom; 887c4762a1bSJed Brown PetscScalar **fx; 888c4762a1bSJed Brown PetscReal *fmin, *fmax, *fint, *ftmp, t; 889c4762a1bSJed Brown PetscInt cStart, cEnd, c, fcount, f, num; 890c4762a1bSJed Brown 891c4762a1bSJed Brown size_t ftableused, ftablealloc; 892c4762a1bSJed Brown 893c4762a1bSJed Brown /* Functionals have indices after registering, this is an upper bound */ 894c4762a1bSJed Brown fcount = user->numMonitorFuncs; 8959566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fint, fcount, &ftmp)); 8969566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fcount, &fdm, fcount, &fv, fcount, &fx)); 897c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 898c4762a1bSJed Brown PetscSection fs; 899c4762a1bSJed Brown const char *name = user->monitorFuncs[f]->name; 900c4762a1bSJed Brown 901c4762a1bSJed Brown fmin[f] = PETSC_MAX_REAL; 902c4762a1bSJed Brown fmax[f] = PETSC_MIN_REAL; 903c4762a1bSJed Brown fint[f] = 0; 904c4762a1bSJed Brown /* Make monitor vecs */ 9059566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &fdm[f])); 9069566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &num, &t)); 9079566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(fdm[f], num, t)); 9089566063dSJacob Faibussowitsch PetscCall(PetscSectionClone(s, &fs)); 9099566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 0, NULL)); 9109566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 1, name)); 9119566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(fdm[f], fs)); 9129566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&fs)); 9139566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm[f], &fv[f])); 9149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fv[f], name)); 9159566063dSJacob Faibussowitsch PetscCall(VecGetArray(fv[f], &fx[f])); 916c4762a1bSJed Brown } 9179566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 9189566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 9199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 9209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 921c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 922c4762a1bSJed Brown PetscFVCellGeom *cg; 923c4762a1bSJed Brown PetscScalar *cx; 924c4762a1bSJed Brown 9259566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 9269566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx)); 927c4762a1bSJed Brown if (!cx) continue; /* not a global cell */ 928c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 929c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 930c4762a1bSJed Brown PetscScalar *fxc; 931c4762a1bSJed Brown 9329566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc)); 933c4762a1bSJed Brown /* I need to make it easier to get interpolated values here */ 9349566063dSJacob Faibussowitsch PetscCall((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx)); 935c4762a1bSJed Brown fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 936c4762a1bSJed Brown } 937c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 938c4762a1bSJed Brown fmin[f] = PetscMin(fmin[f], ftmp[f]); 939c4762a1bSJed Brown fmax[f] = PetscMax(fmax[f], ftmp[f]); 940c4762a1bSJed Brown fint[f] += cg->volume * ftmp[f]; 941c4762a1bSJed Brown } 942c4762a1bSJed Brown } 9439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 9449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 945712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 946712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 947712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 948c4762a1bSJed Brown /* Output functional data */ 949c4762a1bSJed Brown ftablealloc = fcount * 100; 950c4762a1bSJed Brown ftableused = 0; 9519566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ftablealloc, &ftable)); 952c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 953c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 954c4762a1bSJed Brown PetscInt id = func->offset; 955c4762a1bSJed Brown char newline[] = "\n"; 956c4762a1bSJed Brown char buffer[256], *p, *prefix; 957c4762a1bSJed Brown size_t countused, len; 958c4762a1bSJed Brown 959c4762a1bSJed Brown /* Create string with functional outputs */ 960c4762a1bSJed Brown if (f % 3) { 9619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, " ", 2)); 962c4762a1bSJed Brown p = buffer + 2; 963c4762a1bSJed Brown } else if (f) { 9649566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, newline, sizeof(newline) - 1)); 965c4762a1bSJed Brown p = buffer + sizeof(newline) - 1; 966c4762a1bSJed Brown } else { 967c4762a1bSJed Brown p = buffer; 968c4762a1bSJed Brown } 9699566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(p, sizeof buffer - (p - buffer), "%12s [%12.6g,%12.6g] int %12.6g", &countused, func->name, (double)fmin[id], (double)fmax[id], (double)fint[id])); 970c4762a1bSJed Brown countused += p - buffer; 971c4762a1bSJed Brown /* reallocate */ 972c4762a1bSJed Brown if (countused > ftablealloc - ftableused - 1) { 973c4762a1bSJed Brown char *ftablenew; 974c4762a1bSJed Brown 975c4762a1bSJed Brown ftablealloc = 2 * ftablealloc + countused; 9769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ftablealloc, &ftablenew)); 9779566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); 9789566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 979c4762a1bSJed Brown ftable = ftablenew; 980c4762a1bSJed Brown } 9819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused)); 982c4762a1bSJed Brown ftableused += countused; 983c4762a1bSJed Brown ftable[ftableused] = 0; 984c4762a1bSJed Brown /* Output vecs */ 9859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fv[f], &fx[f])); 9869566063dSJacob Faibussowitsch PetscCall(PetscStrlen(func->name, &len)); 9879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len + 2, &prefix)); 988c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(prefix, func->name, len + 2)); 989c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(prefix, "_", len + 2)); 9909566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix)); 9919566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fv[f], NULL, "-vec_view")); 9929566063dSJacob Faibussowitsch PetscCall(PetscFree(prefix)); 9939566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm[f], &fv[f])); 9949566063dSJacob Faibussowitsch PetscCall(DMDestroy(&fdm[f])); 995c4762a1bSJed Brown } 9969566063dSJacob Faibussowitsch PetscCall(PetscFree4(fmin, fmax, fint, ftmp)); 9979566063dSJacob Faibussowitsch PetscCall(PetscFree3(fdm, fv, fx)); 99863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT " time %8.4g |x| (", stepnum, (double)time)); 999c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10009566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ", ")); 10019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%8.4g", (double)xnorms[f * 2 + 0])); 1002c4762a1bSJed Brown } 10039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ") |x|_1 (")); 1004c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10059566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ", ")); 10069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%8.4g", (double)xnorms[f * 2 + 1])); 1007c4762a1bSJed Brown } 10089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ") %s\n", ftable ? ftable : "")); 10099566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 1010c4762a1bSJed Brown } 10119566063dSJacob Faibussowitsch PetscCall(PetscFree(xnorms)); 10123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1013c4762a1bSJed Brown } 1014c4762a1bSJed Brown 1015d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 1016d71ae5a4SJacob Faibussowitsch { 1017c4762a1bSJed Brown MPI_Comm comm; 1018c4762a1bSJed Brown TS ts; 1019c4762a1bSJed Brown DM dm; 1020c4762a1bSJed Brown Vec u; 1021c4762a1bSJed Brown AppCtx user; 1022c4762a1bSJed Brown PetscReal t0, t = 0.0; 1023c4762a1bSJed Brown void *ctxs[2]; 1024c4762a1bSJed Brown 1025c4762a1bSJed Brown ctxs[0] = &t; 1026c4762a1bSJed Brown ctxs[1] = &t; 1027327415f7SBarry Smith PetscFunctionBeginUser; 10289566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1029c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1030c4762a1bSJed Brown user.functionalRegistry = NULL; 1031c4762a1bSJed Brown globalUser = &user; 10329566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 10339566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 10349566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 10359566063dSJacob Faibussowitsch PetscCall(CreateDM(comm, &user, &dm)); 10369566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 10379566063dSJacob Faibussowitsch PetscCall(ProcessMonitorOptions(comm, &user)); 1038c4762a1bSJed Brown 10399566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 10409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 1041c4762a1bSJed Brown if (user.useFV) { 1042c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 1043c4762a1bSJed Brown 10449566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, "", "-use_implicit", &isImplicit)); 1045c4762a1bSJed Brown if (isImplicit) { 10469566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10479566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1048c4762a1bSJed Brown } 10499566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10509566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user)); 1051c4762a1bSJed Brown } else { 10529566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10539566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10549566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1055c4762a1bSJed Brown } 10569566063dSJacob Faibussowitsch if (user.useFV) PetscCall(TSMonitorSet(ts, MonitorFunctionals, &user, NULL)); 10579566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1)); 10589566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2.0)); 10599566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01)); 10609566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 10619566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1062c4762a1bSJed Brown 10639566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u)); 10649566063dSJacob Faibussowitsch if (user.useFV) PetscCall(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1])); 10659566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-init_vec_view")); 10669566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1067c4762a1bSJed Brown t0 = t; 10689566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 10699566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 10709566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 10719566063dSJacob Faibussowitsch if (t > t0) PetscCall(DMTSCheckFromOptions(ts, u)); 10729566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 1073c4762a1bSJed Brown { 1074c4762a1bSJed Brown PetscReal ftime; 1075c4762a1bSJed Brown PetscInt nsteps; 1076c4762a1bSJed Brown TSConvergedReason reason; 1077c4762a1bSJed Brown 10789566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 10799566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &nsteps)); 10809566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 108163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps)); 1082c4762a1bSJed Brown } 1083c4762a1bSJed Brown 10849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 10859566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 10869566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 10879566063dSJacob Faibussowitsch PetscCall(PetscFree(user.monitorFuncs)); 10889566063dSJacob Faibussowitsch PetscCall(FunctionalDestroy(&user.functionalRegistry)); 10899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1090b122ec5aSJacob Faibussowitsch return 0; 1091c4762a1bSJed Brown } 1092c4762a1bSJed Brown 1093c4762a1bSJed Brown /*TEST 1094c4762a1bSJed Brown 109530602db0SMatthew G. Knepley testset: 109630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 109730602db0SMatthew G. Knepley 1098c4762a1bSJed Brown # 2D harmonic velocity, no porosity 1099c4762a1bSJed Brown test: 1100c4762a1bSJed Brown suffix: p1p1 1101c4762a1bSJed Brown requires: !complex !single 1102*b4f8a55aSStefano Zampini args: -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_type ilu -pc_factor_shift_type nonzero -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1103c4762a1bSJed Brown test: 1104c4762a1bSJed Brown suffix: p1p1_xper 1105c4762a1bSJed Brown requires: !complex !single 110630602db0SMatthew G. Knepley 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 1107c4762a1bSJed Brown test: 1108c4762a1bSJed Brown suffix: p1p1_xper_ref 1109c4762a1bSJed Brown requires: !complex !single 111030602db0SMatthew G. Knepley 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 1111c4762a1bSJed Brown test: 1112c4762a1bSJed Brown suffix: p1p1_xyper 1113c4762a1bSJed Brown requires: !complex !single 111430602db0SMatthew G. Knepley 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 1115c4762a1bSJed Brown test: 1116c4762a1bSJed Brown suffix: p1p1_xyper_ref 1117c4762a1bSJed Brown requires: !complex !single 111830602db0SMatthew G. Knepley 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 1119c4762a1bSJed Brown test: 1120c4762a1bSJed Brown suffix: p2p1 1121c4762a1bSJed Brown requires: !complex !single 112230602db0SMatthew G. Knepley 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 1123c4762a1bSJed Brown test: 1124c4762a1bSJed Brown suffix: p2p1_xyper 1125c4762a1bSJed Brown requires: !complex !single 112630602db0SMatthew G. Knepley 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 112730602db0SMatthew G. Knepley 112830602db0SMatthew G. Knepley test: 112930602db0SMatthew G. Knepley suffix: adv_1 113030602db0SMatthew G. Knepley requires: !complex !single 113130602db0SMatthew G. Knepley 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 113230602db0SMatthew G. Knepley 113330602db0SMatthew G. Knepley test: 113430602db0SMatthew G. Knepley suffix: adv_2 113530602db0SMatthew G. Knepley requires: !complex 113630602db0SMatthew G. Knepley TODO: broken memory corruption 113730602db0SMatthew G. Knepley 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 113830602db0SMatthew G. Knepley 113930602db0SMatthew G. Knepley test: 114030602db0SMatthew G. Knepley suffix: adv_3 114130602db0SMatthew G. Knepley requires: !complex 114230602db0SMatthew G. Knepley TODO: broken memory corruption 114330602db0SMatthew G. Knepley 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 114430602db0SMatthew G. Knepley 114530602db0SMatthew G. Knepley test: 114630602db0SMatthew G. Knepley suffix: adv_3_ex 114730602db0SMatthew G. Knepley requires: !complex 114830602db0SMatthew G. Knepley 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 114930602db0SMatthew G. Knepley 115030602db0SMatthew G. Knepley test: 115130602db0SMatthew G. Knepley suffix: adv_4 115230602db0SMatthew G. Knepley requires: !complex 115330602db0SMatthew G. Knepley TODO: broken memory corruption 115430602db0SMatthew G. Knepley 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 115530602db0SMatthew G. Knepley 115630602db0SMatthew G. Knepley # 2D Advection, box, delta 115730602db0SMatthew G. Knepley test: 115830602db0SMatthew G. Knepley suffix: adv_delta_yper_0 115930602db0SMatthew G. Knepley requires: !complex 116030602db0SMatthew G. Knepley TODO: broken 116130602db0SMatthew G. Knepley 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 116230602db0SMatthew G. Knepley 116330602db0SMatthew G. Knepley test: 116430602db0SMatthew G. Knepley suffix: adv_delta_yper_1 116530602db0SMatthew G. Knepley requires: !complex 116630602db0SMatthew G. Knepley TODO: broken 116730602db0SMatthew G. Knepley 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 116830602db0SMatthew G. Knepley 116930602db0SMatthew G. Knepley test: 117030602db0SMatthew G. Knepley suffix: adv_delta_yper_2 117130602db0SMatthew G. Knepley requires: !complex 117230602db0SMatthew G. Knepley TODO: broken 117330602db0SMatthew G. Knepley 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 117430602db0SMatthew G. Knepley 117530602db0SMatthew G. Knepley test: 117630602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_0 117730602db0SMatthew G. Knepley requires: !complex 117830602db0SMatthew G. Knepley TODO: broken 117930602db0SMatthew G. Knepley 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 118030602db0SMatthew G. Knepley 118130602db0SMatthew G. Knepley test: 118230602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_1 118330602db0SMatthew G. Knepley requires: !complex 118430602db0SMatthew G. Knepley TODO: broken 118530602db0SMatthew G. Knepley 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 118630602db0SMatthew G. Knepley 118730602db0SMatthew G. Knepley test: 118830602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_2 118930602db0SMatthew G. Knepley requires: !complex 119030602db0SMatthew G. Knepley TODO: broken 119130602db0SMatthew G. Knepley 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 119230602db0SMatthew G. Knepley 119330602db0SMatthew G. Knepley test: 119430602db0SMatthew G. Knepley suffix: adv_delta_yper_im_0 119530602db0SMatthew G. Knepley requires: !complex 119630602db0SMatthew G. Knepley TODO: broken 119730602db0SMatthew G. Knepley 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 119830602db0SMatthew G. Knepley 119930602db0SMatthew G. Knepley test: 120030602db0SMatthew G. Knepley suffix: adv_delta_yper_im_1 120130602db0SMatthew G. Knepley requires: !complex 120230602db0SMatthew G. Knepley TODO: broken 120330602db0SMatthew G. Knepley 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 120430602db0SMatthew G. Knepley 120530602db0SMatthew G. Knepley test: 120630602db0SMatthew G. Knepley suffix: adv_delta_yper_im_2 120730602db0SMatthew G. Knepley requires: !complex 120830602db0SMatthew G. Knepley TODO: broken 120930602db0SMatthew G. Knepley 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 121030602db0SMatthew G. Knepley 121130602db0SMatthew G. Knepley test: 121230602db0SMatthew G. Knepley suffix: adv_delta_yper_im_3 121330602db0SMatthew G. Knepley requires: !complex 121430602db0SMatthew G. Knepley TODO: broken 121530602db0SMatthew G. Knepley 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 121630602db0SMatthew G. Knepley 121730602db0SMatthew G. Knepley # I believe the nullspace is sin(pi y) 121830602db0SMatthew G. Knepley test: 121930602db0SMatthew G. Knepley suffix: adv_delta_yper_im_4 122030602db0SMatthew G. Knepley requires: !complex 122130602db0SMatthew G. Knepley TODO: broken 122230602db0SMatthew G. Knepley 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 122330602db0SMatthew G. Knepley 122430602db0SMatthew G. Knepley test: 122530602db0SMatthew G. Knepley suffix: adv_delta_yper_im_5 122630602db0SMatthew G. Knepley requires: !complex 122730602db0SMatthew G. Knepley TODO: broken 122830602db0SMatthew G. Knepley 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 122930602db0SMatthew G. Knepley 123030602db0SMatthew G. Knepley test: 123130602db0SMatthew G. Knepley suffix: adv_delta_yper_im_6 123230602db0SMatthew G. Knepley requires: !complex 123330602db0SMatthew G. Knepley TODO: broken 123430602db0SMatthew G. Knepley 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 123530602db0SMatthew G. Knepley # 2D Advection, magma benchmark 1 123630602db0SMatthew G. Knepley 123730602db0SMatthew G. Knepley test: 123830602db0SMatthew G. Knepley suffix: adv_delta_shear_im_0 123930602db0SMatthew G. Knepley requires: !complex 124030602db0SMatthew G. Knepley TODO: broken 124130602db0SMatthew G. Knepley 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 124230602db0SMatthew G. Knepley # 2D Advection, box, gaussian 124330602db0SMatthew G. Knepley 124430602db0SMatthew G. Knepley test: 124530602db0SMatthew G. Knepley suffix: adv_gauss 124630602db0SMatthew G. Knepley requires: !complex 124730602db0SMatthew G. Knepley TODO: broken 124830602db0SMatthew G. Knepley 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 124930602db0SMatthew G. Knepley 125030602db0SMatthew G. Knepley test: 125130602db0SMatthew G. Knepley suffix: adv_gauss_im 125230602db0SMatthew G. Knepley requires: !complex 125330602db0SMatthew G. Knepley TODO: broken 125430602db0SMatthew G. Knepley 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 125530602db0SMatthew G. Knepley 125630602db0SMatthew G. Knepley test: 125730602db0SMatthew G. Knepley suffix: adv_gauss_im_1 125830602db0SMatthew G. Knepley requires: !complex 125930602db0SMatthew G. Knepley TODO: broken 126030602db0SMatthew G. Knepley 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 126130602db0SMatthew G. Knepley 126230602db0SMatthew G. Knepley test: 126330602db0SMatthew G. Knepley suffix: adv_gauss_im_2 126430602db0SMatthew G. Knepley requires: !complex 126530602db0SMatthew G. Knepley TODO: broken 126630602db0SMatthew G. Knepley 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 126730602db0SMatthew G. Knepley 126830602db0SMatthew G. Knepley test: 126930602db0SMatthew G. Knepley suffix: adv_gauss_corner 127030602db0SMatthew G. Knepley requires: !complex 127130602db0SMatthew G. Knepley TODO: broken 127230602db0SMatthew G. Knepley 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 127330602db0SMatthew G. Knepley 127430602db0SMatthew G. Knepley # 2D Advection+Harmonic 12- 127530602db0SMatthew G. Knepley test: 127630602db0SMatthew G. Knepley suffix: adv_harm_0 127730602db0SMatthew G. Knepley requires: !complex 127830602db0SMatthew G. Knepley TODO: broken memory corruption 127930602db0SMatthew G. Knepley 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 128030602db0SMatthew G. Knepley 1281c4762a1bSJed Brown # Must check that FV BCs propagate to coarse meshes 1282c4762a1bSJed Brown # Must check that FV BC ids propagate to coarse meshes 1283c4762a1bSJed Brown # Must check that FE+FV BCs work at the same time 1284c4762a1bSJed Brown # 2D Advection, matching wind in ex11 8-11 1285c4762a1bSJed Brown # NOTE implicit solves are limited by accuracy of FD Jacobian 1286c4762a1bSJed Brown test: 1287c4762a1bSJed Brown suffix: adv_0 1288c4762a1bSJed Brown requires: !complex !single exodusii 128930602db0SMatthew G. Knepley 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 1290c4762a1bSJed Brown 1291c4762a1bSJed Brown test: 1292c4762a1bSJed Brown suffix: adv_0_im 1293c4762a1bSJed Brown requires: !complex exodusii 1294c4762a1bSJed Brown TODO: broken memory corruption 129530602db0SMatthew G. Knepley 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 1296c4762a1bSJed Brown 1297c4762a1bSJed Brown test: 1298c4762a1bSJed Brown suffix: adv_0_im_2 1299c4762a1bSJed Brown requires: !complex exodusii 1300c4762a1bSJed Brown TODO: broken 130130602db0SMatthew G. Knepley 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 1302c4762a1bSJed Brown 1303c4762a1bSJed Brown test: 1304c4762a1bSJed Brown suffix: adv_0_im_3 1305c4762a1bSJed Brown requires: !complex exodusii 1306c4762a1bSJed Brown TODO: broken 130730602db0SMatthew G. Knepley 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 1308c4762a1bSJed Brown 1309c4762a1bSJed Brown test: 1310c4762a1bSJed Brown suffix: adv_0_im_4 1311c4762a1bSJed Brown requires: !complex exodusii 1312c4762a1bSJed Brown TODO: broken 131330602db0SMatthew G. Knepley 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 1314c4762a1bSJed Brown # 2D Advection, misc 1315c4762a1bSJed Brown 1316c4762a1bSJed Brown TEST*/ 1317