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 typedef enum {VEL_ZERO, VEL_CONSTANT, VEL_HARMONIC, VEL_SHEAR} VelocityDistribution; 18c4762a1bSJed Brown 19c4762a1bSJed Brown typedef enum {ZERO, CONSTANT, GAUSSIAN, TILTED, DELTA} PorosityDistribution; 20c4762a1bSJed Brown 21c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* 24c4762a1bSJed Brown FunctionalFunc - Calculates the value of a functional of the solution at a point 25c4762a1bSJed Brown 26c4762a1bSJed Brown Input Parameters: 27c4762a1bSJed Brown + dm - The DM 28c4762a1bSJed Brown . time - The TS time 29c4762a1bSJed Brown . x - The coordinates of the evaluation point 30c4762a1bSJed Brown . u - The field values at point x 31c4762a1bSJed Brown - ctx - A user context, or NULL 32c4762a1bSJed Brown 33c4762a1bSJed Brown Output Parameter: 34c4762a1bSJed Brown . f - The value of the functional at point x 35c4762a1bSJed Brown 36c4762a1bSJed Brown */ 37c4762a1bSJed Brown typedef PetscErrorCode (*FunctionalFunc)(DM, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); 38c4762a1bSJed Brown 39c4762a1bSJed Brown typedef struct _n_Functional *Functional; 40c4762a1bSJed Brown struct _n_Functional { 41c4762a1bSJed Brown char *name; 42c4762a1bSJed Brown FunctionalFunc func; 43c4762a1bSJed Brown void *ctx; 44c4762a1bSJed Brown PetscInt offset; 45c4762a1bSJed Brown Functional next; 46c4762a1bSJed Brown }; 47c4762a1bSJed Brown 48c4762a1bSJed Brown typedef struct { 49c4762a1bSJed Brown /* Problem definition */ 50c4762a1bSJed Brown PetscBool useFV; /* Use a finite volume scheme for advection */ 51c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 52c4762a1bSJed Brown VelocityDistribution velocityDist; 53c4762a1bSJed Brown PorosityDistribution porosityDist; 54c4762a1bSJed Brown PetscReal inflowState; 55c4762a1bSJed Brown PetscReal source[3]; 56c4762a1bSJed Brown /* Monitoring */ 57c4762a1bSJed Brown PetscInt numMonitorFuncs, maxMonitorFunc; 58c4762a1bSJed Brown Functional *monitorFuncs; 59c4762a1bSJed Brown PetscInt errorFunctional; 60c4762a1bSJed Brown Functional functionalRegistry; 61c4762a1bSJed Brown } AppCtx; 62c4762a1bSJed Brown 63c4762a1bSJed Brown static AppCtx *globalUser; 64c4762a1bSJed Brown 65c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 66c4762a1bSJed Brown { 67c4762a1bSJed Brown const char *velocityDist[4] = {"zero", "constant", "harmonic", "shear"}; 68c4762a1bSJed Brown const char *porosityDist[5] = {"zero", "constant", "gaussian", "tilted", "delta"}; 6930602db0SMatthew G. Knepley PetscInt vd, pd, d; 70c4762a1bSJed Brown PetscBool flg; 71c4762a1bSJed Brown PetscErrorCode ierr; 72c4762a1bSJed Brown 73c4762a1bSJed Brown PetscFunctionBeginUser; 74c4762a1bSJed Brown options->useFV = PETSC_FALSE; 75c4762a1bSJed Brown options->velocityDist = VEL_HARMONIC; 76c4762a1bSJed Brown options->porosityDist = ZERO; 77c4762a1bSJed Brown options->inflowState = -2.0; 78c4762a1bSJed Brown options->numMonitorFuncs = 0; 79c4762a1bSJed Brown options->source[0] = 0.5; 80c4762a1bSJed Brown options->source[1] = 0.5; 81c4762a1bSJed Brown options->source[2] = 0.5; 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX");PetscCall(ierr); 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL)); 85c4762a1bSJed Brown vd = options->velocityDist; 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-velocity_dist","Velocity distribution type","ex18.c",velocityDist,4,velocityDist[options->velocityDist],&vd,NULL)); 87c4762a1bSJed Brown options->velocityDist = (VelocityDistribution) vd; 88c4762a1bSJed Brown pd = options->porosityDist; 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-porosity_dist","Initial porosity distribution type","ex18.c",porosityDist,5,porosityDist[options->porosityDist],&pd,NULL)); 90c4762a1bSJed Brown options->porosityDist = (PorosityDistribution) pd; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL)); 9230602db0SMatthew G. Knepley d = 2; 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg)); 943c633725SBarry Smith PetscCheck(!flg || d == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %d", d); 959566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscFunctionReturn(0); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) 101c4762a1bSJed Brown { 102c4762a1bSJed Brown Functional func; 103c4762a1bSJed Brown char *names[256]; 104c4762a1bSJed Brown PetscInt f; 105c4762a1bSJed Brown PetscErrorCode ierr; 106c4762a1bSJed Brown 107c4762a1bSJed Brown PetscFunctionBeginUser; 1089566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX");PetscCall(ierr); 109*dd39110bSPierre Jolivet options->numMonitorFuncs = PETSC_STATIC_ARRAY_LENGTH(names); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL)); 1119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs)); 112c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 113c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 114c4762a1bSJed Brown PetscBool match; 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(names[f], func->name, &match)); 117c4762a1bSJed Brown if (match) break; 118c4762a1bSJed Brown } 1193c633725SBarry Smith PetscCheck(func,comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 120c4762a1bSJed Brown options->monitorFuncs[f] = func; 121c4762a1bSJed Brown /* Jed inserts a de-duplication of functionals here */ 1229566063dSJacob Faibussowitsch PetscCall(PetscFree(names[f])); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 125c4762a1bSJed Brown options->maxMonitorFunc = -1; 126c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 127c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 128c4762a1bSJed Brown Functional call = options->monitorFuncs[f]; 129c4762a1bSJed Brown 130c4762a1bSJed Brown if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown } 1339566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 134c4762a1bSJed Brown PetscFunctionReturn(0); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) 138c4762a1bSJed Brown { 139c4762a1bSJed Brown Functional *ptr, f; 140c4762a1bSJed Brown PetscInt lastoffset = -1; 141c4762a1bSJed Brown 142c4762a1bSJed Brown PetscFunctionBeginUser; 143c4762a1bSJed Brown for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 1449566063dSJacob Faibussowitsch PetscCall(PetscNew(&f)); 1459566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &f->name)); 146c4762a1bSJed Brown f->offset = lastoffset + 1; 147c4762a1bSJed Brown f->func = func; 148c4762a1bSJed Brown f->ctx = ctx; 149c4762a1bSJed Brown f->next = NULL; 150c4762a1bSJed Brown *ptr = f; 151c4762a1bSJed Brown *offset = f->offset; 152c4762a1bSJed Brown PetscFunctionReturn(0); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown static PetscErrorCode FunctionalDestroy(Functional *link) 156c4762a1bSJed Brown { 157c4762a1bSJed Brown Functional next, l; 158c4762a1bSJed Brown 159c4762a1bSJed Brown PetscFunctionBeginUser; 160c4762a1bSJed Brown if (!link) PetscFunctionReturn(0); 161c4762a1bSJed Brown l = *link; 162c4762a1bSJed Brown *link = NULL; 163c4762a1bSJed Brown for (; l; l=next) { 164c4762a1bSJed Brown next = l->next; 1659566063dSJacob Faibussowitsch PetscCall(PetscFree(l->name)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFree(l)); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown PetscFunctionReturn(0); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 172c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 173c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 174c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 175c4762a1bSJed Brown { 176c4762a1bSJed Brown PetscInt comp; 177c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 181c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 182c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 183c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 184c4762a1bSJed Brown { 185c4762a1bSJed Brown PetscScalar wind[3] = {0.0, 0.0, 0.0}; 186c4762a1bSJed Brown PetscInt comp; 187c4762a1bSJed Brown 188c4762a1bSJed Brown constant_u_2d(dim, t, x, Nf, wind, NULL); 189c4762a1bSJed Brown for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp]; 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192c4762a1bSJed Brown static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 193c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 194c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 195c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown PetscInt comp; 198c4762a1bSJed Brown for (comp = 0; comp < dim*dim; ++comp) f1[comp] = 0.0; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 202c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 203c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 204c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 205c4762a1bSJed Brown { 206c4762a1bSJed Brown PetscInt d; 207c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d*dim+d] = 1.0; 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 210c4762a1bSJed Brown static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 211c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 212c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 213c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 214c4762a1bSJed Brown { 215c4762a1bSJed Brown g0[0] = 1.0; 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 219c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 220c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 221c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 222c4762a1bSJed Brown { 223c4762a1bSJed Brown PetscInt comp; 224c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227c4762a1bSJed Brown static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 228c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 229c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 230c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 231c4762a1bSJed Brown { 232c4762a1bSJed Brown PetscInt comp, d; 233c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) { 234c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 235c4762a1bSJed Brown f1[comp*dim+d] = u_x[comp*dim+d]; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown } 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 241c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 242c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 243c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 244c4762a1bSJed Brown { 245c4762a1bSJed Brown f0[0] = -PetscSinReal(2.0*PETSC_PI*x[0]); 246c4762a1bSJed Brown f0[1] = 2.0*PETSC_PI*x[1]*PetscCosReal(2.0*PETSC_PI*x[0]); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249c4762a1bSJed Brown static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 250c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 251c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 252c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 253c4762a1bSJed Brown { 254c4762a1bSJed Brown f0[0] = -2.0*PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1]); 255c4762a1bSJed Brown f0[1] = 2.0*PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0]); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 259c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 260c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 261c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 262c4762a1bSJed Brown { 263c4762a1bSJed Brown const PetscInt Ncomp = dim; 264c4762a1bSJed Brown PetscInt compI, d; 265c4762a1bSJed Brown 266c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 267c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 268c4762a1bSJed Brown g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0; 269c4762a1bSJed Brown } 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 274c4762a1bSJed Brown static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 275c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 276c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 277c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 278c4762a1bSJed Brown { 279c4762a1bSJed Brown PetscInt d; 280c4762a1bSJed Brown f0[0] = u_t[dim]; 281c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += u[dim]*u_x[d*dim+d] + u_x[dim*dim+d]*u[d]; 282c4762a1bSJed Brown } 283c4762a1bSJed Brown 284c4762a1bSJed Brown static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 285c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 286c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 287c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 288c4762a1bSJed Brown { 289c4762a1bSJed Brown PetscInt d; 290c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[0] = 0.0; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 294c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 295c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 296c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 297c4762a1bSJed Brown { 298c4762a1bSJed Brown PetscInt d; 299c4762a1bSJed Brown g0[0] = u_tShift; 300c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[d*dim+d]; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303c4762a1bSJed Brown void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 304c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 305c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 306c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 307c4762a1bSJed Brown { 308c4762a1bSJed Brown PetscInt d; 309c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = u[d]; 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 313c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 314c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 315c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 316c4762a1bSJed Brown { 317c4762a1bSJed Brown PetscInt d; 318c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[dim*dim+d]; 319c4762a1bSJed Brown } 320c4762a1bSJed Brown 321c4762a1bSJed Brown void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 322c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 323c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 324c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 325c4762a1bSJed Brown { 326c4762a1bSJed Brown PetscInt d; 327c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = u[dim]; 328c4762a1bSJed Brown } 329c4762a1bSJed Brown 330c4762a1bSJed 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) 331c4762a1bSJed Brown { 332c4762a1bSJed Brown PetscReal wind[3] = {0.0, 1.0, 0.0}; 333c4762a1bSJed Brown PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim,3), wind, n); 334c4762a1bSJed Brown 335c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 336c4762a1bSJed Brown } 337c4762a1bSJed Brown 338c4762a1bSJed 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) 339c4762a1bSJed Brown { 340c4762a1bSJed Brown PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 341c4762a1bSJed Brown 342c4762a1bSJed Brown #if 1 343c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 344c4762a1bSJed Brown #else 345c4762a1bSJed 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]); */ 346c4762a1bSJed Brown /* Smear it out */ 347c4762a1bSJed Brown flux[0] = 0.5*((uL[dim] + uR[dim]) + (uL[dim] - uR[dim])*tanh(1.0e5*wn)) * wn; 348c4762a1bSJed Brown #endif 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 351c4762a1bSJed Brown static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 352c4762a1bSJed Brown { 353c4762a1bSJed Brown u[0] = 0.0; 354c4762a1bSJed Brown u[1] = 0.0; 355c4762a1bSJed Brown return 0; 356c4762a1bSJed Brown } 357c4762a1bSJed Brown 358c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 359c4762a1bSJed Brown { 360c4762a1bSJed Brown u[0] = 0.0; 361c4762a1bSJed Brown u[1] = 1.0; 362c4762a1bSJed Brown return 0; 363c4762a1bSJed Brown } 364c4762a1bSJed Brown 365c4762a1bSJed Brown /* Coordinates of the point which was at x at t = 0 */ 366c4762a1bSJed Brown static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 367c4762a1bSJed Brown { 368c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 369c4762a1bSJed Brown u[0] = x[0]; 370c4762a1bSJed Brown u[1] = x[1] + t; 37130602db0SMatthew G. Knepley #if 0 3729566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u)); 37330602db0SMatthew G. Knepley #else 37430602db0SMatthew G. Knepley u[1] = u[1] - (int) PetscRealPart(u[1]); 37530602db0SMatthew G. Knepley #endif 376c4762a1bSJed Brown return 0; 377c4762a1bSJed Brown } 378c4762a1bSJed Brown 379c4762a1bSJed Brown /* 380c4762a1bSJed Brown In 2D we use the exact solution: 381c4762a1bSJed Brown 382c4762a1bSJed Brown u = x^2 + y^2 383c4762a1bSJed Brown v = 2 x^2 - 2xy 384c4762a1bSJed Brown phi = h(x + y + (u + v) t) 385c4762a1bSJed Brown f_x = f_y = 4 386c4762a1bSJed Brown 387c4762a1bSJed Brown so that 388c4762a1bSJed Brown 389c4762a1bSJed Brown -\Delta u + f = <-4, -4> + <4, 4> = 0 390c4762a1bSJed Brown {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 391c4762a1bSJed Brown h_t(x + y + (u + v) t) - u . grad phi - phi div u 392c4762a1bSJed Brown = u h' + v h' - u h_x - v h_y 393c4762a1bSJed Brown = 0 394c4762a1bSJed Brown 395c4762a1bSJed Brown We will conserve phi since 396c4762a1bSJed Brown 397c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 398c4762a1bSJed Brown 399c4762a1bSJed Brown Also try h((x + ut)^2 + (y + vt)^2), so that 400c4762a1bSJed Brown 401c4762a1bSJed Brown h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 402c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 403c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 404c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 405c4762a1bSJed Brown = 0 406c4762a1bSJed Brown 407c4762a1bSJed Brown */ 408c4762a1bSJed Brown static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 409c4762a1bSJed Brown { 410c4762a1bSJed Brown u[0] = x[0]*x[0] + x[1]*x[1]; 411c4762a1bSJed Brown u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 412c4762a1bSJed Brown return 0; 413c4762a1bSJed Brown } 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* 416c4762a1bSJed Brown In 2D we use the exact, periodic solution: 417c4762a1bSJed Brown 418c4762a1bSJed Brown u = sin(2 pi x)/4 pi^2 419c4762a1bSJed Brown v = -y cos(2 pi x)/2 pi 420c4762a1bSJed Brown phi = h(x + y + (u + v) t) 421c4762a1bSJed Brown f_x = -sin(2 pi x) 422c4762a1bSJed Brown f_y = 2 pi y cos(2 pi x) 423c4762a1bSJed Brown 424c4762a1bSJed Brown so that 425c4762a1bSJed Brown 426c4762a1bSJed Brown -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 427c4762a1bSJed Brown 428c4762a1bSJed Brown We will conserve phi since 429c4762a1bSJed Brown 430c4762a1bSJed Brown \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 431c4762a1bSJed Brown */ 432c4762a1bSJed Brown static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 433c4762a1bSJed Brown { 434c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 435c4762a1bSJed Brown u[1] = -x[1]*PetscCosReal(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI); 436c4762a1bSJed Brown return 0; 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* 440c4762a1bSJed Brown In 2D we use the exact, doubly periodic solution: 441c4762a1bSJed Brown 442c4762a1bSJed Brown u = sin(2 pi x) cos(2 pi y)/4 pi^2 443c4762a1bSJed Brown v = -sin(2 pi y) cos(2 pi x)/4 pi^2 444c4762a1bSJed Brown phi = h(x + y + (u + v) t) 445c4762a1bSJed Brown f_x = -2sin(2 pi x) cos(2 pi y) 446c4762a1bSJed Brown f_y = 2sin(2 pi y) cos(2 pi x) 447c4762a1bSJed Brown 448c4762a1bSJed Brown so that 449c4762a1bSJed Brown 450c4762a1bSJed 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 451c4762a1bSJed Brown 452c4762a1bSJed Brown We will conserve phi since 453c4762a1bSJed Brown 454c4762a1bSJed Brown \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 455c4762a1bSJed Brown */ 456c4762a1bSJed Brown static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 457c4762a1bSJed Brown { 458c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1])/PetscSqr(2.0*PETSC_PI); 459c4762a1bSJed Brown u[1] = -PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 460c4762a1bSJed Brown return 0; 461c4762a1bSJed Brown } 462c4762a1bSJed Brown 463c4762a1bSJed Brown static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 464c4762a1bSJed Brown { 465c4762a1bSJed Brown u[0] = x[1] - 0.5; 466c4762a1bSJed Brown u[1] = 0.0; 467c4762a1bSJed Brown return 0; 468c4762a1bSJed Brown } 469c4762a1bSJed Brown 470c4762a1bSJed Brown static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 471c4762a1bSJed Brown { 472c4762a1bSJed Brown PetscInt d; 473c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 474c4762a1bSJed Brown return 0; 475c4762a1bSJed Brown } 476c4762a1bSJed Brown 477c4762a1bSJed Brown static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 478c4762a1bSJed Brown { 479c4762a1bSJed Brown u[0] = 0.0; 480c4762a1bSJed Brown return 0; 481c4762a1bSJed Brown } 482c4762a1bSJed Brown 483c4762a1bSJed Brown static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 484c4762a1bSJed Brown { 485c4762a1bSJed Brown u[0] = 1.0; 486c4762a1bSJed Brown return 0; 487c4762a1bSJed Brown } 488c4762a1bSJed Brown 489c4762a1bSJed Brown static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 490c4762a1bSJed Brown { 491c4762a1bSJed Brown PetscReal x0[2]; 492c4762a1bSJed Brown PetscScalar xn[2]; 493c4762a1bSJed Brown 494c4762a1bSJed Brown x0[0] = globalUser->source[0]; 495c4762a1bSJed Brown x0[1] = globalUser->source[1]; 496c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 497c4762a1bSJed Brown { 498c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 499c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 500c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 501c4762a1bSJed Brown 502c4762a1bSJed Brown u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 503c4762a1bSJed Brown } 504c4762a1bSJed Brown return 0; 505c4762a1bSJed Brown } 506c4762a1bSJed Brown 507c4762a1bSJed Brown /* 508c4762a1bSJed Brown Gaussian blob, initially centered on (0.5, 0.5) 509c4762a1bSJed Brown 510c4762a1bSJed Brown xi = x(t) - x0, eta = y(t) - y0 511c4762a1bSJed Brown 512c4762a1bSJed Brown where x(t), y(t) are the integral curves of v(t), 513c4762a1bSJed Brown 514c4762a1bSJed Brown dx/dt . grad f = v . f 515c4762a1bSJed Brown 516c4762a1bSJed Brown Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 517c4762a1bSJed Brown 518c4762a1bSJed Brown v0 f_x + w0 f_y = v . f 519c4762a1bSJed Brown */ 520c4762a1bSJed Brown static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 521c4762a1bSJed Brown { 522c4762a1bSJed Brown const PetscReal x0[2] = {0.5, 0.5}; 523c4762a1bSJed Brown const PetscReal sigma = 1.0/6.0; 524c4762a1bSJed Brown PetscScalar xn[2]; 525c4762a1bSJed Brown 526c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 527c4762a1bSJed Brown { 528c4762a1bSJed Brown /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 529c4762a1bSJed Brown /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 530c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 531c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 532c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 533c4762a1bSJed Brown 534c4762a1bSJed Brown u[0] = PetscExpReal(-r2/(2.0*sigma*sigma))/(sigma*PetscSqrtReal(2.0*PETSC_PI)); 535c4762a1bSJed Brown } 536c4762a1bSJed Brown return 0; 537c4762a1bSJed Brown } 538c4762a1bSJed Brown 539c4762a1bSJed Brown static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 540c4762a1bSJed Brown { 541c4762a1bSJed Brown PetscReal x0[3]; 542c4762a1bSJed Brown const PetscReal wind[3] = {0.0, 1.0, 0.0}; 543c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 544c4762a1bSJed Brown 545c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 546c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 547c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 548c4762a1bSJed Brown return 0; 549c4762a1bSJed Brown } 550c4762a1bSJed Brown 551c4762a1bSJed Brown static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 552c4762a1bSJed Brown { 553c4762a1bSJed Brown PetscReal ur[3]; 554c4762a1bSJed Brown PetscReal x0[3]; 555c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 556c4762a1bSJed Brown 557c4762a1bSJed Brown ur[0] = PetscRealPart(u[0]); ur[1] = PetscRealPart(u[1]); ur[2] = PetscRealPart(u[2]); 558c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 559c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 560c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 561c4762a1bSJed Brown return 0; 562c4762a1bSJed Brown } 563c4762a1bSJed Brown 564c4762a1bSJed Brown static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 565c4762a1bSJed Brown { 566c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 567c4762a1bSJed Brown 568c4762a1bSJed Brown PetscFunctionBeginUser; 569c4762a1bSJed Brown xG[0] = user->inflowState; 570c4762a1bSJed Brown PetscFunctionReturn(0); 571c4762a1bSJed Brown } 572c4762a1bSJed Brown 573c4762a1bSJed Brown static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 574c4762a1bSJed Brown { 575c4762a1bSJed Brown PetscFunctionBeginUser; 57630602db0SMatthew G. Knepley //xG[0] = xI[dim]; 57730602db0SMatthew G. Knepley xG[0] = xI[2]; 578c4762a1bSJed Brown PetscFunctionReturn(0); 579c4762a1bSJed Brown } 580c4762a1bSJed Brown 581c4762a1bSJed Brown static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 582c4762a1bSJed Brown { 583c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 584c4762a1bSJed Brown PetscInt dim; 585c4762a1bSJed Brown 586c4762a1bSJed Brown PetscFunctionBeginUser; 5879566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 588c4762a1bSJed Brown switch (user->porosityDist) { 589c4762a1bSJed Brown case TILTED: 590c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) tilted_phi_2d(dim, time, x, 2, u, (void *) &time); 591c4762a1bSJed Brown else tilted_phi_coupled_2d(dim, time, x, 2, u, (void *) &time); 592c4762a1bSJed Brown break; 593c4762a1bSJed Brown case GAUSSIAN: 594c4762a1bSJed Brown gaussian_phi_2d(dim, time, x, 2, u, (void *) &time); 595c4762a1bSJed Brown break; 596c4762a1bSJed Brown case DELTA: 597c4762a1bSJed Brown delta_phi_2d(dim, time, x, 2, u, (void *) &time); 598c4762a1bSJed Brown break; 599c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 600c4762a1bSJed Brown } 601c4762a1bSJed Brown PetscFunctionReturn(0); 602c4762a1bSJed Brown } 603c4762a1bSJed Brown 604c4762a1bSJed Brown static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 605c4762a1bSJed Brown { 606c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 607c4762a1bSJed Brown PetscScalar yexact[3]={0,0,0}; 608c4762a1bSJed Brown 609c4762a1bSJed Brown PetscFunctionBeginUser; 6109566063dSJacob Faibussowitsch PetscCall(ExactSolution(dm, time, x, yexact, ctx)); 611c4762a1bSJed Brown f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 612c4762a1bSJed Brown PetscFunctionReturn(0); 613c4762a1bSJed Brown } 614c4762a1bSJed Brown 615c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 616c4762a1bSJed Brown { 617c4762a1bSJed Brown PetscFunctionBeginUser; 6189566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 6199566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 62030602db0SMatthew G. Knepley #if 0 62130602db0SMatthew 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; 62230602db0SMatthew G. Knepley const PetscReal L[3] = {1.0, 1.0, 1.0}; 62330602db0SMatthew G. Knepley PetscReal maxCell[3]; 62430602db0SMatthew G. Knepley PetscInt d; 62530602db0SMatthew G. Knepley 6269566063dSJacob Faibussowitsch if (periodic) {for (d = 0; d < 3; ++d) maxCell[d] = 1.1*(L[d]/cells[d]); PetscCall(DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, user->bd));} 62730602db0SMatthew G. Knepley #endif 6289566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 6299566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 630c4762a1bSJed Brown PetscFunctionReturn(0); 631c4762a1bSJed Brown } 632c4762a1bSJed Brown 633c4762a1bSJed Brown static PetscErrorCode SetupBC(DM dm, AppCtx *user) 634c4762a1bSJed Brown { 635348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 63630602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 637c4762a1bSJed Brown PetscDS prob; 638c4762a1bSJed Brown DMLabel label; 639c4762a1bSJed Brown PetscBool check; 64030602db0SMatthew G. Knepley PetscInt dim, n = 3; 64130602db0SMatthew G. Knepley const char *prefix; 642c4762a1bSJed Brown 643c4762a1bSJed Brown PetscFunctionBeginUser; 6449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 6459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 6469566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 647c4762a1bSJed Brown /* Set initial guesses and exact solutions */ 64830602db0SMatthew G. Knepley switch (dim) { 649c4762a1bSJed Brown case 2: 650c4762a1bSJed Brown user->initialGuess[0] = initialVelocity; 651c4762a1bSJed Brown switch(user->porosityDist) { 652c4762a1bSJed Brown case ZERO: user->initialGuess[1] = zero_phi;break; 653c4762a1bSJed Brown case CONSTANT: user->initialGuess[1] = constant_phi;break; 654c4762a1bSJed Brown case GAUSSIAN: user->initialGuess[1] = gaussian_phi_2d;break; 655c4762a1bSJed Brown case DELTA: user->initialGuess[1] = delta_phi_2d;break; 656c4762a1bSJed Brown case TILTED: 657c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 658c4762a1bSJed Brown else user->initialGuess[1] = tilted_phi_coupled_2d; 659c4762a1bSJed Brown break; 660c4762a1bSJed Brown } 661348a1646SMatthew G. Knepley break; 66298921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %D not supported", dim); 663348a1646SMatthew G. Knepley } 664348a1646SMatthew G. Knepley exactFuncs[0] = user->initialGuess[0]; 665348a1646SMatthew G. Knepley exactFuncs[1] = user->initialGuess[1]; 66630602db0SMatthew G. Knepley switch (dim) { 667348a1646SMatthew G. Knepley case 2: 668c4762a1bSJed Brown switch (user->velocityDist) { 669c4762a1bSJed Brown case VEL_ZERO: 670348a1646SMatthew G. Knepley exactFuncs[0] = zero_u_2d; break; 671c4762a1bSJed Brown case VEL_CONSTANT: 672348a1646SMatthew G. Knepley exactFuncs[0] = constant_u_2d; break; 673c4762a1bSJed Brown case VEL_HARMONIC: 67430602db0SMatthew G. Knepley switch (bdt[0]) { 675c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 67630602db0SMatthew G. Knepley switch (bdt[1]) { 677c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 678348a1646SMatthew G. Knepley exactFuncs[0] = doubly_periodic_u_2d; break; 679c4762a1bSJed Brown default: 680348a1646SMatthew G. Knepley exactFuncs[0] = periodic_u_2d; break; 681c4762a1bSJed Brown } 682c4762a1bSJed Brown break; 683c4762a1bSJed Brown default: 684348a1646SMatthew G. Knepley exactFuncs[0] = quadratic_u_2d; break; 685c4762a1bSJed Brown } 686c4762a1bSJed Brown break; 687c4762a1bSJed Brown case VEL_SHEAR: 688348a1646SMatthew G. Knepley exactFuncs[0] = shear_bc; break; 68998921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim); 690c4762a1bSJed Brown } 691348a1646SMatthew G. Knepley break; 69298921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %D not supported", dim); 693c4762a1bSJed Brown } 694c4762a1bSJed Brown { 695c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 696c4762a1bSJed Brown 6979566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 698348a1646SMatthew G. Knepley if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 699c4762a1bSJed Brown } 7009566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-dmts_check", &check)); 701c4762a1bSJed Brown if (check) { 702348a1646SMatthew G. Knepley user->initialGuess[0] = exactFuncs[0]; 703348a1646SMatthew G. Knepley user->initialGuess[1] = exactFuncs[1]; 704c4762a1bSJed Brown } 705c4762a1bSJed Brown /* Set BC */ 7069566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 7079566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 7089566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user)); 7099566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user)); 710c4762a1bSJed Brown if (label) { 711c4762a1bSJed Brown const PetscInt id = 1; 712c4762a1bSJed Brown 7139566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL)); 714c4762a1bSJed Brown } 7159566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &label)); 716c4762a1bSJed Brown if (label && user->useFV) { 717c4762a1bSJed Brown const PetscInt inflowids[] = {100,200,300}, outflowids[] = {101}; 718c4762a1bSJed Brown 719*dd39110bSPierre 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)); 720*dd39110bSPierre 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)); 721c4762a1bSJed Brown } 722c4762a1bSJed Brown PetscFunctionReturn(0); 723c4762a1bSJed Brown } 724c4762a1bSJed Brown 725c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 726c4762a1bSJed Brown { 72730602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 728c4762a1bSJed Brown PetscDS prob; 72930602db0SMatthew G. Knepley PetscInt n = 3; 73030602db0SMatthew G. Knepley const char *prefix; 731c4762a1bSJed Brown 732c4762a1bSJed Brown PetscFunctionBeginUser; 7339566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 7349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 7359566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 736c4762a1bSJed Brown switch (user->velocityDist) { 737c4762a1bSJed Brown case VEL_ZERO: 7389566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u)); 739c4762a1bSJed Brown break; 740c4762a1bSJed Brown case VEL_CONSTANT: 7419566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u)); 7429566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL)); 7439566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL)); 744c4762a1bSJed Brown break; 745c4762a1bSJed Brown case VEL_HARMONIC: 74630602db0SMatthew G. Knepley switch (bdt[0]) { 747c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 74830602db0SMatthew G. Knepley switch (bdt[1]) { 749c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 7509566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u)); 751c4762a1bSJed Brown break; 752c4762a1bSJed Brown default: 7539566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u)); 754c4762a1bSJed Brown break; 755c4762a1bSJed Brown } 756c4762a1bSJed Brown break; 757c4762a1bSJed Brown default: 7589566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u)); 759c4762a1bSJed Brown break; 760c4762a1bSJed Brown } 7619566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 762c4762a1bSJed Brown break; 763c4762a1bSJed Brown case VEL_SHEAR: 7649566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u)); 7659566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 766c4762a1bSJed Brown break; 767c4762a1bSJed Brown } 7689566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_advection, f1_advection)); 7699566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL)); 7709566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL)); 7719566063dSJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_advection)); 7729566063dSJacob Faibussowitsch else PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection)); 773c4762a1bSJed Brown 7749566063dSJacob Faibussowitsch PetscCall(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user)); 775c4762a1bSJed Brown PetscFunctionReturn(0); 776c4762a1bSJed Brown } 777c4762a1bSJed Brown 778c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 779c4762a1bSJed Brown { 780c4762a1bSJed Brown DM cdm = dm; 781c4762a1bSJed Brown PetscQuadrature q; 782c4762a1bSJed Brown PetscFE fe[2]; 783c4762a1bSJed Brown PetscFV fv; 784c4762a1bSJed Brown MPI_Comm comm; 78530602db0SMatthew G. Knepley PetscInt dim; 786c4762a1bSJed Brown 787c4762a1bSJed Brown PetscFunctionBeginUser; 788c4762a1bSJed Brown /* Create finite element */ 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 7909566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7919566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0])); 7929566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity")); 7939566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1])); 7949566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 7959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[1], "porosity")); 796c4762a1bSJed Brown 7979566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(PetscObjectComm((PetscObject) dm), &fv)); 7989566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fv, "porosity")); 7999566063dSJacob Faibussowitsch PetscCall(PetscFVSetFromOptions(fv)); 8009566063dSJacob Faibussowitsch PetscCall(PetscFVSetNumComponents(fv, 1)); 8019566063dSJacob Faibussowitsch PetscCall(PetscFVSetSpatialDimension(fv, dim)); 8029566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe[0], &q)); 8039566063dSJacob Faibussowitsch PetscCall(PetscFVSetQuadrature(fv, q)); 804c4762a1bSJed Brown 8059566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 8069566063dSJacob Faibussowitsch if (user->useFV) PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fv)); 8079566063dSJacob Faibussowitsch else PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 8089566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 8099566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 810c4762a1bSJed Brown 811c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 812c4762a1bSJed Brown while (cdm) { 8139566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 8149566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 815c4762a1bSJed Brown /* Coordinates were never localized for coarse meshes */ 8169566063dSJacob Faibussowitsch if (cdm) PetscCall(DMLocalizeCoordinates(cdm)); 817c4762a1bSJed Brown } 8189566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 8199566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 8209566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fv)); 821c4762a1bSJed Brown PetscFunctionReturn(0); 822c4762a1bSJed Brown } 823c4762a1bSJed Brown 824c4762a1bSJed Brown static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 825c4762a1bSJed Brown { 826c4762a1bSJed Brown PetscFunctionBeginUser; 8279566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, user, dm)); 828c4762a1bSJed Brown /* Handle refinement, etc. */ 8299566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 830c4762a1bSJed Brown /* Construct ghost cells */ 831c4762a1bSJed Brown if (user->useFV) { 832c4762a1bSJed Brown DM gdm; 833c4762a1bSJed Brown 8349566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm)); 8359566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 836c4762a1bSJed Brown *dm = gdm; 837c4762a1bSJed Brown } 838c4762a1bSJed Brown /* Localize coordinates */ 8399566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 8409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)(*dm),"Mesh")); 8419566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 842c4762a1bSJed Brown /* Setup problem */ 8439566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*dm, user)); 844c4762a1bSJed Brown /* Setup BC */ 8459566063dSJacob Faibussowitsch PetscCall(SetupBC(*dm, user)); 846c4762a1bSJed Brown PetscFunctionReturn(0); 847c4762a1bSJed Brown } 848c4762a1bSJed Brown 849c4762a1bSJed Brown static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx) 850c4762a1bSJed Brown { 851c4762a1bSJed Brown PetscDS prob; 852c4762a1bSJed Brown DM dmCell; 853c4762a1bSJed Brown Vec cellgeom; 854c4762a1bSJed Brown const PetscScalar *cgeom; 855c4762a1bSJed Brown PetscScalar *x; 856c4762a1bSJed Brown PetscInt dim, Nf, cStart, cEnd, c; 857c4762a1bSJed Brown 858c4762a1bSJed Brown PetscFunctionBeginUser; 8599566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 8609566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8619566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 8629566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8639566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 8649566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 8659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 8669566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 867c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 868c4762a1bSJed Brown PetscFVCellGeom *cg; 869c4762a1bSJed Brown PetscScalar *xc; 870c4762a1bSJed Brown 8719566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 8729566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc)); 8739566063dSJacob Faibussowitsch if (xc) PetscCall((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx)); 874c4762a1bSJed Brown } 8759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 8769566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 877c4762a1bSJed Brown PetscFunctionReturn(0); 878c4762a1bSJed Brown } 879c4762a1bSJed Brown 880c4762a1bSJed Brown static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 881c4762a1bSJed Brown { 882c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 883c4762a1bSJed Brown char *ftable = NULL; 884c4762a1bSJed Brown DM dm; 885c4762a1bSJed Brown PetscSection s; 886c4762a1bSJed Brown Vec cellgeom; 887c4762a1bSJed Brown const PetscScalar *x; 888c4762a1bSJed Brown PetscScalar *a; 889c4762a1bSJed Brown PetscReal *xnorms; 890c4762a1bSJed Brown PetscInt pStart, pEnd, p, Nf, f; 891c4762a1bSJed Brown 892c4762a1bSJed Brown PetscFunctionBeginUser; 8939566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X, (PetscObject) ts, "-view_solution")); 8949566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 8959566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8969566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 8979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &Nf)); 8989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 8999566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nf*2, &xnorms)); 9009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 901c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 902c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 903c4762a1bSJed Brown PetscInt dof, cdof, d; 904c4762a1bSJed Brown 9059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 9069566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 9079566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, x, &a)); 908c4762a1bSJed Brown /* TODO Use constrained indices here */ 909c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+0] = PetscMax(xnorms[f*2+0], PetscAbsScalar(a[d])); 910c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+1] += PetscAbsScalar(a[d]); 911c4762a1bSJed Brown } 912c4762a1bSJed Brown } 9139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 914c4762a1bSJed Brown if (stepnum >= 0) { /* No summary for final time */ 915c4762a1bSJed Brown DM dmCell, *fdm; 916c4762a1bSJed Brown Vec *fv; 917c4762a1bSJed Brown const PetscScalar *cgeom; 918c4762a1bSJed Brown PetscScalar **fx; 919c4762a1bSJed Brown PetscReal *fmin, *fmax, *fint, *ftmp, t; 920c4762a1bSJed Brown PetscInt cStart, cEnd, c, fcount, f, num; 921c4762a1bSJed Brown 922c4762a1bSJed Brown size_t ftableused,ftablealloc; 923c4762a1bSJed Brown 924c4762a1bSJed Brown /* Functionals have indices after registering, this is an upper bound */ 925c4762a1bSJed Brown fcount = user->numMonitorFuncs; 9269566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fint,fcount,&ftmp)); 9279566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fcount,&fdm,fcount,&fv,fcount,&fx)); 928c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 929c4762a1bSJed Brown PetscSection fs; 930c4762a1bSJed Brown const char *name = user->monitorFuncs[f]->name; 931c4762a1bSJed Brown 932c4762a1bSJed Brown fmin[f] = PETSC_MAX_REAL; 933c4762a1bSJed Brown fmax[f] = PETSC_MIN_REAL; 934c4762a1bSJed Brown fint[f] = 0; 935c4762a1bSJed Brown /* Make monitor vecs */ 9369566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &fdm[f])); 9379566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &num, &t)); 9389566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(fdm[f], num, t)); 9399566063dSJacob Faibussowitsch PetscCall(PetscSectionClone(s, &fs)); 9409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 0, NULL)); 9419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 1, name)); 9429566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(fdm[f], fs)); 9439566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&fs)); 9449566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm[f], &fv[f])); 9459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fv[f], name)); 9469566063dSJacob Faibussowitsch PetscCall(VecGetArray(fv[f], &fx[f])); 947c4762a1bSJed Brown } 9489566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 9499566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 9509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 9519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 952c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 953c4762a1bSJed Brown PetscFVCellGeom *cg; 954c4762a1bSJed Brown PetscScalar *cx; 955c4762a1bSJed Brown 9569566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 9579566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx)); 958c4762a1bSJed Brown if (!cx) continue; /* not a global cell */ 959c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 960c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 961c4762a1bSJed Brown PetscScalar *fxc; 962c4762a1bSJed Brown 9639566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc)); 964c4762a1bSJed Brown /* I need to make it easier to get interpolated values here */ 9659566063dSJacob Faibussowitsch PetscCall((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx)); 966c4762a1bSJed Brown fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 967c4762a1bSJed Brown } 968c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 969c4762a1bSJed Brown fmin[f] = PetscMin(fmin[f], ftmp[f]); 970c4762a1bSJed Brown fmax[f] = PetscMax(fmax[f], ftmp[f]); 971c4762a1bSJed Brown fint[f] += cg->volume * ftmp[f]; 972c4762a1bSJed Brown } 973c4762a1bSJed Brown } 9749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 9759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 9769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 9779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 9789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 979c4762a1bSJed Brown /* Output functional data */ 980c4762a1bSJed Brown ftablealloc = fcount * 100; 981c4762a1bSJed Brown ftableused = 0; 9829566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ftablealloc, &ftable)); 983c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 984c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 985c4762a1bSJed Brown PetscInt id = func->offset; 986c4762a1bSJed Brown char newline[] = "\n"; 987c4762a1bSJed Brown char buffer[256], *p, *prefix; 988c4762a1bSJed Brown size_t countused, len; 989c4762a1bSJed Brown 990c4762a1bSJed Brown /* Create string with functional outputs */ 991c4762a1bSJed Brown if (f % 3) { 9929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, " ", 2)); 993c4762a1bSJed Brown p = buffer + 2; 994c4762a1bSJed Brown } else if (f) { 9959566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, newline, sizeof(newline)-1)); 996c4762a1bSJed Brown p = buffer + sizeof(newline) - 1; 997c4762a1bSJed Brown } else { 998c4762a1bSJed Brown p = buffer; 999c4762a1bSJed Brown } 10009566063dSJacob 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])); 1001c4762a1bSJed Brown countused += p - buffer; 1002c4762a1bSJed Brown /* reallocate */ 1003c4762a1bSJed Brown if (countused > ftablealloc-ftableused-1) { 1004c4762a1bSJed Brown char *ftablenew; 1005c4762a1bSJed Brown 1006c4762a1bSJed Brown ftablealloc = 2*ftablealloc + countused; 10079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ftablealloc, &ftablenew)); 10089566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); 10099566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 1010c4762a1bSJed Brown ftable = ftablenew; 1011c4762a1bSJed Brown } 10129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftable+ftableused, buffer, countused)); 1013c4762a1bSJed Brown ftableused += countused; 1014c4762a1bSJed Brown ftable[ftableused] = 0; 1015c4762a1bSJed Brown /* Output vecs */ 10169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fv[f], &fx[f])); 10179566063dSJacob Faibussowitsch PetscCall(PetscStrlen(func->name, &len)); 10189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len+2,&prefix)); 10199566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(prefix, func->name)); 10209566063dSJacob Faibussowitsch PetscCall(PetscStrcat(prefix, "_")); 10219566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix)); 10229566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fv[f], NULL, "-vec_view")); 10239566063dSJacob Faibussowitsch PetscCall(PetscFree(prefix)); 10249566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm[f], &fv[f])); 10259566063dSJacob Faibussowitsch PetscCall(DMDestroy(&fdm[f])); 1026c4762a1bSJed Brown } 10279566063dSJacob Faibussowitsch PetscCall(PetscFree4(fmin, fmax, fint, ftmp)); 10289566063dSJacob Faibussowitsch PetscCall(PetscFree3(fdm, fv, fx)); 10299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "% 3D time %8.4g |x| (", stepnum, (double) time)); 1030c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10319566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 10329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+0])); 1033c4762a1bSJed Brown } 10349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ") |x|_1 (")); 1035c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10369566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 10379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+1])); 1038c4762a1bSJed Brown } 10399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ") %s\n", ftable ? ftable : "")); 10409566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 1041c4762a1bSJed Brown } 10429566063dSJacob Faibussowitsch PetscCall(PetscFree(xnorms)); 1043c4762a1bSJed Brown PetscFunctionReturn(0); 1044c4762a1bSJed Brown } 1045c4762a1bSJed Brown 1046c4762a1bSJed Brown int main(int argc, char **argv) 1047c4762a1bSJed Brown { 1048c4762a1bSJed Brown MPI_Comm comm; 1049c4762a1bSJed Brown TS ts; 1050c4762a1bSJed Brown DM dm; 1051c4762a1bSJed Brown Vec u; 1052c4762a1bSJed Brown AppCtx user; 1053c4762a1bSJed Brown PetscReal t0, t = 0.0; 1054c4762a1bSJed Brown void *ctxs[2]; 1055c4762a1bSJed Brown 1056c4762a1bSJed Brown ctxs[0] = &t; 1057c4762a1bSJed Brown ctxs[1] = &t; 10589566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char*) 0, help)); 1059c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1060c4762a1bSJed Brown user.functionalRegistry = NULL; 1061c4762a1bSJed Brown globalUser = &user; 10629566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 10639566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 10649566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 10659566063dSJacob Faibussowitsch PetscCall(CreateDM(comm, &user, &dm)); 10669566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 10679566063dSJacob Faibussowitsch PetscCall(ProcessMonitorOptions(comm, &user)); 1068c4762a1bSJed Brown 10699566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 10709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "solution")); 1071c4762a1bSJed Brown if (user.useFV) { 1072c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 1073c4762a1bSJed Brown 10749566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 1075c4762a1bSJed Brown if (isImplicit) { 10769566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10779566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1078c4762a1bSJed Brown } 10799566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10809566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user)); 1081c4762a1bSJed Brown } else { 10829566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10839566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10849566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1085c4762a1bSJed Brown } 10869566063dSJacob Faibussowitsch if (user.useFV) PetscCall(TSMonitorSet(ts, MonitorFunctionals, &user, NULL)); 10879566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1)); 10889566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2.0)); 10899566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,0.01)); 10909566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 10919566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1092c4762a1bSJed Brown 10939566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u)); 10949566063dSJacob Faibussowitsch if (user.useFV) PetscCall(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1])); 10959566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-init_vec_view")); 10969566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1097c4762a1bSJed Brown t0 = t; 10989566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 10999566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 11009566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 11019566063dSJacob Faibussowitsch if (t > t0) PetscCall(DMTSCheckFromOptions(ts, u)); 11029566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 1103c4762a1bSJed Brown { 1104c4762a1bSJed Brown PetscReal ftime; 1105c4762a1bSJed Brown PetscInt nsteps; 1106c4762a1bSJed Brown TSConvergedReason reason; 1107c4762a1bSJed Brown 11089566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 11099566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &nsteps)); 11109566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 11119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %D steps\n", TSConvergedReasons[reason], (double) ftime, nsteps)); 1112c4762a1bSJed Brown } 1113c4762a1bSJed Brown 11149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 11159566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 11169566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 11179566063dSJacob Faibussowitsch PetscCall(PetscFree(user.monitorFuncs)); 11189566063dSJacob Faibussowitsch PetscCall(FunctionalDestroy(&user.functionalRegistry)); 11199566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1120b122ec5aSJacob Faibussowitsch return 0; 1121c4762a1bSJed Brown } 1122c4762a1bSJed Brown 1123c4762a1bSJed Brown /*TEST 1124c4762a1bSJed Brown 112530602db0SMatthew G. Knepley testset: 112630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 112730602db0SMatthew G. Knepley 1128c4762a1bSJed Brown # 2D harmonic velocity, no porosity 1129c4762a1bSJed Brown test: 1130c4762a1bSJed Brown suffix: p1p1 1131c4762a1bSJed Brown requires: !complex !single 113230602db0SMatthew 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 1133c4762a1bSJed Brown test: 1134c4762a1bSJed Brown suffix: p1p1_xper 1135c4762a1bSJed Brown requires: !complex !single 113630602db0SMatthew 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 1137c4762a1bSJed Brown test: 1138c4762a1bSJed Brown suffix: p1p1_xper_ref 1139c4762a1bSJed Brown requires: !complex !single 114030602db0SMatthew 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 1141c4762a1bSJed Brown test: 1142c4762a1bSJed Brown suffix: p1p1_xyper 1143c4762a1bSJed Brown requires: !complex !single 114430602db0SMatthew 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 1145c4762a1bSJed Brown test: 1146c4762a1bSJed Brown suffix: p1p1_xyper_ref 1147c4762a1bSJed Brown requires: !complex !single 114830602db0SMatthew 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 1149c4762a1bSJed Brown test: 1150c4762a1bSJed Brown suffix: p2p1 1151c4762a1bSJed Brown requires: !complex !single 115230602db0SMatthew 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 1153c4762a1bSJed Brown test: 1154c4762a1bSJed Brown suffix: p2p1_xyper 1155c4762a1bSJed Brown requires: !complex !single 115630602db0SMatthew 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 115730602db0SMatthew G. Knepley 115830602db0SMatthew G. Knepley test: 115930602db0SMatthew G. Knepley suffix: adv_1 116030602db0SMatthew G. Knepley requires: !complex !single 116130602db0SMatthew 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 116230602db0SMatthew G. Knepley 116330602db0SMatthew G. Knepley test: 116430602db0SMatthew G. Knepley suffix: adv_2 116530602db0SMatthew G. Knepley requires: !complex 116630602db0SMatthew G. Knepley TODO: broken memory corruption 116730602db0SMatthew 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 116830602db0SMatthew G. Knepley 116930602db0SMatthew G. Knepley test: 117030602db0SMatthew G. Knepley suffix: adv_3 117130602db0SMatthew G. Knepley requires: !complex 117230602db0SMatthew G. Knepley TODO: broken memory corruption 117330602db0SMatthew 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 117430602db0SMatthew G. Knepley 117530602db0SMatthew G. Knepley test: 117630602db0SMatthew G. Knepley suffix: adv_3_ex 117730602db0SMatthew G. Knepley requires: !complex 117830602db0SMatthew 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 117930602db0SMatthew G. Knepley 118030602db0SMatthew G. Knepley test: 118130602db0SMatthew G. Knepley suffix: adv_4 118230602db0SMatthew G. Knepley requires: !complex 118330602db0SMatthew G. Knepley TODO: broken memory corruption 118430602db0SMatthew 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 118530602db0SMatthew G. Knepley 118630602db0SMatthew G. Knepley # 2D Advection, box, delta 118730602db0SMatthew G. Knepley test: 118830602db0SMatthew G. Knepley suffix: adv_delta_yper_0 118930602db0SMatthew G. Knepley requires: !complex 119030602db0SMatthew G. Knepley TODO: broken 119130602db0SMatthew 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 119230602db0SMatthew G. Knepley 119330602db0SMatthew G. Knepley test: 119430602db0SMatthew G. Knepley suffix: adv_delta_yper_1 119530602db0SMatthew G. Knepley requires: !complex 119630602db0SMatthew G. Knepley TODO: broken 119730602db0SMatthew 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 119830602db0SMatthew G. Knepley 119930602db0SMatthew G. Knepley test: 120030602db0SMatthew G. Knepley suffix: adv_delta_yper_2 120130602db0SMatthew G. Knepley requires: !complex 120230602db0SMatthew G. Knepley TODO: broken 120330602db0SMatthew 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 120430602db0SMatthew G. Knepley 120530602db0SMatthew G. Knepley test: 120630602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_0 120730602db0SMatthew G. Knepley requires: !complex 120830602db0SMatthew G. Knepley TODO: broken 120930602db0SMatthew G. Knepley args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_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 121030602db0SMatthew G. Knepley 121130602db0SMatthew G. Knepley test: 121230602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_1 121330602db0SMatthew G. Knepley requires: !complex 121430602db0SMatthew G. Knepley TODO: broken 121530602db0SMatthew G. Knepley args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 1 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_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 121630602db0SMatthew G. Knepley 121730602db0SMatthew G. Knepley test: 121830602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_2 121930602db0SMatthew G. Knepley requires: !complex 122030602db0SMatthew G. Knepley TODO: broken 122130602db0SMatthew 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 122230602db0SMatthew G. Knepley 122330602db0SMatthew G. Knepley test: 122430602db0SMatthew G. Knepley suffix: adv_delta_yper_im_0 122530602db0SMatthew G. Knepley requires: !complex 122630602db0SMatthew G. Knepley TODO: broken 122730602db0SMatthew 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 122830602db0SMatthew G. Knepley 122930602db0SMatthew G. Knepley test: 123030602db0SMatthew G. Knepley suffix: adv_delta_yper_im_1 123130602db0SMatthew G. Knepley requires: !complex 123230602db0SMatthew G. Knepley TODO: broken 123330602db0SMatthew 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 123430602db0SMatthew G. Knepley 123530602db0SMatthew G. Knepley test: 123630602db0SMatthew G. Knepley suffix: adv_delta_yper_im_2 123730602db0SMatthew G. Knepley requires: !complex 123830602db0SMatthew G. Knepley TODO: broken 123930602db0SMatthew 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 124030602db0SMatthew G. Knepley 124130602db0SMatthew G. Knepley test: 124230602db0SMatthew G. Knepley suffix: adv_delta_yper_im_3 124330602db0SMatthew G. Knepley requires: !complex 124430602db0SMatthew G. Knepley TODO: broken 124530602db0SMatthew 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 124630602db0SMatthew G. Knepley 124730602db0SMatthew G. Knepley # I believe the nullspace is sin(pi y) 124830602db0SMatthew G. Knepley test: 124930602db0SMatthew G. Knepley suffix: adv_delta_yper_im_4 125030602db0SMatthew G. Knepley requires: !complex 125130602db0SMatthew G. Knepley TODO: broken 125230602db0SMatthew 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 125330602db0SMatthew G. Knepley 125430602db0SMatthew G. Knepley test: 125530602db0SMatthew G. Knepley suffix: adv_delta_yper_im_5 125630602db0SMatthew G. Knepley requires: !complex 125730602db0SMatthew G. Knepley TODO: broken 125830602db0SMatthew 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 125930602db0SMatthew G. Knepley 126030602db0SMatthew G. Knepley test: 126130602db0SMatthew G. Knepley suffix: adv_delta_yper_im_6 126230602db0SMatthew G. Knepley requires: !complex 126330602db0SMatthew G. Knepley TODO: broken 126430602db0SMatthew 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 126530602db0SMatthew G. Knepley # 2D Advection, magma benchmark 1 126630602db0SMatthew G. Knepley 126730602db0SMatthew G. Knepley test: 126830602db0SMatthew G. Knepley suffix: adv_delta_shear_im_0 126930602db0SMatthew G. Knepley requires: !complex 127030602db0SMatthew G. Knepley TODO: broken 127130602db0SMatthew 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 127230602db0SMatthew G. Knepley # 2D Advection, box, gaussian 127330602db0SMatthew G. Knepley 127430602db0SMatthew G. Knepley test: 127530602db0SMatthew G. Knepley suffix: adv_gauss 127630602db0SMatthew G. Knepley requires: !complex 127730602db0SMatthew G. Knepley TODO: broken 127830602db0SMatthew 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 127930602db0SMatthew G. Knepley 128030602db0SMatthew G. Knepley test: 128130602db0SMatthew G. Knepley suffix: adv_gauss_im 128230602db0SMatthew G. Knepley requires: !complex 128330602db0SMatthew G. Knepley TODO: broken 128430602db0SMatthew 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 128530602db0SMatthew G. Knepley 128630602db0SMatthew G. Knepley test: 128730602db0SMatthew G. Knepley suffix: adv_gauss_im_1 128830602db0SMatthew G. Knepley requires: !complex 128930602db0SMatthew G. Knepley TODO: broken 129030602db0SMatthew 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 129130602db0SMatthew G. Knepley 129230602db0SMatthew G. Knepley test: 129330602db0SMatthew G. Knepley suffix: adv_gauss_im_2 129430602db0SMatthew G. Knepley requires: !complex 129530602db0SMatthew G. Knepley TODO: broken 129630602db0SMatthew 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 129730602db0SMatthew G. Knepley 129830602db0SMatthew G. Knepley test: 129930602db0SMatthew G. Knepley suffix: adv_gauss_corner 130030602db0SMatthew G. Knepley requires: !complex 130130602db0SMatthew G. Knepley TODO: broken 130230602db0SMatthew 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 130330602db0SMatthew G. Knepley 130430602db0SMatthew G. Knepley # 2D Advection+Harmonic 12- 130530602db0SMatthew G. Knepley test: 130630602db0SMatthew G. Knepley suffix: adv_harm_0 130730602db0SMatthew G. Knepley requires: !complex 130830602db0SMatthew G. Knepley TODO: broken memory corruption 130930602db0SMatthew 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 131030602db0SMatthew G. Knepley 1311c4762a1bSJed Brown # Must check that FV BCs propagate to coarse meshes 1312c4762a1bSJed Brown # Must check that FV BC ids propagate to coarse meshes 1313c4762a1bSJed Brown # Must check that FE+FV BCs work at the same time 1314c4762a1bSJed Brown # 2D Advection, matching wind in ex11 8-11 1315c4762a1bSJed Brown # NOTE implicit solves are limited by accuracy of FD Jacobian 1316c4762a1bSJed Brown test: 1317c4762a1bSJed Brown suffix: adv_0 1318c4762a1bSJed Brown requires: !complex !single exodusii 131930602db0SMatthew 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 1320c4762a1bSJed Brown 1321c4762a1bSJed Brown test: 1322c4762a1bSJed Brown suffix: adv_0_im 1323c4762a1bSJed Brown requires: !complex exodusii 1324c4762a1bSJed Brown TODO: broken memory corruption 132530602db0SMatthew 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 1326c4762a1bSJed Brown 1327c4762a1bSJed Brown test: 1328c4762a1bSJed Brown suffix: adv_0_im_2 1329c4762a1bSJed Brown requires: !complex exodusii 1330c4762a1bSJed Brown TODO: broken 133130602db0SMatthew 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 1332c4762a1bSJed Brown 1333c4762a1bSJed Brown test: 1334c4762a1bSJed Brown suffix: adv_0_im_3 1335c4762a1bSJed Brown requires: !complex exodusii 1336c4762a1bSJed Brown TODO: broken 133730602db0SMatthew 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 1338c4762a1bSJed Brown 1339c4762a1bSJed Brown test: 1340c4762a1bSJed Brown suffix: adv_0_im_4 1341c4762a1bSJed Brown requires: !complex exodusii 1342c4762a1bSJed Brown TODO: broken 134330602db0SMatthew 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 1344c4762a1bSJed Brown # 2D Advection, misc 1345c4762a1bSJed Brown 1346c4762a1bSJed Brown TEST*/ 1347