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 17c4762a1bSJed Brown #define ALEN(a) (sizeof(a)/sizeof((a)[0])) 18c4762a1bSJed Brown 19c4762a1bSJed Brown typedef enum {VEL_ZERO, VEL_CONSTANT, VEL_HARMONIC, VEL_SHEAR} VelocityDistribution; 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef enum {ZERO, CONSTANT, GAUSSIAN, TILTED, DELTA} PorosityDistribution; 22c4762a1bSJed Brown 23c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* 26c4762a1bSJed Brown FunctionalFunc - Calculates the value of a functional of the solution at a point 27c4762a1bSJed Brown 28c4762a1bSJed Brown Input Parameters: 29c4762a1bSJed Brown + dm - The DM 30c4762a1bSJed Brown . time - The TS time 31c4762a1bSJed Brown . x - The coordinates of the evaluation point 32c4762a1bSJed Brown . u - The field values at point x 33c4762a1bSJed Brown - ctx - A user context, or NULL 34c4762a1bSJed Brown 35c4762a1bSJed Brown Output Parameter: 36c4762a1bSJed Brown . f - The value of the functional at point x 37c4762a1bSJed Brown 38c4762a1bSJed Brown */ 39c4762a1bSJed Brown typedef PetscErrorCode (*FunctionalFunc)(DM, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); 40c4762a1bSJed Brown 41c4762a1bSJed Brown typedef struct _n_Functional *Functional; 42c4762a1bSJed Brown struct _n_Functional { 43c4762a1bSJed Brown char *name; 44c4762a1bSJed Brown FunctionalFunc func; 45c4762a1bSJed Brown void *ctx; 46c4762a1bSJed Brown PetscInt offset; 47c4762a1bSJed Brown Functional next; 48c4762a1bSJed Brown }; 49c4762a1bSJed Brown 50c4762a1bSJed Brown typedef struct { 51c4762a1bSJed Brown /* Problem definition */ 52c4762a1bSJed Brown PetscBool useFV; /* Use a finite volume scheme for advection */ 53c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 54c4762a1bSJed Brown VelocityDistribution velocityDist; 55c4762a1bSJed Brown PorosityDistribution porosityDist; 56c4762a1bSJed Brown PetscReal inflowState; 57c4762a1bSJed Brown PetscReal source[3]; 58c4762a1bSJed Brown /* Monitoring */ 59c4762a1bSJed Brown PetscInt numMonitorFuncs, maxMonitorFunc; 60c4762a1bSJed Brown Functional *monitorFuncs; 61c4762a1bSJed Brown PetscInt errorFunctional; 62c4762a1bSJed Brown Functional functionalRegistry; 63c4762a1bSJed Brown } AppCtx; 64c4762a1bSJed Brown 65c4762a1bSJed Brown static AppCtx *globalUser; 66c4762a1bSJed Brown 67c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown const char *velocityDist[4] = {"zero", "constant", "harmonic", "shear"}; 70c4762a1bSJed Brown const char *porosityDist[5] = {"zero", "constant", "gaussian", "tilted", "delta"}; 7130602db0SMatthew G. Knepley PetscInt vd, pd, d; 72c4762a1bSJed Brown PetscBool flg; 73c4762a1bSJed Brown PetscErrorCode ierr; 74c4762a1bSJed Brown 75c4762a1bSJed Brown PetscFunctionBeginUser; 76c4762a1bSJed Brown options->useFV = PETSC_FALSE; 77c4762a1bSJed Brown options->velocityDist = VEL_HARMONIC; 78c4762a1bSJed Brown options->porosityDist = ZERO; 79c4762a1bSJed Brown options->inflowState = -2.0; 80c4762a1bSJed Brown options->numMonitorFuncs = 0; 81c4762a1bSJed Brown options->source[0] = 0.5; 82c4762a1bSJed Brown options->source[1] = 0.5; 83c4762a1bSJed Brown options->source[2] = 0.5; 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX");CHKERRQ(ierr); 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL)); 87c4762a1bSJed Brown vd = options->velocityDist; 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-velocity_dist","Velocity distribution type","ex18.c",velocityDist,4,velocityDist[options->velocityDist],&vd,NULL)); 89c4762a1bSJed Brown options->velocityDist = (VelocityDistribution) vd; 90c4762a1bSJed Brown pd = options->porosityDist; 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-porosity_dist","Initial porosity distribution type","ex18.c",porosityDist,5,porosityDist[options->porosityDist],&pd,NULL)); 92c4762a1bSJed Brown options->porosityDist = (PorosityDistribution) pd; 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL)); 9430602db0SMatthew G. Knepley d = 2; 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg)); 963c633725SBarry Smith PetscCheck(!flg || d == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %d", d); 97c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscFunctionReturn(0); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown Functional func; 105c4762a1bSJed Brown char *names[256]; 106c4762a1bSJed Brown PetscInt f; 107c4762a1bSJed Brown PetscErrorCode ierr; 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscFunctionBeginUser; 110c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX");CHKERRQ(ierr); 111c4762a1bSJed Brown options->numMonitorFuncs = ALEN(names); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs)); 114c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 115c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 116c4762a1bSJed Brown PetscBool match; 117c4762a1bSJed Brown 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcasecmp(names[f], func->name, &match)); 119c4762a1bSJed Brown if (match) break; 120c4762a1bSJed Brown } 1213c633725SBarry Smith PetscCheck(func,comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 122c4762a1bSJed Brown options->monitorFuncs[f] = func; 123c4762a1bSJed Brown /* Jed inserts a de-duplication of functionals here */ 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(names[f])); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 127c4762a1bSJed Brown options->maxMonitorFunc = -1; 128c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 129c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 130c4762a1bSJed Brown Functional call = options->monitorFuncs[f]; 131c4762a1bSJed Brown 132c4762a1bSJed Brown if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 135c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 136c4762a1bSJed Brown PetscFunctionReturn(0); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) 140c4762a1bSJed Brown { 141c4762a1bSJed Brown Functional *ptr, f; 142c4762a1bSJed Brown PetscInt lastoffset = -1; 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscFunctionBeginUser; 145c4762a1bSJed Brown for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&f)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(name, &f->name)); 148c4762a1bSJed Brown f->offset = lastoffset + 1; 149c4762a1bSJed Brown f->func = func; 150c4762a1bSJed Brown f->ctx = ctx; 151c4762a1bSJed Brown f->next = NULL; 152c4762a1bSJed Brown *ptr = f; 153c4762a1bSJed Brown *offset = f->offset; 154c4762a1bSJed Brown PetscFunctionReturn(0); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 157c4762a1bSJed Brown static PetscErrorCode FunctionalDestroy(Functional *link) 158c4762a1bSJed Brown { 159c4762a1bSJed Brown Functional next, l; 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscFunctionBeginUser; 162c4762a1bSJed Brown if (!link) PetscFunctionReturn(0); 163c4762a1bSJed Brown l = *link; 164c4762a1bSJed Brown *link = NULL; 165c4762a1bSJed Brown for (; l; l=next) { 166c4762a1bSJed Brown next = l->next; 1675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(l->name)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(l)); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown PetscFunctionReturn(0); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 174c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 175c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 176c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown PetscInt comp; 179c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 183c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 184c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 185c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 186c4762a1bSJed Brown { 187c4762a1bSJed Brown PetscScalar wind[3] = {0.0, 0.0, 0.0}; 188c4762a1bSJed Brown PetscInt comp; 189c4762a1bSJed Brown 190c4762a1bSJed Brown 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 194c4762a1bSJed Brown static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 195c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 196c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 197c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 198c4762a1bSJed Brown { 199c4762a1bSJed Brown PetscInt comp; 200c4762a1bSJed Brown for (comp = 0; comp < dim*dim; ++comp) f1[comp] = 0.0; 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203c4762a1bSJed Brown static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 204c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 205c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 206c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 207c4762a1bSJed Brown { 208c4762a1bSJed Brown PetscInt d; 209c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d*dim+d] = 1.0; 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 213c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 214c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 215c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 216c4762a1bSJed Brown { 217c4762a1bSJed Brown g0[0] = 1.0; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 221c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 222c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 223c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 224c4762a1bSJed Brown { 225c4762a1bSJed Brown PetscInt comp; 226c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 230c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 231c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 232c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 233c4762a1bSJed Brown { 234c4762a1bSJed Brown PetscInt comp, d; 235c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) { 236c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 237c4762a1bSJed Brown f1[comp*dim+d] = u_x[comp*dim+d]; 238c4762a1bSJed Brown } 239c4762a1bSJed Brown } 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 243c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 244c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 245c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 246c4762a1bSJed Brown { 247c4762a1bSJed Brown f0[0] = -PetscSinReal(2.0*PETSC_PI*x[0]); 248c4762a1bSJed Brown f0[1] = 2.0*PETSC_PI*x[1]*PetscCosReal(2.0*PETSC_PI*x[0]); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 252c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 253c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 254c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 255c4762a1bSJed Brown { 256c4762a1bSJed Brown f0[0] = -2.0*PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1]); 257c4762a1bSJed Brown f0[1] = 2.0*PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0]); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 261c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 262c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 263c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 264c4762a1bSJed Brown { 265c4762a1bSJed Brown const PetscInt Ncomp = dim; 266c4762a1bSJed Brown PetscInt compI, d; 267c4762a1bSJed Brown 268c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 269c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 270c4762a1bSJed Brown g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0; 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 276c4762a1bSJed Brown static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 277c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 278c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 279c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 280c4762a1bSJed Brown { 281c4762a1bSJed Brown PetscInt d; 282c4762a1bSJed Brown f0[0] = u_t[dim]; 283c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += u[dim]*u_x[d*dim+d] + u_x[dim*dim+d]*u[d]; 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 287c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 288c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 289c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 290c4762a1bSJed Brown { 291c4762a1bSJed Brown PetscInt d; 292c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[0] = 0.0; 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 296c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 297c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 298c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 299c4762a1bSJed Brown { 300c4762a1bSJed Brown PetscInt d; 301c4762a1bSJed Brown g0[0] = u_tShift; 302c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[d*dim+d]; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305c4762a1bSJed Brown void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 306c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 307c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 308c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 309c4762a1bSJed Brown { 310c4762a1bSJed Brown PetscInt d; 311c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = u[d]; 312c4762a1bSJed Brown } 313c4762a1bSJed Brown 314c4762a1bSJed Brown void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 315c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 316c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 317c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 318c4762a1bSJed Brown { 319c4762a1bSJed Brown PetscInt d; 320c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[dim*dim+d]; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 323c4762a1bSJed Brown void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 324c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 325c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 326c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 327c4762a1bSJed Brown { 328c4762a1bSJed Brown PetscInt d; 329c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = u[dim]; 330c4762a1bSJed Brown } 331c4762a1bSJed Brown 332c4762a1bSJed Brown 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) 333c4762a1bSJed Brown { 334c4762a1bSJed Brown PetscReal wind[3] = {0.0, 1.0, 0.0}; 335c4762a1bSJed Brown PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim,3), wind, n); 336c4762a1bSJed Brown 337c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340c4762a1bSJed Brown 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) 341c4762a1bSJed Brown { 342c4762a1bSJed Brown PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 343c4762a1bSJed Brown 344c4762a1bSJed Brown #if 1 345c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 346c4762a1bSJed Brown #else 347c4762a1bSJed 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]); */ 348c4762a1bSJed Brown /* Smear it out */ 349c4762a1bSJed Brown flux[0] = 0.5*((uL[dim] + uR[dim]) + (uL[dim] - uR[dim])*tanh(1.0e5*wn)) * wn; 350c4762a1bSJed Brown #endif 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 353c4762a1bSJed Brown static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 354c4762a1bSJed Brown { 355c4762a1bSJed Brown u[0] = 0.0; 356c4762a1bSJed Brown u[1] = 0.0; 357c4762a1bSJed Brown return 0; 358c4762a1bSJed Brown } 359c4762a1bSJed Brown 360c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 361c4762a1bSJed Brown { 362c4762a1bSJed Brown u[0] = 0.0; 363c4762a1bSJed Brown u[1] = 1.0; 364c4762a1bSJed Brown return 0; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* Coordinates of the point which was at x at t = 0 */ 368c4762a1bSJed Brown static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 369c4762a1bSJed Brown { 370c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 371c4762a1bSJed Brown u[0] = x[0]; 372c4762a1bSJed Brown u[1] = x[1] + t; 37330602db0SMatthew G. Knepley #if 0 37430602db0SMatthew G. Knepley PetscErrorCode ierr; 3755f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u)); 37630602db0SMatthew G. Knepley #else 37730602db0SMatthew G. Knepley u[1] = u[1] - (int) PetscRealPart(u[1]); 37830602db0SMatthew G. Knepley #endif 379c4762a1bSJed Brown return 0; 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 382c4762a1bSJed Brown /* 383c4762a1bSJed Brown In 2D we use the exact solution: 384c4762a1bSJed Brown 385c4762a1bSJed Brown u = x^2 + y^2 386c4762a1bSJed Brown v = 2 x^2 - 2xy 387c4762a1bSJed Brown phi = h(x + y + (u + v) t) 388c4762a1bSJed Brown f_x = f_y = 4 389c4762a1bSJed Brown 390c4762a1bSJed Brown so that 391c4762a1bSJed Brown 392c4762a1bSJed Brown -\Delta u + f = <-4, -4> + <4, 4> = 0 393c4762a1bSJed Brown {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 394c4762a1bSJed Brown h_t(x + y + (u + v) t) - u . grad phi - phi div u 395c4762a1bSJed Brown = u h' + v h' - u h_x - v h_y 396c4762a1bSJed Brown = 0 397c4762a1bSJed Brown 398c4762a1bSJed Brown We will conserve phi since 399c4762a1bSJed Brown 400c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 401c4762a1bSJed Brown 402c4762a1bSJed Brown Also try h((x + ut)^2 + (y + vt)^2), so that 403c4762a1bSJed Brown 404c4762a1bSJed Brown h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 405c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 406c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 407c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 408c4762a1bSJed Brown = 0 409c4762a1bSJed Brown 410c4762a1bSJed Brown */ 411c4762a1bSJed Brown static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 412c4762a1bSJed Brown { 413c4762a1bSJed Brown u[0] = x[0]*x[0] + x[1]*x[1]; 414c4762a1bSJed Brown u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 415c4762a1bSJed Brown return 0; 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 418c4762a1bSJed Brown /* 419c4762a1bSJed Brown In 2D we use the exact, periodic solution: 420c4762a1bSJed Brown 421c4762a1bSJed Brown u = sin(2 pi x)/4 pi^2 422c4762a1bSJed Brown v = -y cos(2 pi x)/2 pi 423c4762a1bSJed Brown phi = h(x + y + (u + v) t) 424c4762a1bSJed Brown f_x = -sin(2 pi x) 425c4762a1bSJed Brown f_y = 2 pi y cos(2 pi x) 426c4762a1bSJed Brown 427c4762a1bSJed Brown so that 428c4762a1bSJed Brown 429c4762a1bSJed Brown -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 430c4762a1bSJed Brown 431c4762a1bSJed Brown We will conserve phi since 432c4762a1bSJed Brown 433c4762a1bSJed Brown \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 434c4762a1bSJed Brown */ 435c4762a1bSJed Brown static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 436c4762a1bSJed Brown { 437c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 438c4762a1bSJed Brown u[1] = -x[1]*PetscCosReal(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI); 439c4762a1bSJed Brown return 0; 440c4762a1bSJed Brown } 441c4762a1bSJed Brown 442c4762a1bSJed Brown /* 443c4762a1bSJed Brown In 2D we use the exact, doubly periodic solution: 444c4762a1bSJed Brown 445c4762a1bSJed Brown u = sin(2 pi x) cos(2 pi y)/4 pi^2 446c4762a1bSJed Brown v = -sin(2 pi y) cos(2 pi x)/4 pi^2 447c4762a1bSJed Brown phi = h(x + y + (u + v) t) 448c4762a1bSJed Brown f_x = -2sin(2 pi x) cos(2 pi y) 449c4762a1bSJed Brown f_y = 2sin(2 pi y) cos(2 pi x) 450c4762a1bSJed Brown 451c4762a1bSJed Brown so that 452c4762a1bSJed Brown 453c4762a1bSJed 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 454c4762a1bSJed Brown 455c4762a1bSJed Brown We will conserve phi since 456c4762a1bSJed Brown 457c4762a1bSJed Brown \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 458c4762a1bSJed Brown */ 459c4762a1bSJed Brown static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 460c4762a1bSJed Brown { 461c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1])/PetscSqr(2.0*PETSC_PI); 462c4762a1bSJed Brown u[1] = -PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 463c4762a1bSJed Brown return 0; 464c4762a1bSJed Brown } 465c4762a1bSJed Brown 466c4762a1bSJed Brown static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 467c4762a1bSJed Brown { 468c4762a1bSJed Brown u[0] = x[1] - 0.5; 469c4762a1bSJed Brown u[1] = 0.0; 470c4762a1bSJed Brown return 0; 471c4762a1bSJed Brown } 472c4762a1bSJed Brown 473c4762a1bSJed Brown static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 474c4762a1bSJed Brown { 475c4762a1bSJed Brown PetscInt d; 476c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 477c4762a1bSJed Brown return 0; 478c4762a1bSJed Brown } 479c4762a1bSJed Brown 480c4762a1bSJed Brown static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 481c4762a1bSJed Brown { 482c4762a1bSJed Brown u[0] = 0.0; 483c4762a1bSJed Brown return 0; 484c4762a1bSJed Brown } 485c4762a1bSJed Brown 486c4762a1bSJed Brown static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 487c4762a1bSJed Brown { 488c4762a1bSJed Brown u[0] = 1.0; 489c4762a1bSJed Brown return 0; 490c4762a1bSJed Brown } 491c4762a1bSJed Brown 492c4762a1bSJed Brown static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 493c4762a1bSJed Brown { 494c4762a1bSJed Brown PetscReal x0[2]; 495c4762a1bSJed Brown PetscScalar xn[2]; 496c4762a1bSJed Brown 497c4762a1bSJed Brown x0[0] = globalUser->source[0]; 498c4762a1bSJed Brown x0[1] = globalUser->source[1]; 499c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 500c4762a1bSJed Brown { 501c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 502c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 503c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 504c4762a1bSJed Brown 505c4762a1bSJed Brown u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 506c4762a1bSJed Brown } 507c4762a1bSJed Brown return 0; 508c4762a1bSJed Brown } 509c4762a1bSJed Brown 510c4762a1bSJed Brown /* 511c4762a1bSJed Brown Gaussian blob, initially centered on (0.5, 0.5) 512c4762a1bSJed Brown 513c4762a1bSJed Brown xi = x(t) - x0, eta = y(t) - y0 514c4762a1bSJed Brown 515c4762a1bSJed Brown where x(t), y(t) are the integral curves of v(t), 516c4762a1bSJed Brown 517c4762a1bSJed Brown dx/dt . grad f = v . f 518c4762a1bSJed Brown 519c4762a1bSJed Brown Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 520c4762a1bSJed Brown 521c4762a1bSJed Brown v0 f_x + w0 f_y = v . f 522c4762a1bSJed Brown */ 523c4762a1bSJed Brown static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 524c4762a1bSJed Brown { 525c4762a1bSJed Brown const PetscReal x0[2] = {0.5, 0.5}; 526c4762a1bSJed Brown const PetscReal sigma = 1.0/6.0; 527c4762a1bSJed Brown PetscScalar xn[2]; 528c4762a1bSJed Brown 529c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 530c4762a1bSJed Brown { 531c4762a1bSJed Brown /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 532c4762a1bSJed Brown /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 533c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 534c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 535c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 536c4762a1bSJed Brown 537c4762a1bSJed Brown u[0] = PetscExpReal(-r2/(2.0*sigma*sigma))/(sigma*PetscSqrtReal(2.0*PETSC_PI)); 538c4762a1bSJed Brown } 539c4762a1bSJed Brown return 0; 540c4762a1bSJed Brown } 541c4762a1bSJed Brown 542c4762a1bSJed Brown static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 543c4762a1bSJed Brown { 544c4762a1bSJed Brown PetscReal x0[3]; 545c4762a1bSJed Brown const PetscReal wind[3] = {0.0, 1.0, 0.0}; 546c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 547c4762a1bSJed Brown 548c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 549c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 550c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 551c4762a1bSJed Brown return 0; 552c4762a1bSJed Brown } 553c4762a1bSJed Brown 554c4762a1bSJed Brown static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 555c4762a1bSJed Brown { 556c4762a1bSJed Brown PetscReal ur[3]; 557c4762a1bSJed Brown PetscReal x0[3]; 558c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 559c4762a1bSJed Brown 560c4762a1bSJed Brown ur[0] = PetscRealPart(u[0]); ur[1] = PetscRealPart(u[1]); ur[2] = PetscRealPart(u[2]); 561c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 562c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 563c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 564c4762a1bSJed Brown return 0; 565c4762a1bSJed Brown } 566c4762a1bSJed Brown 567c4762a1bSJed Brown static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 568c4762a1bSJed Brown { 569c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 570c4762a1bSJed Brown 571c4762a1bSJed Brown PetscFunctionBeginUser; 572c4762a1bSJed Brown xG[0] = user->inflowState; 573c4762a1bSJed Brown PetscFunctionReturn(0); 574c4762a1bSJed Brown } 575c4762a1bSJed Brown 576c4762a1bSJed Brown static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 577c4762a1bSJed Brown { 578c4762a1bSJed Brown PetscFunctionBeginUser; 57930602db0SMatthew G. Knepley //xG[0] = xI[dim]; 58030602db0SMatthew G. Knepley xG[0] = xI[2]; 581c4762a1bSJed Brown PetscFunctionReturn(0); 582c4762a1bSJed Brown } 583c4762a1bSJed Brown 584c4762a1bSJed Brown static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 585c4762a1bSJed Brown { 586c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 587c4762a1bSJed Brown PetscInt dim; 588c4762a1bSJed Brown 589c4762a1bSJed Brown PetscFunctionBeginUser; 5905f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 591c4762a1bSJed Brown switch (user->porosityDist) { 592c4762a1bSJed Brown case TILTED: 593c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) tilted_phi_2d(dim, time, x, 2, u, (void *) &time); 594c4762a1bSJed Brown else tilted_phi_coupled_2d(dim, time, x, 2, u, (void *) &time); 595c4762a1bSJed Brown break; 596c4762a1bSJed Brown case GAUSSIAN: 597c4762a1bSJed Brown gaussian_phi_2d(dim, time, x, 2, u, (void *) &time); 598c4762a1bSJed Brown break; 599c4762a1bSJed Brown case DELTA: 600c4762a1bSJed Brown delta_phi_2d(dim, time, x, 2, u, (void *) &time); 601c4762a1bSJed Brown break; 602c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 603c4762a1bSJed Brown } 604c4762a1bSJed Brown PetscFunctionReturn(0); 605c4762a1bSJed Brown } 606c4762a1bSJed Brown 607c4762a1bSJed Brown static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 608c4762a1bSJed Brown { 609c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 610c4762a1bSJed Brown PetscScalar yexact[3]={0,0,0}; 611c4762a1bSJed Brown 612c4762a1bSJed Brown PetscFunctionBeginUser; 6135f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(dm, time, x, yexact, ctx)); 614c4762a1bSJed Brown f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 615c4762a1bSJed Brown PetscFunctionReturn(0); 616c4762a1bSJed Brown } 617c4762a1bSJed Brown 618c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 619c4762a1bSJed Brown { 620c4762a1bSJed Brown PetscFunctionBeginUser; 6215f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 6225f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 62330602db0SMatthew G. Knepley #if 0 62430602db0SMatthew G. Knepley PetscBool periodic = user->bd[0] == DM_BOUNDARY_PERIODIC || user->bd[0] == DM_BOUNDARY_TWIST || user->bd[1] == DM_BOUNDARY_PERIODIC || user->bd[1] == DM_BOUNDARY_TWIST ? PETSC_TRUE : PETSC_FALSE; 62530602db0SMatthew G. Knepley const PetscReal L[3] = {1.0, 1.0, 1.0}; 62630602db0SMatthew G. Knepley PetscReal maxCell[3]; 62730602db0SMatthew G. Knepley PetscInt d; 62830602db0SMatthew G. Knepley 6295f80ce2aSJacob Faibussowitsch if (periodic) {for (d = 0; d < 3; ++d) maxCell[d] = 1.1*(L[d]/cells[d]); CHKERRQ(DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, user->bd));} 63030602db0SMatthew G. Knepley #endif 6315f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 6325f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 633c4762a1bSJed Brown PetscFunctionReturn(0); 634c4762a1bSJed Brown } 635c4762a1bSJed Brown 636c4762a1bSJed Brown static PetscErrorCode SetupBC(DM dm, AppCtx *user) 637c4762a1bSJed Brown { 638348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 63930602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 640c4762a1bSJed Brown PetscDS prob; 641c4762a1bSJed Brown DMLabel label; 642c4762a1bSJed Brown PetscBool check; 64330602db0SMatthew G. Knepley PetscInt dim, n = 3; 64430602db0SMatthew G. Knepley const char *prefix; 645c4762a1bSJed Brown 646c4762a1bSJed Brown PetscFunctionBeginUser; 6475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 6485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 6495f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 650c4762a1bSJed Brown /* Set initial guesses and exact solutions */ 65130602db0SMatthew G. Knepley switch (dim) { 652c4762a1bSJed Brown case 2: 653c4762a1bSJed Brown user->initialGuess[0] = initialVelocity; 654c4762a1bSJed Brown switch(user->porosityDist) { 655c4762a1bSJed Brown case ZERO: user->initialGuess[1] = zero_phi;break; 656c4762a1bSJed Brown case CONSTANT: user->initialGuess[1] = constant_phi;break; 657c4762a1bSJed Brown case GAUSSIAN: user->initialGuess[1] = gaussian_phi_2d;break; 658c4762a1bSJed Brown case DELTA: user->initialGuess[1] = delta_phi_2d;break; 659c4762a1bSJed Brown case TILTED: 660c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 661c4762a1bSJed Brown else user->initialGuess[1] = tilted_phi_coupled_2d; 662c4762a1bSJed Brown break; 663c4762a1bSJed Brown } 664348a1646SMatthew G. Knepley break; 66598921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %D not supported", dim); 666348a1646SMatthew G. Knepley } 667348a1646SMatthew G. Knepley exactFuncs[0] = user->initialGuess[0]; 668348a1646SMatthew G. Knepley exactFuncs[1] = user->initialGuess[1]; 66930602db0SMatthew G. Knepley switch (dim) { 670348a1646SMatthew G. Knepley case 2: 671c4762a1bSJed Brown switch (user->velocityDist) { 672c4762a1bSJed Brown case VEL_ZERO: 673348a1646SMatthew G. Knepley exactFuncs[0] = zero_u_2d; break; 674c4762a1bSJed Brown case VEL_CONSTANT: 675348a1646SMatthew G. Knepley exactFuncs[0] = constant_u_2d; break; 676c4762a1bSJed Brown case VEL_HARMONIC: 67730602db0SMatthew G. Knepley switch (bdt[0]) { 678c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 67930602db0SMatthew G. Knepley switch (bdt[1]) { 680c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 681348a1646SMatthew G. Knepley exactFuncs[0] = doubly_periodic_u_2d; break; 682c4762a1bSJed Brown default: 683348a1646SMatthew G. Knepley exactFuncs[0] = periodic_u_2d; break; 684c4762a1bSJed Brown } 685c4762a1bSJed Brown break; 686c4762a1bSJed Brown default: 687348a1646SMatthew G. Knepley exactFuncs[0] = quadratic_u_2d; break; 688c4762a1bSJed Brown } 689c4762a1bSJed Brown break; 690c4762a1bSJed Brown case VEL_SHEAR: 691348a1646SMatthew G. Knepley exactFuncs[0] = shear_bc; break; 69298921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim); 693c4762a1bSJed Brown } 694348a1646SMatthew G. Knepley break; 69598921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %D not supported", dim); 696c4762a1bSJed Brown } 697c4762a1bSJed Brown { 698c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 699c4762a1bSJed Brown 7005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 701348a1646SMatthew G. Knepley if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 702c4762a1bSJed Brown } 7035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-dmts_check", &check)); 704c4762a1bSJed Brown if (check) { 705348a1646SMatthew G. Knepley user->initialGuess[0] = exactFuncs[0]; 706348a1646SMatthew G. Knepley user->initialGuess[1] = exactFuncs[1]; 707c4762a1bSJed Brown } 708c4762a1bSJed Brown /* Set BC */ 7095f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &prob)); 7105f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 7115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user)); 7125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user)); 713c4762a1bSJed Brown if (label) { 714c4762a1bSJed Brown const PetscInt id = 1; 715c4762a1bSJed Brown 7165f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL)); 717c4762a1bSJed Brown } 7185f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "Face Sets", &label)); 719c4762a1bSJed Brown if (label && user->useFV) { 720c4762a1bSJed Brown const PetscInt inflowids[] = {100,200,300}, outflowids[] = {101}; 721c4762a1bSJed Brown 7225f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "inflow", label, ALEN(inflowids), inflowids, 1, 0, NULL, (void (*)(void)) advect_inflow, NULL, user, NULL)); 7235f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL_RIEMANN, "outflow", label, ALEN(outflowids), outflowids, 1, 0, NULL, (void (*)(void)) advect_outflow, NULL, user, NULL)); 724c4762a1bSJed Brown } 725c4762a1bSJed Brown PetscFunctionReturn(0); 726c4762a1bSJed Brown } 727c4762a1bSJed Brown 728c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 729c4762a1bSJed Brown { 73030602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 731c4762a1bSJed Brown PetscDS prob; 73230602db0SMatthew G. Knepley PetscInt n = 3; 73330602db0SMatthew G. Knepley const char *prefix; 734c4762a1bSJed Brown 735c4762a1bSJed Brown PetscFunctionBeginUser; 7365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 7375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 7385f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &prob)); 739c4762a1bSJed Brown switch (user->velocityDist) { 740c4762a1bSJed Brown case VEL_ZERO: 7415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u)); 742c4762a1bSJed Brown break; 743c4762a1bSJed Brown case VEL_CONSTANT: 7445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u)); 7455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL)); 7465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL)); 747c4762a1bSJed Brown break; 748c4762a1bSJed Brown case VEL_HARMONIC: 74930602db0SMatthew G. Knepley switch (bdt[0]) { 750c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 75130602db0SMatthew G. Knepley switch (bdt[1]) { 752c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 7535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u)); 754c4762a1bSJed Brown break; 755c4762a1bSJed Brown default: 7565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u)); 757c4762a1bSJed Brown break; 758c4762a1bSJed Brown } 759c4762a1bSJed Brown break; 760c4762a1bSJed Brown default: 7615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u)); 762c4762a1bSJed Brown break; 763c4762a1bSJed Brown } 7645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 765c4762a1bSJed Brown break; 766c4762a1bSJed Brown case VEL_SHEAR: 7675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u)); 7685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 769c4762a1bSJed Brown break; 770c4762a1bSJed Brown } 7715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 1, f0_advection, f1_advection)); 7725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL)); 7735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL)); 7745f80ce2aSJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) CHKERRQ(PetscDSSetRiemannSolver(prob, 1, riemann_advection)); 7755f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection)); 776c4762a1bSJed Brown 7775f80ce2aSJacob Faibussowitsch CHKERRQ(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user)); 778c4762a1bSJed Brown PetscFunctionReturn(0); 779c4762a1bSJed Brown } 780c4762a1bSJed Brown 781c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 782c4762a1bSJed Brown { 783c4762a1bSJed Brown DM cdm = dm; 784c4762a1bSJed Brown PetscQuadrature q; 785c4762a1bSJed Brown PetscFE fe[2]; 786c4762a1bSJed Brown PetscFV fv; 787c4762a1bSJed Brown MPI_Comm comm; 78830602db0SMatthew G. Knepley PetscInt dim; 789c4762a1bSJed Brown 790c4762a1bSJed Brown PetscFunctionBeginUser; 791c4762a1bSJed Brown /* Create finite element */ 7925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 7935f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 7945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0])); 7955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[0], "velocity")); 7965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1])); 7975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe[0], fe[1])); 7985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[1], "porosity")); 799c4762a1bSJed Brown 8005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVCreate(PetscObjectComm((PetscObject) dm), &fv)); 8015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fv, "porosity")); 8025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetFromOptions(fv)); 8035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetNumComponents(fv, 1)); 8045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetSpatialDimension(fv, dim)); 8055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetQuadrature(fe[0], &q)); 8065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetQuadrature(fv, q)); 807c4762a1bSJed Brown 8085f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 8095f80ce2aSJacob Faibussowitsch if (user->useFV) CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fv)); 8105f80ce2aSJacob Faibussowitsch else CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 8115f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 8125f80ce2aSJacob Faibussowitsch CHKERRQ(SetupProblem(dm, user)); 813c4762a1bSJed Brown 814c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 815c4762a1bSJed Brown while (cdm) { 8165f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 8175f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 818c4762a1bSJed Brown /* Coordinates were never localized for coarse meshes */ 8195f80ce2aSJacob Faibussowitsch if (cdm) CHKERRQ(DMLocalizeCoordinates(cdm)); 820c4762a1bSJed Brown } 8215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[0])); 8225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe[1])); 8235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVDestroy(&fv)); 824c4762a1bSJed Brown PetscFunctionReturn(0); 825c4762a1bSJed Brown } 826c4762a1bSJed Brown 827c4762a1bSJed Brown static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 828c4762a1bSJed Brown { 829c4762a1bSJed Brown PetscFunctionBeginUser; 8305f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(comm, user, dm)); 831c4762a1bSJed Brown /* Handle refinement, etc. */ 8325f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 833c4762a1bSJed Brown /* Construct ghost cells */ 834c4762a1bSJed Brown if (user->useFV) { 835c4762a1bSJed Brown DM gdm; 836c4762a1bSJed Brown 8375f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm)); 8385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 839c4762a1bSJed Brown *dm = gdm; 840c4762a1bSJed Brown } 841c4762a1bSJed Brown /* Localize coordinates */ 8425f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalizeCoordinates(*dm)); 8435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)(*dm),"Mesh")); 8445f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 845c4762a1bSJed Brown /* Setup problem */ 8465f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(*dm, user)); 847c4762a1bSJed Brown /* Setup BC */ 8485f80ce2aSJacob Faibussowitsch CHKERRQ(SetupBC(*dm, user)); 849c4762a1bSJed Brown PetscFunctionReturn(0); 850c4762a1bSJed Brown } 851c4762a1bSJed Brown 852c4762a1bSJed Brown static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx) 853c4762a1bSJed Brown { 854c4762a1bSJed Brown PetscDS prob; 855c4762a1bSJed Brown DM dmCell; 856c4762a1bSJed Brown Vec cellgeom; 857c4762a1bSJed Brown const PetscScalar *cgeom; 858c4762a1bSJed Brown PetscScalar *x; 859c4762a1bSJed Brown PetscInt dim, Nf, cStart, cEnd, c; 860c4762a1bSJed Brown 861c4762a1bSJed Brown PetscFunctionBeginUser; 8625f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &prob)); 8635f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 8645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetNumFields(prob, &Nf)); 8655f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8665f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(cellgeom, &dmCell)); 8675f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 8685f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(cellgeom, &cgeom)); 8695f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(X, &x)); 870c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 871c4762a1bSJed Brown PetscFVCellGeom *cg; 872c4762a1bSJed Brown PetscScalar *xc; 873c4762a1bSJed Brown 8745f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 8755f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc)); 8765f80ce2aSJacob Faibussowitsch if (xc) CHKERRQ((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx)); 877c4762a1bSJed Brown } 8785f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(cellgeom, &cgeom)); 8795f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(X, &x)); 880c4762a1bSJed Brown PetscFunctionReturn(0); 881c4762a1bSJed Brown } 882c4762a1bSJed Brown 883c4762a1bSJed Brown static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 884c4762a1bSJed Brown { 885c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 886c4762a1bSJed Brown char *ftable = NULL; 887c4762a1bSJed Brown DM dm; 888c4762a1bSJed Brown PetscSection s; 889c4762a1bSJed Brown Vec cellgeom; 890c4762a1bSJed Brown const PetscScalar *x; 891c4762a1bSJed Brown PetscScalar *a; 892c4762a1bSJed Brown PetscReal *xnorms; 893c4762a1bSJed Brown PetscInt pStart, pEnd, p, Nf, f; 894c4762a1bSJed Brown 895c4762a1bSJed Brown PetscFunctionBeginUser; 8965f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(X, (PetscObject) ts, "-view_solution")); 8975f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(X, &dm)); 8985f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8995f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(dm, &s)); 9005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetNumFields(s, &Nf)); 9015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetChart(s, &pStart, &pEnd)); 9025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(Nf*2, &xnorms)); 9035f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X, &x)); 904c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 905c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 906c4762a1bSJed Brown PetscInt dof, cdof, d; 907c4762a1bSJed Brown 9085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetFieldDof(s, p, f, &dof)); 9095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 9105f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointGlobalFieldRead(dm, p, f, x, &a)); 911c4762a1bSJed Brown /* TODO Use constrained indices here */ 912c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+0] = PetscMax(xnorms[f*2+0], PetscAbsScalar(a[d])); 913c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+1] += PetscAbsScalar(a[d]); 914c4762a1bSJed Brown } 915c4762a1bSJed Brown } 9165f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X, &x)); 917c4762a1bSJed Brown if (stepnum >= 0) { /* No summary for final time */ 918c4762a1bSJed Brown DM dmCell, *fdm; 919c4762a1bSJed Brown Vec *fv; 920c4762a1bSJed Brown const PetscScalar *cgeom; 921c4762a1bSJed Brown PetscScalar **fx; 922c4762a1bSJed Brown PetscReal *fmin, *fmax, *fint, *ftmp, t; 923c4762a1bSJed Brown PetscInt cStart, cEnd, c, fcount, f, num; 924c4762a1bSJed Brown 925c4762a1bSJed Brown size_t ftableused,ftablealloc; 926c4762a1bSJed Brown 927c4762a1bSJed Brown /* Functionals have indices after registering, this is an upper bound */ 928c4762a1bSJed Brown fcount = user->numMonitorFuncs; 9295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fint,fcount,&ftmp)); 9305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(fcount,&fdm,fcount,&fv,fcount,&fx)); 931c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 932c4762a1bSJed Brown PetscSection fs; 933c4762a1bSJed Brown const char *name = user->monitorFuncs[f]->name; 934c4762a1bSJed Brown 935c4762a1bSJed Brown fmin[f] = PETSC_MAX_REAL; 936c4762a1bSJed Brown fmax[f] = PETSC_MIN_REAL; 937c4762a1bSJed Brown fint[f] = 0; 938c4762a1bSJed Brown /* Make monitor vecs */ 9395f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, &fdm[f])); 9405f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetOutputSequenceNumber(dm, &num, &t)); 9415f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(fdm[f], num, t)); 9425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionClone(s, &fs)); 9435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldName(fs, 0, NULL)); 9445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldName(fs, 1, name)); 9455f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(fdm[f], fs)); 9465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&fs)); 9475f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(fdm[f], &fv[f])); 9485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fv[f], name)); 9495f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(fv[f], &fx[f])); 950c4762a1bSJed Brown } 9515f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 9525f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(cellgeom, &dmCell)); 9535f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(cellgeom, &cgeom)); 9545f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X, &x)); 955c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 956c4762a1bSJed Brown PetscFVCellGeom *cg; 957c4762a1bSJed Brown PetscScalar *cx; 958c4762a1bSJed Brown 9595f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 9605f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx)); 961c4762a1bSJed Brown if (!cx) continue; /* not a global cell */ 962c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 963c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 964c4762a1bSJed Brown PetscScalar *fxc; 965c4762a1bSJed Brown 9665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc)); 967c4762a1bSJed Brown /* I need to make it easier to get interpolated values here */ 9685f80ce2aSJacob Faibussowitsch CHKERRQ((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx)); 969c4762a1bSJed Brown fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 970c4762a1bSJed Brown } 971c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 972c4762a1bSJed Brown fmin[f] = PetscMin(fmin[f], ftmp[f]); 973c4762a1bSJed Brown fmax[f] = PetscMax(fmax[f], ftmp[f]); 974c4762a1bSJed Brown fint[f] += cg->volume * ftmp[f]; 975c4762a1bSJed Brown } 976c4762a1bSJed Brown } 9775f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(cellgeom, &cgeom)); 9785f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X, &x)); 9795f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 9805f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 9815f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 982c4762a1bSJed Brown /* Output functional data */ 983c4762a1bSJed Brown ftablealloc = fcount * 100; 984c4762a1bSJed Brown ftableused = 0; 9855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(ftablealloc, &ftable)); 986c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 987c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 988c4762a1bSJed Brown PetscInt id = func->offset; 989c4762a1bSJed Brown char newline[] = "\n"; 990c4762a1bSJed Brown char buffer[256], *p, *prefix; 991c4762a1bSJed Brown size_t countused, len; 992c4762a1bSJed Brown 993c4762a1bSJed Brown /* Create string with functional outputs */ 994c4762a1bSJed Brown if (f % 3) { 9955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(buffer, " ", 2)); 996c4762a1bSJed Brown p = buffer + 2; 997c4762a1bSJed Brown } else if (f) { 9985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(buffer, newline, sizeof(newline)-1)); 999c4762a1bSJed Brown p = buffer + sizeof(newline) - 1; 1000c4762a1bSJed Brown } else { 1001c4762a1bSJed Brown p = buffer; 1002c4762a1bSJed Brown } 10035f80ce2aSJacob Faibussowitsch CHKERRQ(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])); 1004c4762a1bSJed Brown countused += p - buffer; 1005c4762a1bSJed Brown /* reallocate */ 1006c4762a1bSJed Brown if (countused > ftablealloc-ftableused-1) { 1007c4762a1bSJed Brown char *ftablenew; 1008c4762a1bSJed Brown 1009c4762a1bSJed Brown ftablealloc = 2*ftablealloc + countused; 10105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(ftablealloc, &ftablenew)); 10115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(ftablenew, ftable, ftableused)); 10125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ftable)); 1013c4762a1bSJed Brown ftable = ftablenew; 1014c4762a1bSJed Brown } 10155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(ftable+ftableused, buffer, countused)); 1016c4762a1bSJed Brown ftableused += countused; 1017c4762a1bSJed Brown ftable[ftableused] = 0; 1018c4762a1bSJed Brown /* Output vecs */ 10195f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(fv[f], &fx[f])); 10205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlen(func->name, &len)); 10215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(len+2,&prefix)); 10225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(prefix, func->name)); 10235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(prefix, "_")); 10245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix)); 10255f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(fv[f], NULL, "-vec_view")); 10265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(prefix)); 10275f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(fdm[f], &fv[f])); 10285f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&fdm[f])); 1029c4762a1bSJed Brown } 10305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree4(fmin, fmax, fint, ftmp)); 10315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(fdm, fv, fx)); 10325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), "% 3D time %8.4g |x| (", stepnum, (double) time)); 1033c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10345f80ce2aSJacob Faibussowitsch if (f > 0) CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 10355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+0])); 1036c4762a1bSJed Brown } 10375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), ") |x|_1 (")); 1038c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10395f80ce2aSJacob Faibussowitsch if (f > 0) CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 10405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+1])); 1041c4762a1bSJed Brown } 10425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), ") %s\n", ftable ? ftable : "")); 10435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ftable)); 1044c4762a1bSJed Brown } 10455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(xnorms)); 1046c4762a1bSJed Brown PetscFunctionReturn(0); 1047c4762a1bSJed Brown } 1048c4762a1bSJed Brown 1049c4762a1bSJed Brown int main(int argc, char **argv) 1050c4762a1bSJed Brown { 1051c4762a1bSJed Brown MPI_Comm comm; 1052c4762a1bSJed Brown TS ts; 1053c4762a1bSJed Brown DM dm; 1054c4762a1bSJed Brown Vec u; 1055c4762a1bSJed Brown AppCtx user; 1056c4762a1bSJed Brown PetscReal t0, t = 0.0; 1057c4762a1bSJed Brown void *ctxs[2]; 1058c4762a1bSJed Brown 1059c4762a1bSJed Brown ctxs[0] = &t; 1060c4762a1bSJed Brown ctxs[1] = &t; 1061*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, (char*) 0, help)); 1062c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1063c4762a1bSJed Brown user.functionalRegistry = NULL; 1064c4762a1bSJed Brown globalUser = &user; 10655f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(comm, &user)); 10665f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(comm, &ts)); 10675f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts, TSBEULER)); 10685f80ce2aSJacob Faibussowitsch CHKERRQ(CreateDM(comm, &user, &dm)); 10695f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, dm)); 10705f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessMonitorOptions(comm, &user)); 1071c4762a1bSJed Brown 10725f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 10735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "solution")); 1074c4762a1bSJed Brown if (user.useFV) { 1075c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 1076c4762a1bSJed Brown 10775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 1078c4762a1bSJed Brown if (isImplicit) { 10795f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10805f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1081c4762a1bSJed Brown } 10825f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10835f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user)); 1084c4762a1bSJed Brown } else { 10855f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10865f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10875f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1088c4762a1bSJed Brown } 10895f80ce2aSJacob Faibussowitsch if (user.useFV) CHKERRQ(TSMonitorSet(ts, MonitorFunctionals, &user, NULL)); 10905f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts, 1)); 10915f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts, 2.0)); 10925f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,0.01)); 10935f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 10945f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 1095c4762a1bSJed Brown 10965f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u)); 10975f80ce2aSJacob Faibussowitsch if (user.useFV) CHKERRQ(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1])); 10985f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-init_vec_view")); 10995f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts, &t)); 1100c4762a1bSJed Brown t0 = t; 11015f80ce2aSJacob Faibussowitsch CHKERRQ(DMTSCheckFromOptions(ts, u)); 11025f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts, u)); 11035f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts, &t)); 11045f80ce2aSJacob Faibussowitsch if (t > t0) CHKERRQ(DMTSCheckFromOptions(ts, u)); 11055f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view")); 1106c4762a1bSJed Brown { 1107c4762a1bSJed Brown PetscReal ftime; 1108c4762a1bSJed Brown PetscInt nsteps; 1109c4762a1bSJed Brown TSConvergedReason reason; 1110c4762a1bSJed Brown 11115f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts, &ftime)); 11125f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts, &nsteps)); 11135f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetConvergedReason(ts, &reason)); 11145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %D steps\n", TSConvergedReasons[reason], (double) ftime, nsteps)); 1115c4762a1bSJed Brown } 1116c4762a1bSJed Brown 11175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 11185f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 11195f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 11205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.monitorFuncs)); 11215f80ce2aSJacob Faibussowitsch CHKERRQ(FunctionalDestroy(&user.functionalRegistry)); 1122*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 1123*b122ec5aSJacob Faibussowitsch return 0; 1124c4762a1bSJed Brown } 1125c4762a1bSJed Brown 1126c4762a1bSJed Brown /*TEST 1127c4762a1bSJed Brown 112830602db0SMatthew G. Knepley testset: 112930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 113030602db0SMatthew G. Knepley 1131c4762a1bSJed Brown # 2D harmonic velocity, no porosity 1132c4762a1bSJed Brown test: 1133c4762a1bSJed Brown suffix: p1p1 1134c4762a1bSJed Brown requires: !complex !single 113530602db0SMatthew G. Knepley args: -velocity_petscspace_degree 1 -porosity_petscspace_degree 1 -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type greedy -pc_factor_shift_type nonzero -ts_monitor -snes_error_if_not_converged -ksp_error_if_not_converged -dmts_check 1136c4762a1bSJed Brown test: 1137c4762a1bSJed Brown suffix: p1p1_xper 1138c4762a1bSJed Brown requires: !complex !single 113930602db0SMatthew 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 1140c4762a1bSJed Brown test: 1141c4762a1bSJed Brown suffix: p1p1_xper_ref 1142c4762a1bSJed Brown requires: !complex !single 114330602db0SMatthew 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 1144c4762a1bSJed Brown test: 1145c4762a1bSJed Brown suffix: p1p1_xyper 1146c4762a1bSJed Brown requires: !complex !single 114730602db0SMatthew 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 1148c4762a1bSJed Brown test: 1149c4762a1bSJed Brown suffix: p1p1_xyper_ref 1150c4762a1bSJed Brown requires: !complex !single 115130602db0SMatthew 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 1152c4762a1bSJed Brown test: 1153c4762a1bSJed Brown suffix: p2p1 1154c4762a1bSJed Brown requires: !complex !single 115530602db0SMatthew 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 1156c4762a1bSJed Brown test: 1157c4762a1bSJed Brown suffix: p2p1_xyper 1158c4762a1bSJed Brown requires: !complex !single 115930602db0SMatthew 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 116030602db0SMatthew G. Knepley 116130602db0SMatthew G. Knepley test: 116230602db0SMatthew G. Knepley suffix: adv_1 116330602db0SMatthew G. Knepley requires: !complex !single 116430602db0SMatthew 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 116530602db0SMatthew G. Knepley 116630602db0SMatthew G. Knepley test: 116730602db0SMatthew G. Knepley suffix: adv_2 116830602db0SMatthew G. Knepley requires: !complex 116930602db0SMatthew G. Knepley TODO: broken memory corruption 117030602db0SMatthew 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 117130602db0SMatthew G. Knepley 117230602db0SMatthew G. Knepley test: 117330602db0SMatthew G. Knepley suffix: adv_3 117430602db0SMatthew G. Knepley requires: !complex 117530602db0SMatthew G. Knepley TODO: broken memory corruption 117630602db0SMatthew 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 117730602db0SMatthew G. Knepley 117830602db0SMatthew G. Knepley test: 117930602db0SMatthew G. Knepley suffix: adv_3_ex 118030602db0SMatthew G. Knepley requires: !complex 118130602db0SMatthew 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 118230602db0SMatthew G. Knepley 118330602db0SMatthew G. Knepley test: 118430602db0SMatthew G. Knepley suffix: adv_4 118530602db0SMatthew G. Knepley requires: !complex 118630602db0SMatthew G. Knepley TODO: broken memory corruption 118730602db0SMatthew 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 118830602db0SMatthew G. Knepley 118930602db0SMatthew G. Knepley # 2D Advection, box, delta 119030602db0SMatthew G. Knepley test: 119130602db0SMatthew G. Knepley suffix: adv_delta_yper_0 119230602db0SMatthew G. Knepley requires: !complex 119330602db0SMatthew G. Knepley TODO: broken 119430602db0SMatthew 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 119530602db0SMatthew G. Knepley 119630602db0SMatthew G. Knepley test: 119730602db0SMatthew G. Knepley suffix: adv_delta_yper_1 119830602db0SMatthew G. Knepley requires: !complex 119930602db0SMatthew G. Knepley TODO: broken 120030602db0SMatthew 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 120130602db0SMatthew G. Knepley 120230602db0SMatthew G. Knepley test: 120330602db0SMatthew G. Knepley suffix: adv_delta_yper_2 120430602db0SMatthew G. Knepley requires: !complex 120530602db0SMatthew G. Knepley TODO: broken 120630602db0SMatthew 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 120730602db0SMatthew G. Knepley 120830602db0SMatthew G. Knepley test: 120930602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_0 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 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 121330602db0SMatthew G. Knepley 121430602db0SMatthew G. Knepley test: 121530602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_1 121630602db0SMatthew G. Knepley requires: !complex 121730602db0SMatthew G. Knepley TODO: broken 121830602db0SMatthew 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 121930602db0SMatthew G. Knepley 122030602db0SMatthew G. Knepley test: 122130602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_2 122230602db0SMatthew G. Knepley requires: !complex 122330602db0SMatthew G. Knepley TODO: broken 122430602db0SMatthew 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 122530602db0SMatthew G. Knepley 122630602db0SMatthew G. Knepley test: 122730602db0SMatthew G. Knepley suffix: adv_delta_yper_im_0 122830602db0SMatthew G. Knepley requires: !complex 122930602db0SMatthew G. Knepley TODO: broken 123030602db0SMatthew 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 123130602db0SMatthew G. Knepley 123230602db0SMatthew G. Knepley test: 123330602db0SMatthew G. Knepley suffix: adv_delta_yper_im_1 123430602db0SMatthew G. Knepley requires: !complex 123530602db0SMatthew G. Knepley TODO: broken 123630602db0SMatthew 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 123730602db0SMatthew G. Knepley 123830602db0SMatthew G. Knepley test: 123930602db0SMatthew G. Knepley suffix: adv_delta_yper_im_2 124030602db0SMatthew G. Knepley requires: !complex 124130602db0SMatthew G. Knepley TODO: broken 124230602db0SMatthew 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 124330602db0SMatthew G. Knepley 124430602db0SMatthew G. Knepley test: 124530602db0SMatthew G. Knepley suffix: adv_delta_yper_im_3 124630602db0SMatthew G. Knepley requires: !complex 124730602db0SMatthew G. Knepley TODO: broken 124830602db0SMatthew 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 124930602db0SMatthew G. Knepley 125030602db0SMatthew G. Knepley # I believe the nullspace is sin(pi y) 125130602db0SMatthew G. Knepley test: 125230602db0SMatthew G. Knepley suffix: adv_delta_yper_im_4 125330602db0SMatthew G. Knepley requires: !complex 125430602db0SMatthew G. Knepley TODO: broken 125530602db0SMatthew 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 125630602db0SMatthew G. Knepley 125730602db0SMatthew G. Knepley test: 125830602db0SMatthew G. Knepley suffix: adv_delta_yper_im_5 125930602db0SMatthew G. Knepley requires: !complex 126030602db0SMatthew G. Knepley TODO: broken 126130602db0SMatthew 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 126230602db0SMatthew G. Knepley 126330602db0SMatthew G. Knepley test: 126430602db0SMatthew G. Knepley suffix: adv_delta_yper_im_6 126530602db0SMatthew G. Knepley requires: !complex 126630602db0SMatthew G. Knepley TODO: broken 126730602db0SMatthew 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 126830602db0SMatthew G. Knepley # 2D Advection, magma benchmark 1 126930602db0SMatthew G. Knepley 127030602db0SMatthew G. Knepley test: 127130602db0SMatthew G. Knepley suffix: adv_delta_shear_im_0 127230602db0SMatthew G. Knepley requires: !complex 127330602db0SMatthew G. Knepley TODO: broken 127430602db0SMatthew 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 127530602db0SMatthew G. Knepley # 2D Advection, box, gaussian 127630602db0SMatthew G. Knepley 127730602db0SMatthew G. Knepley test: 127830602db0SMatthew G. Knepley suffix: adv_gauss 127930602db0SMatthew G. Knepley requires: !complex 128030602db0SMatthew G. Knepley TODO: broken 128130602db0SMatthew 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 128230602db0SMatthew G. Knepley 128330602db0SMatthew G. Knepley test: 128430602db0SMatthew G. Knepley suffix: adv_gauss_im 128530602db0SMatthew G. Knepley requires: !complex 128630602db0SMatthew G. Knepley TODO: broken 128730602db0SMatthew 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 128830602db0SMatthew G. Knepley 128930602db0SMatthew G. Knepley test: 129030602db0SMatthew G. Knepley suffix: adv_gauss_im_1 129130602db0SMatthew G. Knepley requires: !complex 129230602db0SMatthew G. Knepley TODO: broken 129330602db0SMatthew 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 129430602db0SMatthew G. Knepley 129530602db0SMatthew G. Knepley test: 129630602db0SMatthew G. Knepley suffix: adv_gauss_im_2 129730602db0SMatthew G. Knepley requires: !complex 129830602db0SMatthew G. Knepley TODO: broken 129930602db0SMatthew 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 130030602db0SMatthew G. Knepley 130130602db0SMatthew G. Knepley test: 130230602db0SMatthew G. Knepley suffix: adv_gauss_corner 130330602db0SMatthew G. Knepley requires: !complex 130430602db0SMatthew G. Knepley TODO: broken 130530602db0SMatthew 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 130630602db0SMatthew G. Knepley 130730602db0SMatthew G. Knepley # 2D Advection+Harmonic 12- 130830602db0SMatthew G. Knepley test: 130930602db0SMatthew G. Knepley suffix: adv_harm_0 131030602db0SMatthew G. Knepley requires: !complex 131130602db0SMatthew G. Knepley TODO: broken memory corruption 131230602db0SMatthew 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 131330602db0SMatthew G. Knepley 1314c4762a1bSJed Brown # Must check that FV BCs propagate to coarse meshes 1315c4762a1bSJed Brown # Must check that FV BC ids propagate to coarse meshes 1316c4762a1bSJed Brown # Must check that FE+FV BCs work at the same time 1317c4762a1bSJed Brown # 2D Advection, matching wind in ex11 8-11 1318c4762a1bSJed Brown # NOTE implicit solves are limited by accuracy of FD Jacobian 1319c4762a1bSJed Brown test: 1320c4762a1bSJed Brown suffix: adv_0 1321c4762a1bSJed Brown requires: !complex !single exodusii 132230602db0SMatthew 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 1323c4762a1bSJed Brown 1324c4762a1bSJed Brown test: 1325c4762a1bSJed Brown suffix: adv_0_im 1326c4762a1bSJed Brown requires: !complex exodusii 1327c4762a1bSJed Brown TODO: broken memory corruption 132830602db0SMatthew 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 1329c4762a1bSJed Brown 1330c4762a1bSJed Brown test: 1331c4762a1bSJed Brown suffix: adv_0_im_2 1332c4762a1bSJed Brown requires: !complex exodusii 1333c4762a1bSJed Brown TODO: broken 133430602db0SMatthew 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 1335c4762a1bSJed Brown 1336c4762a1bSJed Brown test: 1337c4762a1bSJed Brown suffix: adv_0_im_3 1338c4762a1bSJed Brown requires: !complex exodusii 1339c4762a1bSJed Brown TODO: broken 134030602db0SMatthew 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 1341c4762a1bSJed Brown 1342c4762a1bSJed Brown test: 1343c4762a1bSJed Brown suffix: adv_0_im_4 1344c4762a1bSJed Brown requires: !complex exodusii 1345c4762a1bSJed Brown TODO: broken 134630602db0SMatthew 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 1347c4762a1bSJed Brown # 2D Advection, misc 1348c4762a1bSJed Brown 1349c4762a1bSJed Brown TEST*/ 1350