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(); 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) 110d71ae5a4SJacob Faibussowitsch { 111c4762a1bSJed Brown Functional func; 112c4762a1bSJed Brown char *names[256]; 113c4762a1bSJed Brown PetscInt f; 114c4762a1bSJed Brown 115c4762a1bSJed Brown PetscFunctionBeginUser; 116d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX"); 117dd39110bSPierre Jolivet options->numMonitorFuncs = PETSC_STATIC_ARRAY_LENGTH(names); 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL)); 1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs)); 120c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 121c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 122c4762a1bSJed Brown PetscBool match; 123c4762a1bSJed Brown 1249566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(names[f], func->name, &match)); 125c4762a1bSJed Brown if (match) break; 126c4762a1bSJed Brown } 1273c633725SBarry Smith PetscCheck(func, comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 128c4762a1bSJed Brown options->monitorFuncs[f] = func; 129c4762a1bSJed Brown /* Jed inserts a de-duplication of functionals here */ 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(names[f])); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 133c4762a1bSJed Brown options->maxMonitorFunc = -1; 134c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 135c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 136c4762a1bSJed Brown Functional call = options->monitorFuncs[f]; 137c4762a1bSJed Brown 138c4762a1bSJed Brown if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown } 141d0609cedSBarry Smith PetscOptionsEnd(); 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145d71ae5a4SJacob Faibussowitsch static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) 146d71ae5a4SJacob Faibussowitsch { 147c4762a1bSJed Brown Functional *ptr, f; 148c4762a1bSJed Brown PetscInt lastoffset = -1; 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscFunctionBeginUser; 151c4762a1bSJed Brown for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 1529566063dSJacob Faibussowitsch PetscCall(PetscNew(&f)); 1539566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &f->name)); 154c4762a1bSJed Brown f->offset = lastoffset + 1; 155c4762a1bSJed Brown f->func = func; 156c4762a1bSJed Brown f->ctx = ctx; 157c4762a1bSJed Brown f->next = NULL; 158c4762a1bSJed Brown *ptr = f; 159c4762a1bSJed Brown *offset = f->offset; 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163d71ae5a4SJacob Faibussowitsch static PetscErrorCode FunctionalDestroy(Functional *link) 164d71ae5a4SJacob Faibussowitsch { 165c4762a1bSJed Brown Functional next, l; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBeginUser; 1683ba16761SJacob Faibussowitsch if (!link) PetscFunctionReturn(PETSC_SUCCESS); 169c4762a1bSJed Brown l = *link; 170c4762a1bSJed Brown *link = NULL; 171c4762a1bSJed Brown for (; l; l = next) { 172c4762a1bSJed Brown next = l->next; 1739566063dSJacob Faibussowitsch PetscCall(PetscFree(l->name)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFree(l)); 175c4762a1bSJed Brown } 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179d71ae5a4SJacob 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[]) 180d71ae5a4SJacob Faibussowitsch { 181c4762a1bSJed Brown PetscInt comp; 182c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185d71ae5a4SJacob 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[]) 186d71ae5a4SJacob Faibussowitsch { 187c4762a1bSJed Brown PetscScalar wind[3] = {0.0, 0.0, 0.0}; 188c4762a1bSJed Brown PetscInt comp; 189c4762a1bSJed Brown 1903ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, constant_u_2d(dim, t, x, Nf, wind, NULL)); 191c4762a1bSJed Brown for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp]; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194d71ae5a4SJacob 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[]) 195d71ae5a4SJacob Faibussowitsch { 196c4762a1bSJed Brown PetscInt comp; 197c4762a1bSJed Brown for (comp = 0; comp < dim * dim; ++comp) f1[comp] = 0.0; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown 200d71ae5a4SJacob 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[]) 201d71ae5a4SJacob Faibussowitsch { 202c4762a1bSJed Brown PetscInt d; 203c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d * dim + d] = 1.0; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206d71ae5a4SJacob 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[]) 207d71ae5a4SJacob Faibussowitsch { 208c4762a1bSJed Brown g0[0] = 1.0; 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211d71ae5a4SJacob 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[]) 212d71ae5a4SJacob Faibussowitsch { 213c4762a1bSJed Brown PetscInt comp; 214c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217d71ae5a4SJacob 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[]) 218d71ae5a4SJacob Faibussowitsch { 219c4762a1bSJed Brown PetscInt comp, d; 220c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) { 221ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f1[comp * dim + d] = u_x[comp * dim + d]; 222c4762a1bSJed Brown } 223c4762a1bSJed Brown } 224c4762a1bSJed Brown 225d71ae5a4SJacob 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[]) 226d71ae5a4SJacob Faibussowitsch { 227c4762a1bSJed Brown f0[0] = -PetscSinReal(2.0 * PETSC_PI * x[0]); 228c4762a1bSJed Brown f0[1] = 2.0 * PETSC_PI * x[1] * PetscCosReal(2.0 * PETSC_PI * x[0]); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231d71ae5a4SJacob 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[]) 232d71ae5a4SJacob Faibussowitsch { 233c4762a1bSJed Brown f0[0] = -2.0 * PetscSinReal(2.0 * PETSC_PI * x[0]) * PetscCosReal(2.0 * PETSC_PI * x[1]); 234c4762a1bSJed Brown f0[1] = 2.0 * PetscSinReal(2.0 * PETSC_PI * x[1]) * PetscCosReal(2.0 * PETSC_PI * x[0]); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237d71ae5a4SJacob 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[]) 238d71ae5a4SJacob Faibussowitsch { 239c4762a1bSJed Brown const PetscInt Ncomp = dim; 240c4762a1bSJed Brown PetscInt compI, d; 241c4762a1bSJed Brown 242c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 243ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g3[((compI * Ncomp + compI) * dim + d) * dim + d] = 1.0; 244c4762a1bSJed Brown } 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 248d71ae5a4SJacob 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[]) 249d71ae5a4SJacob Faibussowitsch { 250c4762a1bSJed Brown PetscInt d; 251c4762a1bSJed Brown f0[0] = u_t[dim]; 252c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += u[dim] * u_x[d * dim + d] + u_x[dim * dim + d] * u[d]; 253c4762a1bSJed Brown } 254c4762a1bSJed Brown 255d71ae5a4SJacob 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[]) 256d71ae5a4SJacob Faibussowitsch { 257c4762a1bSJed Brown PetscInt d; 258c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[0] = 0.0; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261d71ae5a4SJacob 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[]) 262d71ae5a4SJacob Faibussowitsch { 263c4762a1bSJed Brown PetscInt d; 264c4762a1bSJed Brown g0[0] = u_tShift; 265c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[d * dim + d]; 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268d71ae5a4SJacob 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[]) 269d71ae5a4SJacob Faibussowitsch { 270c4762a1bSJed Brown PetscInt d; 271c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = u[d]; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274d71ae5a4SJacob 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[]) 275d71ae5a4SJacob Faibussowitsch { 276c4762a1bSJed Brown PetscInt d; 277c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[dim * dim + d]; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280d71ae5a4SJacob 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[]) 281d71ae5a4SJacob Faibussowitsch { 282c4762a1bSJed Brown PetscInt d; 283c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = u[dim]; 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286d71ae5a4SJacob 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) 287d71ae5a4SJacob Faibussowitsch { 288c4762a1bSJed Brown PetscReal wind[3] = {0.0, 1.0, 0.0}; 289c4762a1bSJed Brown PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim, 3), wind, n); 290c4762a1bSJed Brown 291c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 292c4762a1bSJed Brown } 293c4762a1bSJed Brown 294d71ae5a4SJacob 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) 295d71ae5a4SJacob Faibussowitsch { 296c4762a1bSJed Brown PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 297c4762a1bSJed Brown 298c4762a1bSJed Brown #if 1 299c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 300c4762a1bSJed Brown #else 301c4762a1bSJed 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]); */ 302c4762a1bSJed Brown /* Smear it out */ 303c4762a1bSJed Brown flux[0] = 0.5 * ((uL[dim] + uR[dim]) + (uL[dim] - uR[dim]) * tanh(1.0e5 * wn)) * wn; 304c4762a1bSJed Brown #endif 305c4762a1bSJed Brown } 306c4762a1bSJed Brown 307d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 308d71ae5a4SJacob Faibussowitsch { 309c4762a1bSJed Brown u[0] = 0.0; 310c4762a1bSJed Brown u[1] = 0.0; 3113ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 312c4762a1bSJed Brown } 313c4762a1bSJed Brown 314d71ae5a4SJacob Faibussowitsch static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 315d71ae5a4SJacob Faibussowitsch { 316c4762a1bSJed Brown u[0] = 0.0; 317c4762a1bSJed Brown u[1] = 1.0; 3183ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 319c4762a1bSJed Brown } 320c4762a1bSJed Brown 321c4762a1bSJed Brown /* Coordinates of the point which was at x at t = 0 */ 322d71ae5a4SJacob Faibussowitsch static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 323d71ae5a4SJacob Faibussowitsch { 324c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 325c4762a1bSJed Brown u[0] = x[0]; 326c4762a1bSJed Brown u[1] = x[1] + t; 32730602db0SMatthew G. Knepley #if 0 3289566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u)); 32930602db0SMatthew G. Knepley #else 33030602db0SMatthew G. Knepley u[1] = u[1] - (int)PetscRealPart(u[1]); 33130602db0SMatthew G. Knepley #endif 3323ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 335c4762a1bSJed Brown /* 336c4762a1bSJed Brown In 2D we use the exact solution: 337c4762a1bSJed Brown 338c4762a1bSJed Brown u = x^2 + y^2 339c4762a1bSJed Brown v = 2 x^2 - 2xy 340c4762a1bSJed Brown phi = h(x + y + (u + v) t) 341c4762a1bSJed Brown f_x = f_y = 4 342c4762a1bSJed Brown 343c4762a1bSJed Brown so that 344c4762a1bSJed Brown 345c4762a1bSJed Brown -\Delta u + f = <-4, -4> + <4, 4> = 0 346c4762a1bSJed Brown {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 347c4762a1bSJed Brown h_t(x + y + (u + v) t) - u . grad phi - phi div u 348c4762a1bSJed Brown = u h' + v h' - u h_x - v h_y 349c4762a1bSJed Brown = 0 350c4762a1bSJed Brown 351c4762a1bSJed Brown We will conserve phi since 352c4762a1bSJed Brown 353c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 354c4762a1bSJed Brown 355c4762a1bSJed Brown Also try h((x + ut)^2 + (y + vt)^2), so that 356c4762a1bSJed Brown 357c4762a1bSJed Brown h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 358c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 359c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 360c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 361c4762a1bSJed Brown = 0 362c4762a1bSJed Brown 363c4762a1bSJed Brown */ 364d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 365d71ae5a4SJacob Faibussowitsch { 366c4762a1bSJed Brown u[0] = x[0] * x[0] + x[1] * x[1]; 367c4762a1bSJed Brown u[1] = 2.0 * x[0] * x[0] - 2.0 * x[0] * x[1]; 3683ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown 371c4762a1bSJed Brown /* 372c4762a1bSJed Brown In 2D we use the exact, periodic solution: 373c4762a1bSJed Brown 374c4762a1bSJed Brown u = sin(2 pi x)/4 pi^2 375c4762a1bSJed Brown v = -y cos(2 pi x)/2 pi 376c4762a1bSJed Brown phi = h(x + y + (u + v) t) 377c4762a1bSJed Brown f_x = -sin(2 pi x) 378c4762a1bSJed Brown f_y = 2 pi y cos(2 pi x) 379c4762a1bSJed Brown 380c4762a1bSJed Brown so that 381c4762a1bSJed Brown 382c4762a1bSJed Brown -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 383c4762a1bSJed Brown 384c4762a1bSJed Brown We will conserve phi since 385c4762a1bSJed Brown 386c4762a1bSJed Brown \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 387c4762a1bSJed Brown */ 388d71ae5a4SJacob Faibussowitsch static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 389d71ae5a4SJacob Faibussowitsch { 390c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]) / PetscSqr(2.0 * PETSC_PI); 391c4762a1bSJed Brown u[1] = -x[1] * PetscCosReal(2.0 * PETSC_PI * x[0]) / (2.0 * PETSC_PI); 3923ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 393c4762a1bSJed Brown } 394c4762a1bSJed Brown 395c4762a1bSJed Brown /* 396c4762a1bSJed Brown In 2D we use the exact, doubly periodic solution: 397c4762a1bSJed Brown 398c4762a1bSJed Brown u = sin(2 pi x) cos(2 pi y)/4 pi^2 399c4762a1bSJed Brown v = -sin(2 pi y) cos(2 pi x)/4 pi^2 400c4762a1bSJed Brown phi = h(x + y + (u + v) t) 401c4762a1bSJed Brown f_x = -2sin(2 pi x) cos(2 pi y) 402c4762a1bSJed Brown f_y = 2sin(2 pi y) cos(2 pi x) 403c4762a1bSJed Brown 404c4762a1bSJed Brown so that 405c4762a1bSJed Brown 406c4762a1bSJed 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 407c4762a1bSJed Brown 408c4762a1bSJed Brown We will conserve phi since 409c4762a1bSJed Brown 410c4762a1bSJed Brown \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 411c4762a1bSJed Brown */ 412d71ae5a4SJacob Faibussowitsch static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 413d71ae5a4SJacob Faibussowitsch { 414c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]) * PetscCosReal(2.0 * PETSC_PI * x[1]) / PetscSqr(2.0 * PETSC_PI); 415c4762a1bSJed Brown u[1] = -PetscSinReal(2.0 * PETSC_PI * x[1]) * PetscCosReal(2.0 * PETSC_PI * x[0]) / PetscSqr(2.0 * PETSC_PI); 4163ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 417c4762a1bSJed Brown } 418c4762a1bSJed Brown 419d71ae5a4SJacob Faibussowitsch static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 420d71ae5a4SJacob Faibussowitsch { 421c4762a1bSJed Brown u[0] = x[1] - 0.5; 422c4762a1bSJed Brown u[1] = 0.0; 4233ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 424c4762a1bSJed Brown } 425c4762a1bSJed Brown 426d71ae5a4SJacob Faibussowitsch static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 427d71ae5a4SJacob Faibussowitsch { 428c4762a1bSJed Brown PetscInt d; 429c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 4303ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 431c4762a1bSJed Brown } 432c4762a1bSJed Brown 433d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 434d71ae5a4SJacob Faibussowitsch { 435c4762a1bSJed Brown u[0] = 0.0; 4363ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 439d71ae5a4SJacob Faibussowitsch static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 440d71ae5a4SJacob Faibussowitsch { 441c4762a1bSJed Brown u[0] = 1.0; 4423ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 443c4762a1bSJed Brown } 444c4762a1bSJed Brown 445d71ae5a4SJacob Faibussowitsch static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 446d71ae5a4SJacob Faibussowitsch { 447c4762a1bSJed Brown PetscReal x0[2]; 448c4762a1bSJed Brown PetscScalar xn[2]; 449c4762a1bSJed Brown 450c4762a1bSJed Brown x0[0] = globalUser->source[0]; 451c4762a1bSJed Brown x0[1] = globalUser->source[1]; 4523ba16761SJacob Faibussowitsch PetscCall(constant_x_2d(dim, time, x0, Nf, xn, ctx)); 453c4762a1bSJed Brown { 454c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 455c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 456c4762a1bSJed Brown const PetscReal r2 = xi * xi + eta * eta; 457c4762a1bSJed Brown 458c4762a1bSJed Brown u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 459c4762a1bSJed Brown } 4603ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 461c4762a1bSJed Brown } 462c4762a1bSJed Brown 463c4762a1bSJed Brown /* 464c4762a1bSJed Brown Gaussian blob, initially centered on (0.5, 0.5) 465c4762a1bSJed Brown 466c4762a1bSJed Brown xi = x(t) - x0, eta = y(t) - y0 467c4762a1bSJed Brown 468c4762a1bSJed Brown where x(t), y(t) are the integral curves of v(t), 469c4762a1bSJed Brown 470c4762a1bSJed Brown dx/dt . grad f = v . f 471c4762a1bSJed Brown 472c4762a1bSJed Brown Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 473c4762a1bSJed Brown 474c4762a1bSJed Brown v0 f_x + w0 f_y = v . f 475c4762a1bSJed Brown */ 476d71ae5a4SJacob Faibussowitsch static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 477d71ae5a4SJacob Faibussowitsch { 478c4762a1bSJed Brown const PetscReal x0[2] = {0.5, 0.5}; 479c4762a1bSJed Brown const PetscReal sigma = 1.0 / 6.0; 480c4762a1bSJed Brown PetscScalar xn[2]; 481c4762a1bSJed Brown 4823ba16761SJacob Faibussowitsch PetscCall(constant_x_2d(dim, time, x0, Nf, xn, ctx)); 483c4762a1bSJed Brown { 484c4762a1bSJed Brown /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 485c4762a1bSJed Brown /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 486c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 487c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 488c4762a1bSJed Brown const PetscReal r2 = xi * xi + eta * eta; 489c4762a1bSJed Brown 490c4762a1bSJed Brown u[0] = PetscExpReal(-r2 / (2.0 * sigma * sigma)) / (sigma * PetscSqrtReal(2.0 * PETSC_PI)); 491c4762a1bSJed Brown } 4923ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 493c4762a1bSJed Brown } 494c4762a1bSJed Brown 495d71ae5a4SJacob Faibussowitsch static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 496d71ae5a4SJacob Faibussowitsch { 497c4762a1bSJed Brown PetscReal x0[3]; 498c4762a1bSJed Brown const PetscReal wind[3] = {0.0, 1.0, 0.0}; 499c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 500c4762a1bSJed Brown 501c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 502c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0 * x[0] + 3.0 * x[1]; 503c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 5043ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 505c4762a1bSJed Brown } 506c4762a1bSJed Brown 507d71ae5a4SJacob Faibussowitsch static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 508d71ae5a4SJacob Faibussowitsch { 509c4762a1bSJed Brown PetscReal ur[3]; 510c4762a1bSJed Brown PetscReal x0[3]; 511c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 512c4762a1bSJed Brown 5139371c9d4SSatish Balay ur[0] = PetscRealPart(u[0]); 5149371c9d4SSatish Balay ur[1] = PetscRealPart(u[1]); 5159371c9d4SSatish Balay ur[2] = PetscRealPart(u[2]); 516c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 517c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0 * x[0] + 3.0 * x[1]; 518c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 5193ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 520c4762a1bSJed Brown } 521c4762a1bSJed Brown 522d71ae5a4SJacob Faibussowitsch static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 523d71ae5a4SJacob Faibussowitsch { 524c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 525c4762a1bSJed Brown 526c4762a1bSJed Brown PetscFunctionBeginUser; 527c4762a1bSJed Brown xG[0] = user->inflowState; 5283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 529c4762a1bSJed Brown } 530c4762a1bSJed Brown 531d71ae5a4SJacob Faibussowitsch static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 532d71ae5a4SJacob Faibussowitsch { 533c4762a1bSJed Brown PetscFunctionBeginUser; 53430602db0SMatthew G. Knepley //xG[0] = xI[dim]; 53530602db0SMatthew G. Knepley xG[0] = xI[2]; 5363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 537c4762a1bSJed Brown } 538c4762a1bSJed Brown 539d71ae5a4SJacob Faibussowitsch static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 540d71ae5a4SJacob Faibussowitsch { 541c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 542c4762a1bSJed Brown PetscInt dim; 543c4762a1bSJed Brown 544c4762a1bSJed Brown PetscFunctionBeginUser; 5459566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 546c4762a1bSJed Brown switch (user->porosityDist) { 547c4762a1bSJed Brown case TILTED: 5483ba16761SJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) PetscCall(tilted_phi_2d(dim, time, x, 2, u, (void *)&time)); 5493ba16761SJacob Faibussowitsch else PetscCall(tilted_phi_coupled_2d(dim, time, x, 2, u, (void *)&time)); 550c4762a1bSJed Brown break; 551d71ae5a4SJacob Faibussowitsch case GAUSSIAN: 5523ba16761SJacob Faibussowitsch PetscCall(gaussian_phi_2d(dim, time, x, 2, u, (void *)&time)); 553d71ae5a4SJacob Faibussowitsch break; 554d71ae5a4SJacob Faibussowitsch case DELTA: 5553ba16761SJacob Faibussowitsch PetscCall(delta_phi_2d(dim, time, x, 2, u, (void *)&time)); 556d71ae5a4SJacob Faibussowitsch break; 557d71ae5a4SJacob Faibussowitsch default: 558d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 559c4762a1bSJed Brown } 5603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 561c4762a1bSJed Brown } 562c4762a1bSJed Brown 563d71ae5a4SJacob Faibussowitsch static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 564d71ae5a4SJacob Faibussowitsch { 565c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 566c4762a1bSJed Brown PetscScalar yexact[3] = {0, 0, 0}; 567c4762a1bSJed Brown 568c4762a1bSJed Brown PetscFunctionBeginUser; 5699566063dSJacob Faibussowitsch PetscCall(ExactSolution(dm, time, x, yexact, ctx)); 570c4762a1bSJed Brown f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 5713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 572c4762a1bSJed Brown } 573c4762a1bSJed Brown 574d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 575d71ae5a4SJacob Faibussowitsch { 576c4762a1bSJed Brown PetscFunctionBeginUser; 5779566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 5789566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 5799566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 5809566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 5813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 582c4762a1bSJed Brown } 583c4762a1bSJed Brown 584d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupBC(DM dm, AppCtx *user) 585d71ae5a4SJacob Faibussowitsch { 586348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 58730602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 588c4762a1bSJed Brown PetscDS prob; 589c4762a1bSJed Brown DMLabel label; 590c4762a1bSJed Brown PetscBool check; 59130602db0SMatthew G. Knepley PetscInt dim, n = 3; 59230602db0SMatthew G. Knepley const char *prefix; 593c4762a1bSJed Brown 594c4762a1bSJed Brown PetscFunctionBeginUser; 5959566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 5969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *)bdt, &n, NULL)); 5979566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 598c4762a1bSJed Brown /* Set initial guesses and exact solutions */ 59930602db0SMatthew G. Knepley switch (dim) { 600c4762a1bSJed Brown case 2: 601c4762a1bSJed Brown user->initialGuess[0] = initialVelocity; 602c4762a1bSJed Brown switch (user->porosityDist) { 603d71ae5a4SJacob Faibussowitsch case ZERO: 604d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = zero_phi; 605d71ae5a4SJacob Faibussowitsch break; 606d71ae5a4SJacob Faibussowitsch case CONSTANT: 607d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = constant_phi; 608d71ae5a4SJacob Faibussowitsch break; 609d71ae5a4SJacob Faibussowitsch case GAUSSIAN: 610d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = gaussian_phi_2d; 611d71ae5a4SJacob Faibussowitsch break; 612d71ae5a4SJacob Faibussowitsch case DELTA: 613d71ae5a4SJacob Faibussowitsch user->initialGuess[1] = delta_phi_2d; 614d71ae5a4SJacob Faibussowitsch break; 615c4762a1bSJed Brown case TILTED: 616c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 617c4762a1bSJed Brown else user->initialGuess[1] = tilted_phi_coupled_2d; 618c4762a1bSJed Brown break; 619c4762a1bSJed Brown } 620348a1646SMatthew G. Knepley break; 621d71ae5a4SJacob Faibussowitsch default: 622d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 623348a1646SMatthew G. Knepley } 624348a1646SMatthew G. Knepley exactFuncs[0] = user->initialGuess[0]; 625348a1646SMatthew G. Knepley exactFuncs[1] = user->initialGuess[1]; 62630602db0SMatthew G. Knepley switch (dim) { 627348a1646SMatthew G. Knepley case 2: 628c4762a1bSJed Brown switch (user->velocityDist) { 629d71ae5a4SJacob Faibussowitsch case VEL_ZERO: 630d71ae5a4SJacob Faibussowitsch exactFuncs[0] = zero_u_2d; 631d71ae5a4SJacob Faibussowitsch break; 632d71ae5a4SJacob Faibussowitsch case VEL_CONSTANT: 633d71ae5a4SJacob Faibussowitsch exactFuncs[0] = constant_u_2d; 634d71ae5a4SJacob Faibussowitsch break; 635c4762a1bSJed Brown case VEL_HARMONIC: 63630602db0SMatthew G. Knepley switch (bdt[0]) { 637c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 63830602db0SMatthew G. Knepley switch (bdt[1]) { 639d71ae5a4SJacob Faibussowitsch case DM_BOUNDARY_PERIODIC: 640d71ae5a4SJacob Faibussowitsch exactFuncs[0] = doubly_periodic_u_2d; 641d71ae5a4SJacob Faibussowitsch break; 642d71ae5a4SJacob Faibussowitsch default: 643d71ae5a4SJacob Faibussowitsch exactFuncs[0] = periodic_u_2d; 644d71ae5a4SJacob Faibussowitsch break; 645c4762a1bSJed Brown } 646c4762a1bSJed Brown break; 647d71ae5a4SJacob Faibussowitsch default: 648d71ae5a4SJacob Faibussowitsch exactFuncs[0] = quadratic_u_2d; 649d71ae5a4SJacob Faibussowitsch break; 650c4762a1bSJed Brown } 651c4762a1bSJed Brown break; 652d71ae5a4SJacob Faibussowitsch case VEL_SHEAR: 653d71ae5a4SJacob Faibussowitsch exactFuncs[0] = shear_bc; 654d71ae5a4SJacob Faibussowitsch break; 655d71ae5a4SJacob Faibussowitsch default: 656d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 657c4762a1bSJed Brown } 658348a1646SMatthew G. Knepley break; 659d71ae5a4SJacob Faibussowitsch default: 660d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 661c4762a1bSJed Brown } 662c4762a1bSJed Brown { 663c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 664c4762a1bSJed Brown 6659566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, "", "-use_implicit", &isImplicit)); 666348a1646SMatthew G. Knepley if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 667c4762a1bSJed Brown } 6689566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-dmts_check", &check)); 669c4762a1bSJed Brown if (check) { 670348a1646SMatthew G. Knepley user->initialGuess[0] = exactFuncs[0]; 671348a1646SMatthew G. Knepley user->initialGuess[1] = exactFuncs[1]; 672c4762a1bSJed Brown } 673c4762a1bSJed Brown /* Set BC */ 6749566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 6759566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 6769566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user)); 6779566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user)); 678c4762a1bSJed Brown if (label) { 679c4762a1bSJed Brown const PetscInt id = 1; 680c4762a1bSJed Brown 6819566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, user, NULL)); 682c4762a1bSJed Brown } 6839566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &label)); 684c4762a1bSJed Brown if (label && user->useFV) { 685c4762a1bSJed Brown const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101}; 686c4762a1bSJed Brown 687dd39110bSPierre 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)); 688dd39110bSPierre 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)); 689c4762a1bSJed Brown } 6903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 691c4762a1bSJed Brown } 692c4762a1bSJed Brown 693d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 694d71ae5a4SJacob Faibussowitsch { 69530602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 696c4762a1bSJed Brown PetscDS prob; 69730602db0SMatthew G. Knepley PetscInt n = 3; 69830602db0SMatthew G. Knepley const char *prefix; 699c4762a1bSJed Brown 700c4762a1bSJed Brown PetscFunctionBeginUser; 7019566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 7029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *)bdt, &n, NULL)); 7039566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 704c4762a1bSJed Brown switch (user->velocityDist) { 705d71ae5a4SJacob Faibussowitsch case VEL_ZERO: 706d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u)); 707d71ae5a4SJacob Faibussowitsch break; 708c4762a1bSJed Brown case VEL_CONSTANT: 7099566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u)); 7109566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL)); 7119566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL)); 712c4762a1bSJed Brown break; 713c4762a1bSJed Brown case VEL_HARMONIC: 71430602db0SMatthew G. Knepley switch (bdt[0]) { 715c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 71630602db0SMatthew G. Knepley switch (bdt[1]) { 717d71ae5a4SJacob Faibussowitsch case DM_BOUNDARY_PERIODIC: 718d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u)); 719d71ae5a4SJacob Faibussowitsch break; 720d71ae5a4SJacob Faibussowitsch default: 721d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u)); 722d71ae5a4SJacob Faibussowitsch break; 723c4762a1bSJed Brown } 724c4762a1bSJed Brown break; 725d71ae5a4SJacob Faibussowitsch default: 726d71ae5a4SJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u)); 727d71ae5a4SJacob Faibussowitsch break; 728c4762a1bSJed Brown } 7299566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 730c4762a1bSJed Brown break; 731c4762a1bSJed Brown case VEL_SHEAR: 7329566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u)); 7339566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 734c4762a1bSJed Brown break; 735c4762a1bSJed Brown } 7369566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_advection, f1_advection)); 7379566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL)); 7389566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL)); 7399566063dSJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_advection)); 7409566063dSJacob Faibussowitsch else PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection)); 741c4762a1bSJed Brown 7429566063dSJacob Faibussowitsch PetscCall(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user)); 7433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 744c4762a1bSJed Brown } 745c4762a1bSJed Brown 746d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 747d71ae5a4SJacob Faibussowitsch { 748c4762a1bSJed Brown DM cdm = dm; 749c4762a1bSJed Brown PetscQuadrature q; 750c4762a1bSJed Brown PetscFE fe[2]; 751c4762a1bSJed Brown PetscFV fv; 752c4762a1bSJed Brown MPI_Comm comm; 75330602db0SMatthew G. Knepley PetscInt dim; 754c4762a1bSJed Brown 755c4762a1bSJed Brown PetscFunctionBeginUser; 756c4762a1bSJed Brown /* Create finite element */ 7579566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 7589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7599566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0])); 7609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 7619566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1])); 7629566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 7639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "porosity")); 764c4762a1bSJed Brown 7659566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)dm), &fv)); 7669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fv, "porosity")); 7679566063dSJacob Faibussowitsch PetscCall(PetscFVSetFromOptions(fv)); 7689566063dSJacob Faibussowitsch PetscCall(PetscFVSetNumComponents(fv, 1)); 7699566063dSJacob Faibussowitsch PetscCall(PetscFVSetSpatialDimension(fv, dim)); 7709566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe[0], &q)); 7719566063dSJacob Faibussowitsch PetscCall(PetscFVSetQuadrature(fv, q)); 772c4762a1bSJed Brown 7739566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 7749566063dSJacob Faibussowitsch if (user->useFV) PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fv)); 7759566063dSJacob Faibussowitsch else PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 7769566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 7779566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 778c4762a1bSJed Brown 779c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 780c4762a1bSJed Brown while (cdm) { 7819566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 7829566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 783c4762a1bSJed Brown /* Coordinates were never localized for coarse meshes */ 7849566063dSJacob Faibussowitsch if (cdm) PetscCall(DMLocalizeCoordinates(cdm)); 785c4762a1bSJed Brown } 7869566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 7879566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 7889566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fv)); 7893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 790c4762a1bSJed Brown } 791c4762a1bSJed Brown 792d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 793d71ae5a4SJacob Faibussowitsch { 794c4762a1bSJed Brown PetscFunctionBeginUser; 7959566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, user, dm)); 796c4762a1bSJed Brown /* Handle refinement, etc. */ 7979566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 798c4762a1bSJed Brown /* Construct ghost cells */ 799c4762a1bSJed Brown if (user->useFV) { 800c4762a1bSJed Brown DM gdm; 801c4762a1bSJed Brown 8029566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm)); 8039566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 804c4762a1bSJed Brown *dm = gdm; 805c4762a1bSJed Brown } 806c4762a1bSJed Brown /* Localize coordinates */ 8079566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 8089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)(*dm), "Mesh")); 8099566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 810c4762a1bSJed Brown /* Setup problem */ 8119566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*dm, user)); 812c4762a1bSJed Brown /* Setup BC */ 8139566063dSJacob Faibussowitsch PetscCall(SetupBC(*dm, user)); 8143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 815c4762a1bSJed Brown } 816c4762a1bSJed Brown 817d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx) 818d71ae5a4SJacob Faibussowitsch { 819c4762a1bSJed Brown PetscDS prob; 820c4762a1bSJed Brown DM dmCell; 821c4762a1bSJed Brown Vec cellgeom; 822c4762a1bSJed Brown const PetscScalar *cgeom; 823c4762a1bSJed Brown PetscScalar *x; 824c4762a1bSJed Brown PetscInt dim, Nf, cStart, cEnd, c; 825c4762a1bSJed Brown 826c4762a1bSJed Brown PetscFunctionBeginUser; 8279566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 8289566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8299566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 8309566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8319566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 8329566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 8339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 8349566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 835c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 836c4762a1bSJed Brown PetscFVCellGeom *cg; 837c4762a1bSJed Brown PetscScalar *xc; 838c4762a1bSJed Brown 8399566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 8409566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc)); 8419566063dSJacob Faibussowitsch if (xc) PetscCall((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx)); 842c4762a1bSJed Brown } 8439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 8449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 8453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 846c4762a1bSJed Brown } 847c4762a1bSJed Brown 848d71ae5a4SJacob Faibussowitsch static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 849d71ae5a4SJacob Faibussowitsch { 850c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 851c4762a1bSJed Brown char *ftable = NULL; 852c4762a1bSJed Brown DM dm; 853c4762a1bSJed Brown PetscSection s; 854c4762a1bSJed Brown Vec cellgeom; 855c4762a1bSJed Brown const PetscScalar *x; 856c4762a1bSJed Brown PetscScalar *a; 857c4762a1bSJed Brown PetscReal *xnorms; 858c4762a1bSJed Brown PetscInt pStart, pEnd, p, Nf, f; 859c4762a1bSJed Brown 860c4762a1bSJed Brown PetscFunctionBeginUser; 8619566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X, (PetscObject)ts, "-view_solution")); 8629566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 8639566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8649566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 8659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &Nf)); 8669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 8679566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nf * 2, &xnorms)); 8689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 869c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 870c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 871c4762a1bSJed Brown PetscInt dof, cdof, d; 872c4762a1bSJed Brown 8739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 8749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 8759566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, x, &a)); 876c4762a1bSJed Brown /* TODO Use constrained indices here */ 877c4762a1bSJed Brown for (d = 0; d < dof - cdof; ++d) xnorms[f * 2 + 0] = PetscMax(xnorms[f * 2 + 0], PetscAbsScalar(a[d])); 878c4762a1bSJed Brown for (d = 0; d < dof - cdof; ++d) xnorms[f * 2 + 1] += PetscAbsScalar(a[d]); 879c4762a1bSJed Brown } 880c4762a1bSJed Brown } 8819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 882c4762a1bSJed Brown if (stepnum >= 0) { /* No summary for final time */ 883c4762a1bSJed Brown DM dmCell, *fdm; 884c4762a1bSJed Brown Vec *fv; 885c4762a1bSJed Brown const PetscScalar *cgeom; 886c4762a1bSJed Brown PetscScalar **fx; 887c4762a1bSJed Brown PetscReal *fmin, *fmax, *fint, *ftmp, t; 888c4762a1bSJed Brown PetscInt cStart, cEnd, c, fcount, f, num; 889c4762a1bSJed Brown 890c4762a1bSJed Brown size_t ftableused, ftablealloc; 891c4762a1bSJed Brown 892c4762a1bSJed Brown /* Functionals have indices after registering, this is an upper bound */ 893c4762a1bSJed Brown fcount = user->numMonitorFuncs; 8949566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fint, fcount, &ftmp)); 8959566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fcount, &fdm, fcount, &fv, fcount, &fx)); 896c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 897c4762a1bSJed Brown PetscSection fs; 898c4762a1bSJed Brown const char *name = user->monitorFuncs[f]->name; 899c4762a1bSJed Brown 900c4762a1bSJed Brown fmin[f] = PETSC_MAX_REAL; 901c4762a1bSJed Brown fmax[f] = PETSC_MIN_REAL; 902c4762a1bSJed Brown fint[f] = 0; 903c4762a1bSJed Brown /* Make monitor vecs */ 9049566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &fdm[f])); 9059566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &num, &t)); 9069566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(fdm[f], num, t)); 9079566063dSJacob Faibussowitsch PetscCall(PetscSectionClone(s, &fs)); 9089566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 0, NULL)); 9099566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 1, name)); 9109566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(fdm[f], fs)); 9119566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&fs)); 9129566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm[f], &fv[f])); 9139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fv[f], name)); 9149566063dSJacob Faibussowitsch PetscCall(VecGetArray(fv[f], &fx[f])); 915c4762a1bSJed Brown } 9169566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 9179566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 9189566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 9199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 920c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 921c4762a1bSJed Brown PetscFVCellGeom *cg; 922c4762a1bSJed Brown PetscScalar *cx; 923c4762a1bSJed Brown 9249566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 9259566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx)); 926c4762a1bSJed Brown if (!cx) continue; /* not a global cell */ 927c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 928c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 929c4762a1bSJed Brown PetscScalar *fxc; 930c4762a1bSJed Brown 9319566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc)); 932c4762a1bSJed Brown /* I need to make it easier to get interpolated values here */ 9339566063dSJacob Faibussowitsch PetscCall((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx)); 934c4762a1bSJed Brown fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 935c4762a1bSJed Brown } 936c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 937c4762a1bSJed Brown fmin[f] = PetscMin(fmin[f], ftmp[f]); 938c4762a1bSJed Brown fmax[f] = PetscMax(fmax[f], ftmp[f]); 939c4762a1bSJed Brown fint[f] += cg->volume * ftmp[f]; 940c4762a1bSJed Brown } 941c4762a1bSJed Brown } 9429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 9439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 944712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 945712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 946712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 947c4762a1bSJed Brown /* Output functional data */ 948c4762a1bSJed Brown ftablealloc = fcount * 100; 949c4762a1bSJed Brown ftableused = 0; 9509566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ftablealloc, &ftable)); 951c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 952c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 953c4762a1bSJed Brown PetscInt id = func->offset; 954c4762a1bSJed Brown char newline[] = "\n"; 955c4762a1bSJed Brown char buffer[256], *p, *prefix; 956c4762a1bSJed Brown size_t countused, len; 957c4762a1bSJed Brown 958c4762a1bSJed Brown /* Create string with functional outputs */ 959c4762a1bSJed Brown if (f % 3) { 9609566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, " ", 2)); 961c4762a1bSJed Brown p = buffer + 2; 962c4762a1bSJed Brown } else if (f) { 9639566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, newline, sizeof(newline) - 1)); 964c4762a1bSJed Brown p = buffer + sizeof(newline) - 1; 965c4762a1bSJed Brown } else { 966c4762a1bSJed Brown p = buffer; 967c4762a1bSJed Brown } 9689566063dSJacob 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])); 969c4762a1bSJed Brown countused += p - buffer; 970c4762a1bSJed Brown /* reallocate */ 971c4762a1bSJed Brown if (countused > ftablealloc - ftableused - 1) { 972c4762a1bSJed Brown char *ftablenew; 973c4762a1bSJed Brown 974c4762a1bSJed Brown ftablealloc = 2 * ftablealloc + countused; 9759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ftablealloc, &ftablenew)); 9769566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); 9779566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 978c4762a1bSJed Brown ftable = ftablenew; 979c4762a1bSJed Brown } 9809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused)); 981c4762a1bSJed Brown ftableused += countused; 982c4762a1bSJed Brown ftable[ftableused] = 0; 983c4762a1bSJed Brown /* Output vecs */ 9849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fv[f], &fx[f])); 9859566063dSJacob Faibussowitsch PetscCall(PetscStrlen(func->name, &len)); 9869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len + 2, &prefix)); 987c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(prefix, func->name, len + 2)); 988c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(prefix, "_", len + 2)); 9899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix)); 9909566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fv[f], NULL, "-vec_view")); 9919566063dSJacob Faibussowitsch PetscCall(PetscFree(prefix)); 9929566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm[f], &fv[f])); 9939566063dSJacob Faibussowitsch PetscCall(DMDestroy(&fdm[f])); 994c4762a1bSJed Brown } 9959566063dSJacob Faibussowitsch PetscCall(PetscFree4(fmin, fmax, fint, ftmp)); 9969566063dSJacob Faibussowitsch PetscCall(PetscFree3(fdm, fv, fx)); 99763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT " time %8.4g |x| (", stepnum, (double)time)); 998c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 9999566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ", ")); 10009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%8.4g", (double)xnorms[f * 2 + 0])); 1001c4762a1bSJed Brown } 10029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ") |x|_1 (")); 1003c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10049566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ", ")); 10059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%8.4g", (double)xnorms[f * 2 + 1])); 1006c4762a1bSJed Brown } 10079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ") %s\n", ftable ? ftable : "")); 10089566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 1009c4762a1bSJed Brown } 10109566063dSJacob Faibussowitsch PetscCall(PetscFree(xnorms)); 10113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1012c4762a1bSJed Brown } 1013c4762a1bSJed Brown 1014d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 1015d71ae5a4SJacob Faibussowitsch { 1016c4762a1bSJed Brown MPI_Comm comm; 1017c4762a1bSJed Brown TS ts; 1018c4762a1bSJed Brown DM dm; 1019c4762a1bSJed Brown Vec u; 1020c4762a1bSJed Brown AppCtx user; 1021c4762a1bSJed Brown PetscReal t0, t = 0.0; 1022*4d86920dSPierre Jolivet void *ctxs[2] = {&t, &t}; 1023c4762a1bSJed Brown 1024327415f7SBarry Smith PetscFunctionBeginUser; 10259566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1026c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1027c4762a1bSJed Brown user.functionalRegistry = NULL; 1028c4762a1bSJed Brown globalUser = &user; 10299566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 10309566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 10319566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 10329566063dSJacob Faibussowitsch PetscCall(CreateDM(comm, &user, &dm)); 10339566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 10349566063dSJacob Faibussowitsch PetscCall(ProcessMonitorOptions(comm, &user)); 1035c4762a1bSJed Brown 10369566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 10379566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 1038c4762a1bSJed Brown if (user.useFV) { 1039c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 1040c4762a1bSJed Brown 10419566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, "", "-use_implicit", &isImplicit)); 1042c4762a1bSJed Brown if (isImplicit) { 10439566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10449566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1045c4762a1bSJed Brown } 10469566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10479566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user)); 1048c4762a1bSJed Brown } else { 10499566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10509566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10519566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1052c4762a1bSJed Brown } 10539566063dSJacob Faibussowitsch if (user.useFV) PetscCall(TSMonitorSet(ts, MonitorFunctionals, &user, NULL)); 10549566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1)); 10559566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2.0)); 10569566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01)); 10579566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 10589566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1059c4762a1bSJed Brown 10609566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u)); 10619566063dSJacob Faibussowitsch if (user.useFV) PetscCall(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1])); 10629566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-init_vec_view")); 10639566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1064c4762a1bSJed Brown t0 = t; 10659566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 10669566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 10679566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 10689566063dSJacob Faibussowitsch if (t > t0) PetscCall(DMTSCheckFromOptions(ts, u)); 10699566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 1070c4762a1bSJed Brown { 1071c4762a1bSJed Brown PetscReal ftime; 1072c4762a1bSJed Brown PetscInt nsteps; 1073c4762a1bSJed Brown TSConvergedReason reason; 1074c4762a1bSJed Brown 10759566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 10769566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &nsteps)); 10779566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 107863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps)); 1079c4762a1bSJed Brown } 1080c4762a1bSJed Brown 10819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 10829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 10839566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 10849566063dSJacob Faibussowitsch PetscCall(PetscFree(user.monitorFuncs)); 10859566063dSJacob Faibussowitsch PetscCall(FunctionalDestroy(&user.functionalRegistry)); 10869566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1087b122ec5aSJacob Faibussowitsch return 0; 1088c4762a1bSJed Brown } 1089c4762a1bSJed Brown 1090c4762a1bSJed Brown /*TEST 1091c4762a1bSJed Brown 109230602db0SMatthew G. Knepley testset: 109330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 109430602db0SMatthew G. Knepley 1095c4762a1bSJed Brown # 2D harmonic velocity, no porosity 1096c4762a1bSJed Brown test: 1097c4762a1bSJed Brown suffix: p1p1 1098c4762a1bSJed Brown requires: !complex !single 1099b4f8a55aSStefano 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 1100c4762a1bSJed Brown test: 1101c4762a1bSJed Brown suffix: p1p1_xper 1102c4762a1bSJed Brown requires: !complex !single 110330602db0SMatthew 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 1104c4762a1bSJed Brown test: 1105c4762a1bSJed Brown suffix: p1p1_xper_ref 1106c4762a1bSJed Brown requires: !complex !single 110730602db0SMatthew 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 1108c4762a1bSJed Brown test: 1109c4762a1bSJed Brown suffix: p1p1_xyper 1110c4762a1bSJed Brown requires: !complex !single 111130602db0SMatthew 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 1112c4762a1bSJed Brown test: 1113c4762a1bSJed Brown suffix: p1p1_xyper_ref 1114c4762a1bSJed Brown requires: !complex !single 111530602db0SMatthew 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 1116c4762a1bSJed Brown test: 1117c4762a1bSJed Brown suffix: p2p1 1118c4762a1bSJed Brown requires: !complex !single 111930602db0SMatthew 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 1120c4762a1bSJed Brown test: 1121c4762a1bSJed Brown suffix: p2p1_xyper 1122c4762a1bSJed Brown requires: !complex !single 112330602db0SMatthew 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 112430602db0SMatthew G. Knepley 112530602db0SMatthew G. Knepley test: 112630602db0SMatthew G. Knepley suffix: adv_1 112730602db0SMatthew G. Knepley requires: !complex !single 112830602db0SMatthew 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 112930602db0SMatthew G. Knepley 113030602db0SMatthew G. Knepley test: 113130602db0SMatthew G. Knepley suffix: adv_2 113230602db0SMatthew G. Knepley requires: !complex 113330602db0SMatthew G. Knepley TODO: broken memory corruption 113430602db0SMatthew 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 113530602db0SMatthew G. Knepley 113630602db0SMatthew G. Knepley test: 113730602db0SMatthew G. Knepley suffix: adv_3 113830602db0SMatthew G. Knepley requires: !complex 113930602db0SMatthew G. Knepley TODO: broken memory corruption 114030602db0SMatthew 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 114130602db0SMatthew G. Knepley 114230602db0SMatthew G. Knepley test: 114330602db0SMatthew G. Knepley suffix: adv_3_ex 114430602db0SMatthew G. Knepley requires: !complex 114530602db0SMatthew 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 114630602db0SMatthew G. Knepley 114730602db0SMatthew G. Knepley test: 114830602db0SMatthew G. Knepley suffix: adv_4 114930602db0SMatthew G. Knepley requires: !complex 115030602db0SMatthew G. Knepley TODO: broken memory corruption 115130602db0SMatthew 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 115230602db0SMatthew G. Knepley 115330602db0SMatthew G. Knepley # 2D Advection, box, delta 115430602db0SMatthew G. Knepley test: 115530602db0SMatthew G. Knepley suffix: adv_delta_yper_0 115630602db0SMatthew G. Knepley requires: !complex 115730602db0SMatthew G. Knepley TODO: broken 115830602db0SMatthew 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 115930602db0SMatthew G. Knepley 116030602db0SMatthew G. Knepley test: 116130602db0SMatthew G. Knepley suffix: adv_delta_yper_1 116230602db0SMatthew G. Knepley requires: !complex 116330602db0SMatthew G. Knepley TODO: broken 116430602db0SMatthew 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 116530602db0SMatthew G. Knepley 116630602db0SMatthew G. Knepley test: 116730602db0SMatthew G. Knepley suffix: adv_delta_yper_2 116830602db0SMatthew G. Knepley requires: !complex 116930602db0SMatthew G. Knepley TODO: broken 117030602db0SMatthew 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 117130602db0SMatthew G. Knepley 117230602db0SMatthew G. Knepley test: 117330602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_0 117430602db0SMatthew G. Knepley requires: !complex 117530602db0SMatthew G. Knepley TODO: broken 117630602db0SMatthew 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 117730602db0SMatthew G. Knepley 117830602db0SMatthew G. Knepley test: 117930602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_1 118030602db0SMatthew G. Knepley requires: !complex 118130602db0SMatthew G. Knepley TODO: broken 118230602db0SMatthew 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 118330602db0SMatthew G. Knepley 118430602db0SMatthew G. Knepley test: 118530602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_2 118630602db0SMatthew G. Knepley requires: !complex 118730602db0SMatthew G. Knepley TODO: broken 118830602db0SMatthew 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 118930602db0SMatthew G. Knepley 119030602db0SMatthew G. Knepley test: 119130602db0SMatthew G. Knepley suffix: adv_delta_yper_im_0 119230602db0SMatthew G. Knepley requires: !complex 119330602db0SMatthew G. Knepley TODO: broken 119430602db0SMatthew 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 119530602db0SMatthew G. Knepley 119630602db0SMatthew G. Knepley test: 119730602db0SMatthew G. Knepley suffix: adv_delta_yper_im_1 119830602db0SMatthew G. Knepley requires: !complex 119930602db0SMatthew G. Knepley TODO: broken 120030602db0SMatthew 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 120130602db0SMatthew G. Knepley 120230602db0SMatthew G. Knepley test: 120330602db0SMatthew G. Knepley suffix: adv_delta_yper_im_2 120430602db0SMatthew G. Knepley requires: !complex 120530602db0SMatthew G. Knepley TODO: broken 120630602db0SMatthew 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 120730602db0SMatthew G. Knepley 120830602db0SMatthew G. Knepley test: 120930602db0SMatthew G. Knepley suffix: adv_delta_yper_im_3 121030602db0SMatthew G. Knepley requires: !complex 121130602db0SMatthew G. Knepley TODO: broken 121230602db0SMatthew 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 121330602db0SMatthew G. Knepley 121430602db0SMatthew G. Knepley # I believe the nullspace is sin(pi y) 121530602db0SMatthew G. Knepley test: 121630602db0SMatthew G. Knepley suffix: adv_delta_yper_im_4 121730602db0SMatthew G. Knepley requires: !complex 121830602db0SMatthew G. Knepley TODO: broken 121930602db0SMatthew 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 122030602db0SMatthew G. Knepley 122130602db0SMatthew G. Knepley test: 122230602db0SMatthew G. Knepley suffix: adv_delta_yper_im_5 122330602db0SMatthew G. Knepley requires: !complex 122430602db0SMatthew G. Knepley TODO: broken 122530602db0SMatthew 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 122630602db0SMatthew G. Knepley 122730602db0SMatthew G. Knepley test: 122830602db0SMatthew G. Knepley suffix: adv_delta_yper_im_6 122930602db0SMatthew G. Knepley requires: !complex 123030602db0SMatthew G. Knepley TODO: broken 123130602db0SMatthew 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 123230602db0SMatthew G. Knepley # 2D Advection, magma benchmark 1 123330602db0SMatthew G. Knepley 123430602db0SMatthew G. Knepley test: 123530602db0SMatthew G. Knepley suffix: adv_delta_shear_im_0 123630602db0SMatthew G. Knepley requires: !complex 123730602db0SMatthew G. Knepley TODO: broken 123830602db0SMatthew 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 123930602db0SMatthew G. Knepley # 2D Advection, box, gaussian 124030602db0SMatthew G. Knepley 124130602db0SMatthew G. Knepley test: 124230602db0SMatthew G. Knepley suffix: adv_gauss 124330602db0SMatthew G. Knepley requires: !complex 124430602db0SMatthew G. Knepley TODO: broken 124530602db0SMatthew 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 124630602db0SMatthew G. Knepley 124730602db0SMatthew G. Knepley test: 124830602db0SMatthew G. Knepley suffix: adv_gauss_im 124930602db0SMatthew G. Knepley requires: !complex 125030602db0SMatthew G. Knepley TODO: broken 125130602db0SMatthew 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 125230602db0SMatthew G. Knepley 125330602db0SMatthew G. Knepley test: 125430602db0SMatthew G. Knepley suffix: adv_gauss_im_1 125530602db0SMatthew G. Knepley requires: !complex 125630602db0SMatthew G. Knepley TODO: broken 125730602db0SMatthew 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 125830602db0SMatthew G. Knepley 125930602db0SMatthew G. Knepley test: 126030602db0SMatthew G. Knepley suffix: adv_gauss_im_2 126130602db0SMatthew G. Knepley requires: !complex 126230602db0SMatthew G. Knepley TODO: broken 126330602db0SMatthew 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 126430602db0SMatthew G. Knepley 126530602db0SMatthew G. Knepley test: 126630602db0SMatthew G. Knepley suffix: adv_gauss_corner 126730602db0SMatthew G. Knepley requires: !complex 126830602db0SMatthew G. Knepley TODO: broken 126930602db0SMatthew 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 127030602db0SMatthew G. Knepley 127130602db0SMatthew G. Knepley # 2D Advection+Harmonic 12- 127230602db0SMatthew G. Knepley test: 127330602db0SMatthew G. Knepley suffix: adv_harm_0 127430602db0SMatthew G. Knepley requires: !complex 127530602db0SMatthew G. Knepley TODO: broken memory corruption 127630602db0SMatthew 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 127730602db0SMatthew G. Knepley 1278c4762a1bSJed Brown # Must check that FV BCs propagate to coarse meshes 1279c4762a1bSJed Brown # Must check that FV BC ids propagate to coarse meshes 1280c4762a1bSJed Brown # Must check that FE+FV BCs work at the same time 1281c4762a1bSJed Brown # 2D Advection, matching wind in ex11 8-11 1282c4762a1bSJed Brown # NOTE implicit solves are limited by accuracy of FD Jacobian 1283c4762a1bSJed Brown test: 1284c4762a1bSJed Brown suffix: adv_0 1285c4762a1bSJed Brown requires: !complex !single exodusii 128630602db0SMatthew 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 1287c4762a1bSJed Brown 1288c4762a1bSJed Brown test: 1289c4762a1bSJed Brown suffix: adv_0_im 1290c4762a1bSJed Brown requires: !complex exodusii 1291c4762a1bSJed Brown TODO: broken memory corruption 129230602db0SMatthew 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 1293c4762a1bSJed Brown 1294c4762a1bSJed Brown test: 1295c4762a1bSJed Brown suffix: adv_0_im_2 1296c4762a1bSJed Brown requires: !complex exodusii 1297c4762a1bSJed Brown TODO: broken 129830602db0SMatthew 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 1299c4762a1bSJed Brown 1300c4762a1bSJed Brown test: 1301c4762a1bSJed Brown suffix: adv_0_im_3 1302c4762a1bSJed Brown requires: !complex exodusii 1303c4762a1bSJed Brown TODO: broken 130430602db0SMatthew 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 1305c4762a1bSJed Brown 1306c4762a1bSJed Brown test: 1307c4762a1bSJed Brown suffix: adv_0_im_4 1308c4762a1bSJed Brown requires: !complex exodusii 1309c4762a1bSJed Brown TODO: broken 131030602db0SMatthew 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 1311c4762a1bSJed Brown # 2D Advection, misc 1312c4762a1bSJed Brown 1313c4762a1bSJed Brown TEST*/ 1314