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 72c4762a1bSJed Brown PetscFunctionBeginUser; 73c4762a1bSJed Brown options->useFV = PETSC_FALSE; 74c4762a1bSJed Brown options->velocityDist = VEL_HARMONIC; 75c4762a1bSJed Brown options->porosityDist = ZERO; 76c4762a1bSJed Brown options->inflowState = -2.0; 77c4762a1bSJed Brown options->numMonitorFuncs = 0; 78c4762a1bSJed Brown options->source[0] = 0.5; 79c4762a1bSJed Brown options->source[1] = 0.5; 80c4762a1bSJed Brown options->source[2] = 0.5; 81c4762a1bSJed Brown 82d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX"); 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL)); 84c4762a1bSJed Brown vd = options->velocityDist; 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-velocity_dist","Velocity distribution type","ex18.c",velocityDist,4,velocityDist[options->velocityDist],&vd,NULL)); 86c4762a1bSJed Brown options->velocityDist = (VelocityDistribution) vd; 87c4762a1bSJed Brown pd = options->porosityDist; 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-porosity_dist","Initial porosity distribution type","ex18.c",porosityDist,5,porosityDist[options->porosityDist],&pd,NULL)); 89c4762a1bSJed Brown options->porosityDist = (PorosityDistribution) pd; 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL)); 9130602db0SMatthew G. Knepley d = 2; 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg)); 9363a3b9bcSJacob Faibussowitsch PetscCheck(!flg || d == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %" PetscInt_FMT, d); 94d0609cedSBarry Smith PetscOptionsEnd(); 95c4762a1bSJed Brown 96c4762a1bSJed Brown PetscFunctionReturn(0); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99c4762a1bSJed Brown static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) 100c4762a1bSJed Brown { 101c4762a1bSJed Brown Functional func; 102c4762a1bSJed Brown char *names[256]; 103c4762a1bSJed Brown PetscInt f; 104c4762a1bSJed Brown 105c4762a1bSJed Brown PetscFunctionBeginUser; 106d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX"); 107dd39110bSPierre Jolivet options->numMonitorFuncs = PETSC_STATIC_ARRAY_LENGTH(names); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL)); 1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs)); 110c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 111c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 112c4762a1bSJed Brown PetscBool match; 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(names[f], func->name, &match)); 115c4762a1bSJed Brown if (match) break; 116c4762a1bSJed Brown } 1173c633725SBarry Smith PetscCheck(func,comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 118c4762a1bSJed Brown options->monitorFuncs[f] = func; 119c4762a1bSJed Brown /* Jed inserts a de-duplication of functionals here */ 1209566063dSJacob Faibussowitsch PetscCall(PetscFree(names[f])); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 123c4762a1bSJed Brown options->maxMonitorFunc = -1; 124c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 125c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 126c4762a1bSJed Brown Functional call = options->monitorFuncs[f]; 127c4762a1bSJed Brown 128c4762a1bSJed Brown if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown } 131d0609cedSBarry Smith PetscOptionsEnd(); 132c4762a1bSJed Brown PetscFunctionReturn(0); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) 136c4762a1bSJed Brown { 137c4762a1bSJed Brown Functional *ptr, f; 138c4762a1bSJed Brown PetscInt lastoffset = -1; 139c4762a1bSJed Brown 140c4762a1bSJed Brown PetscFunctionBeginUser; 141c4762a1bSJed Brown for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 1429566063dSJacob Faibussowitsch PetscCall(PetscNew(&f)); 1439566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &f->name)); 144c4762a1bSJed Brown f->offset = lastoffset + 1; 145c4762a1bSJed Brown f->func = func; 146c4762a1bSJed Brown f->ctx = ctx; 147c4762a1bSJed Brown f->next = NULL; 148c4762a1bSJed Brown *ptr = f; 149c4762a1bSJed Brown *offset = f->offset; 150c4762a1bSJed Brown PetscFunctionReturn(0); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown static PetscErrorCode FunctionalDestroy(Functional *link) 154c4762a1bSJed Brown { 155c4762a1bSJed Brown Functional next, l; 156c4762a1bSJed Brown 157c4762a1bSJed Brown PetscFunctionBeginUser; 158c4762a1bSJed Brown if (!link) PetscFunctionReturn(0); 159c4762a1bSJed Brown l = *link; 160c4762a1bSJed Brown *link = NULL; 161c4762a1bSJed Brown for (; l; l=next) { 162c4762a1bSJed Brown next = l->next; 1639566063dSJacob Faibussowitsch PetscCall(PetscFree(l->name)); 1649566063dSJacob Faibussowitsch PetscCall(PetscFree(l)); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown PetscFunctionReturn(0); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 170c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 171c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 172c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 173c4762a1bSJed Brown { 174c4762a1bSJed Brown PetscInt comp; 175c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 179c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 180c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 181c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 182c4762a1bSJed Brown { 183c4762a1bSJed Brown PetscScalar wind[3] = {0.0, 0.0, 0.0}; 184c4762a1bSJed Brown PetscInt comp; 185c4762a1bSJed Brown 186c4762a1bSJed Brown constant_u_2d(dim, t, x, Nf, wind, NULL); 187c4762a1bSJed Brown for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp]; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 191c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 192c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 193c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 194c4762a1bSJed Brown { 195c4762a1bSJed Brown PetscInt comp; 196c4762a1bSJed Brown for (comp = 0; comp < dim*dim; ++comp) f1[comp] = 0.0; 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 200c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 201c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 202c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 203c4762a1bSJed Brown { 204c4762a1bSJed Brown PetscInt d; 205c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d*dim+d] = 1.0; 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 209c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 210c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 211c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 212c4762a1bSJed Brown { 213c4762a1bSJed Brown g0[0] = 1.0; 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 217c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 218c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 219c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown PetscInt comp; 222c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown 225c4762a1bSJed Brown static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 226c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 227c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 228c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 229c4762a1bSJed Brown { 230c4762a1bSJed Brown PetscInt comp, d; 231c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) { 232c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 233c4762a1bSJed Brown f1[comp*dim+d] = u_x[comp*dim+d]; 234c4762a1bSJed Brown } 235c4762a1bSJed Brown } 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 239c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 240c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 241c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 242c4762a1bSJed Brown { 243c4762a1bSJed Brown f0[0] = -PetscSinReal(2.0*PETSC_PI*x[0]); 244c4762a1bSJed Brown f0[1] = 2.0*PETSC_PI*x[1]*PetscCosReal(2.0*PETSC_PI*x[0]); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 248c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 249c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 250c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 251c4762a1bSJed Brown { 252c4762a1bSJed Brown f0[0] = -2.0*PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1]); 253c4762a1bSJed Brown f0[1] = 2.0*PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0]); 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256c4762a1bSJed Brown void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 257c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 258c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 259c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 260c4762a1bSJed Brown { 261c4762a1bSJed Brown const PetscInt Ncomp = dim; 262c4762a1bSJed Brown PetscInt compI, d; 263c4762a1bSJed Brown 264c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 265c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 266c4762a1bSJed Brown g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0; 267c4762a1bSJed Brown } 268c4762a1bSJed Brown } 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 272c4762a1bSJed Brown static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 273c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 274c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 275c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 276c4762a1bSJed Brown { 277c4762a1bSJed Brown PetscInt d; 278c4762a1bSJed Brown f0[0] = u_t[dim]; 279c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += u[dim]*u_x[d*dim+d] + u_x[dim*dim+d]*u[d]; 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 283c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 284c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 285c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 286c4762a1bSJed Brown { 287c4762a1bSJed Brown PetscInt d; 288c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[0] = 0.0; 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 292c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 293c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 294c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 295c4762a1bSJed Brown { 296c4762a1bSJed Brown PetscInt d; 297c4762a1bSJed Brown g0[0] = u_tShift; 298c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[d*dim+d]; 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 301c4762a1bSJed Brown void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 302c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 303c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 304c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 305c4762a1bSJed Brown { 306c4762a1bSJed Brown PetscInt d; 307c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = u[d]; 308c4762a1bSJed Brown } 309c4762a1bSJed Brown 310c4762a1bSJed Brown void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 311c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 312c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 313c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 314c4762a1bSJed Brown { 315c4762a1bSJed Brown PetscInt d; 316c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[dim*dim+d]; 317c4762a1bSJed Brown } 318c4762a1bSJed Brown 319c4762a1bSJed Brown void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 320c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 321c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 322c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 323c4762a1bSJed Brown { 324c4762a1bSJed Brown PetscInt d; 325c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = u[dim]; 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 328c4762a1bSJed 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) 329c4762a1bSJed Brown { 330c4762a1bSJed Brown PetscReal wind[3] = {0.0, 1.0, 0.0}; 331c4762a1bSJed Brown PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim,3), wind, n); 332c4762a1bSJed Brown 333c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown 336c4762a1bSJed 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) 337c4762a1bSJed Brown { 338c4762a1bSJed Brown PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 339c4762a1bSJed Brown 340c4762a1bSJed Brown #if 1 341c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 342c4762a1bSJed Brown #else 343c4762a1bSJed 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]); */ 344c4762a1bSJed Brown /* Smear it out */ 345c4762a1bSJed Brown flux[0] = 0.5*((uL[dim] + uR[dim]) + (uL[dim] - uR[dim])*tanh(1.0e5*wn)) * wn; 346c4762a1bSJed Brown #endif 347c4762a1bSJed Brown } 348c4762a1bSJed Brown 349c4762a1bSJed Brown static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 350c4762a1bSJed Brown { 351c4762a1bSJed Brown u[0] = 0.0; 352c4762a1bSJed Brown u[1] = 0.0; 353c4762a1bSJed Brown return 0; 354c4762a1bSJed Brown } 355c4762a1bSJed Brown 356c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 357c4762a1bSJed Brown { 358c4762a1bSJed Brown u[0] = 0.0; 359c4762a1bSJed Brown u[1] = 1.0; 360c4762a1bSJed Brown return 0; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363c4762a1bSJed Brown /* Coordinates of the point which was at x at t = 0 */ 364c4762a1bSJed Brown static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 365c4762a1bSJed Brown { 366c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 367c4762a1bSJed Brown u[0] = x[0]; 368c4762a1bSJed Brown u[1] = x[1] + t; 36930602db0SMatthew G. Knepley #if 0 3709566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u)); 37130602db0SMatthew G. Knepley #else 37230602db0SMatthew G. Knepley u[1] = u[1] - (int) PetscRealPart(u[1]); 37330602db0SMatthew G. Knepley #endif 374c4762a1bSJed Brown return 0; 375c4762a1bSJed Brown } 376c4762a1bSJed Brown 377c4762a1bSJed Brown /* 378c4762a1bSJed Brown In 2D we use the exact solution: 379c4762a1bSJed Brown 380c4762a1bSJed Brown u = x^2 + y^2 381c4762a1bSJed Brown v = 2 x^2 - 2xy 382c4762a1bSJed Brown phi = h(x + y + (u + v) t) 383c4762a1bSJed Brown f_x = f_y = 4 384c4762a1bSJed Brown 385c4762a1bSJed Brown so that 386c4762a1bSJed Brown 387c4762a1bSJed Brown -\Delta u + f = <-4, -4> + <4, 4> = 0 388c4762a1bSJed Brown {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 389c4762a1bSJed Brown h_t(x + y + (u + v) t) - u . grad phi - phi div u 390c4762a1bSJed Brown = u h' + v h' - u h_x - v h_y 391c4762a1bSJed Brown = 0 392c4762a1bSJed Brown 393c4762a1bSJed Brown We will conserve phi since 394c4762a1bSJed Brown 395c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 396c4762a1bSJed Brown 397c4762a1bSJed Brown Also try h((x + ut)^2 + (y + vt)^2), so that 398c4762a1bSJed Brown 399c4762a1bSJed Brown h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 400c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 401c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 402c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 403c4762a1bSJed Brown = 0 404c4762a1bSJed Brown 405c4762a1bSJed Brown */ 406c4762a1bSJed Brown static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 407c4762a1bSJed Brown { 408c4762a1bSJed Brown u[0] = x[0]*x[0] + x[1]*x[1]; 409c4762a1bSJed Brown u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 410c4762a1bSJed Brown return 0; 411c4762a1bSJed Brown } 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* 414c4762a1bSJed Brown In 2D we use the exact, periodic solution: 415c4762a1bSJed Brown 416c4762a1bSJed Brown u = sin(2 pi x)/4 pi^2 417c4762a1bSJed Brown v = -y cos(2 pi x)/2 pi 418c4762a1bSJed Brown phi = h(x + y + (u + v) t) 419c4762a1bSJed Brown f_x = -sin(2 pi x) 420c4762a1bSJed Brown f_y = 2 pi y cos(2 pi x) 421c4762a1bSJed Brown 422c4762a1bSJed Brown so that 423c4762a1bSJed Brown 424c4762a1bSJed Brown -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 425c4762a1bSJed Brown 426c4762a1bSJed Brown We will conserve phi since 427c4762a1bSJed Brown 428c4762a1bSJed Brown \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 429c4762a1bSJed Brown */ 430c4762a1bSJed Brown static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 431c4762a1bSJed Brown { 432c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 433c4762a1bSJed Brown u[1] = -x[1]*PetscCosReal(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI); 434c4762a1bSJed Brown return 0; 435c4762a1bSJed Brown } 436c4762a1bSJed Brown 437c4762a1bSJed Brown /* 438c4762a1bSJed Brown In 2D we use the exact, doubly periodic solution: 439c4762a1bSJed Brown 440c4762a1bSJed Brown u = sin(2 pi x) cos(2 pi y)/4 pi^2 441c4762a1bSJed Brown v = -sin(2 pi y) cos(2 pi x)/4 pi^2 442c4762a1bSJed Brown phi = h(x + y + (u + v) t) 443c4762a1bSJed Brown f_x = -2sin(2 pi x) cos(2 pi y) 444c4762a1bSJed Brown f_y = 2sin(2 pi y) cos(2 pi x) 445c4762a1bSJed Brown 446c4762a1bSJed Brown so that 447c4762a1bSJed Brown 448c4762a1bSJed 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 449c4762a1bSJed Brown 450c4762a1bSJed Brown We will conserve phi since 451c4762a1bSJed Brown 452c4762a1bSJed Brown \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 453c4762a1bSJed Brown */ 454c4762a1bSJed Brown static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 455c4762a1bSJed Brown { 456c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1])/PetscSqr(2.0*PETSC_PI); 457c4762a1bSJed Brown u[1] = -PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 458c4762a1bSJed Brown return 0; 459c4762a1bSJed Brown } 460c4762a1bSJed Brown 461c4762a1bSJed Brown static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 462c4762a1bSJed Brown { 463c4762a1bSJed Brown u[0] = x[1] - 0.5; 464c4762a1bSJed Brown u[1] = 0.0; 465c4762a1bSJed Brown return 0; 466c4762a1bSJed Brown } 467c4762a1bSJed Brown 468c4762a1bSJed Brown static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 469c4762a1bSJed Brown { 470c4762a1bSJed Brown PetscInt d; 471c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 472c4762a1bSJed Brown return 0; 473c4762a1bSJed Brown } 474c4762a1bSJed Brown 475c4762a1bSJed Brown static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 476c4762a1bSJed Brown { 477c4762a1bSJed Brown u[0] = 0.0; 478c4762a1bSJed Brown return 0; 479c4762a1bSJed Brown } 480c4762a1bSJed Brown 481c4762a1bSJed Brown static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 482c4762a1bSJed Brown { 483c4762a1bSJed Brown u[0] = 1.0; 484c4762a1bSJed Brown return 0; 485c4762a1bSJed Brown } 486c4762a1bSJed Brown 487c4762a1bSJed Brown static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 488c4762a1bSJed Brown { 489c4762a1bSJed Brown PetscReal x0[2]; 490c4762a1bSJed Brown PetscScalar xn[2]; 491c4762a1bSJed Brown 492c4762a1bSJed Brown x0[0] = globalUser->source[0]; 493c4762a1bSJed Brown x0[1] = globalUser->source[1]; 494c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 495c4762a1bSJed Brown { 496c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 497c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 498c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 499c4762a1bSJed Brown 500c4762a1bSJed Brown u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 501c4762a1bSJed Brown } 502c4762a1bSJed Brown return 0; 503c4762a1bSJed Brown } 504c4762a1bSJed Brown 505c4762a1bSJed Brown /* 506c4762a1bSJed Brown Gaussian blob, initially centered on (0.5, 0.5) 507c4762a1bSJed Brown 508c4762a1bSJed Brown xi = x(t) - x0, eta = y(t) - y0 509c4762a1bSJed Brown 510c4762a1bSJed Brown where x(t), y(t) are the integral curves of v(t), 511c4762a1bSJed Brown 512c4762a1bSJed Brown dx/dt . grad f = v . f 513c4762a1bSJed Brown 514c4762a1bSJed Brown Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 515c4762a1bSJed Brown 516c4762a1bSJed Brown v0 f_x + w0 f_y = v . f 517c4762a1bSJed Brown */ 518c4762a1bSJed Brown static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 519c4762a1bSJed Brown { 520c4762a1bSJed Brown const PetscReal x0[2] = {0.5, 0.5}; 521c4762a1bSJed Brown const PetscReal sigma = 1.0/6.0; 522c4762a1bSJed Brown PetscScalar xn[2]; 523c4762a1bSJed Brown 524c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 525c4762a1bSJed Brown { 526c4762a1bSJed Brown /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 527c4762a1bSJed Brown /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 528c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 529c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 530c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 531c4762a1bSJed Brown 532c4762a1bSJed Brown u[0] = PetscExpReal(-r2/(2.0*sigma*sigma))/(sigma*PetscSqrtReal(2.0*PETSC_PI)); 533c4762a1bSJed Brown } 534c4762a1bSJed Brown return 0; 535c4762a1bSJed Brown } 536c4762a1bSJed Brown 537c4762a1bSJed Brown static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 538c4762a1bSJed Brown { 539c4762a1bSJed Brown PetscReal x0[3]; 540c4762a1bSJed Brown const PetscReal wind[3] = {0.0, 1.0, 0.0}; 541c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 542c4762a1bSJed Brown 543c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 544c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 545c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 546c4762a1bSJed Brown return 0; 547c4762a1bSJed Brown } 548c4762a1bSJed Brown 549c4762a1bSJed Brown static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 550c4762a1bSJed Brown { 551c4762a1bSJed Brown PetscReal ur[3]; 552c4762a1bSJed Brown PetscReal x0[3]; 553c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 554c4762a1bSJed Brown 555c4762a1bSJed Brown ur[0] = PetscRealPart(u[0]); ur[1] = PetscRealPart(u[1]); ur[2] = PetscRealPart(u[2]); 556c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 557c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 558c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 559c4762a1bSJed Brown return 0; 560c4762a1bSJed Brown } 561c4762a1bSJed Brown 562c4762a1bSJed Brown static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 563c4762a1bSJed Brown { 564c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 565c4762a1bSJed Brown 566c4762a1bSJed Brown PetscFunctionBeginUser; 567c4762a1bSJed Brown xG[0] = user->inflowState; 568c4762a1bSJed Brown PetscFunctionReturn(0); 569c4762a1bSJed Brown } 570c4762a1bSJed Brown 571c4762a1bSJed Brown static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 572c4762a1bSJed Brown { 573c4762a1bSJed Brown PetscFunctionBeginUser; 57430602db0SMatthew G. Knepley //xG[0] = xI[dim]; 57530602db0SMatthew G. Knepley xG[0] = xI[2]; 576c4762a1bSJed Brown PetscFunctionReturn(0); 577c4762a1bSJed Brown } 578c4762a1bSJed Brown 579c4762a1bSJed Brown static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 580c4762a1bSJed Brown { 581c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 582c4762a1bSJed Brown PetscInt dim; 583c4762a1bSJed Brown 584c4762a1bSJed Brown PetscFunctionBeginUser; 5859566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 586c4762a1bSJed Brown switch (user->porosityDist) { 587c4762a1bSJed Brown case TILTED: 588c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) tilted_phi_2d(dim, time, x, 2, u, (void *) &time); 589c4762a1bSJed Brown else tilted_phi_coupled_2d(dim, time, x, 2, u, (void *) &time); 590c4762a1bSJed Brown break; 591c4762a1bSJed Brown case GAUSSIAN: 592c4762a1bSJed Brown gaussian_phi_2d(dim, time, x, 2, u, (void *) &time); 593c4762a1bSJed Brown break; 594c4762a1bSJed Brown case DELTA: 595c4762a1bSJed Brown delta_phi_2d(dim, time, x, 2, u, (void *) &time); 596c4762a1bSJed Brown break; 597c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 598c4762a1bSJed Brown } 599c4762a1bSJed Brown PetscFunctionReturn(0); 600c4762a1bSJed Brown } 601c4762a1bSJed Brown 602c4762a1bSJed Brown static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 603c4762a1bSJed Brown { 604c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 605c4762a1bSJed Brown PetscScalar yexact[3]={0,0,0}; 606c4762a1bSJed Brown 607c4762a1bSJed Brown PetscFunctionBeginUser; 6089566063dSJacob Faibussowitsch PetscCall(ExactSolution(dm, time, x, yexact, ctx)); 609c4762a1bSJed Brown f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 610c4762a1bSJed Brown PetscFunctionReturn(0); 611c4762a1bSJed Brown } 612c4762a1bSJed Brown 613c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 614c4762a1bSJed Brown { 615c4762a1bSJed Brown PetscFunctionBeginUser; 6169566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 6179566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 6189566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 6199566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 620c4762a1bSJed Brown PetscFunctionReturn(0); 621c4762a1bSJed Brown } 622c4762a1bSJed Brown 623c4762a1bSJed Brown static PetscErrorCode SetupBC(DM dm, AppCtx *user) 624c4762a1bSJed Brown { 625348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 62630602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 627c4762a1bSJed Brown PetscDS prob; 628c4762a1bSJed Brown DMLabel label; 629c4762a1bSJed Brown PetscBool check; 63030602db0SMatthew G. Knepley PetscInt dim, n = 3; 63130602db0SMatthew G. Knepley const char *prefix; 632c4762a1bSJed Brown 633c4762a1bSJed Brown PetscFunctionBeginUser; 6349566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 6359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 6369566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 637c4762a1bSJed Brown /* Set initial guesses and exact solutions */ 63830602db0SMatthew G. Knepley switch (dim) { 639c4762a1bSJed Brown case 2: 640c4762a1bSJed Brown user->initialGuess[0] = initialVelocity; 641c4762a1bSJed Brown switch(user->porosityDist) { 642c4762a1bSJed Brown case ZERO: user->initialGuess[1] = zero_phi;break; 643c4762a1bSJed Brown case CONSTANT: user->initialGuess[1] = constant_phi;break; 644c4762a1bSJed Brown case GAUSSIAN: user->initialGuess[1] = gaussian_phi_2d;break; 645c4762a1bSJed Brown case DELTA: user->initialGuess[1] = delta_phi_2d;break; 646c4762a1bSJed Brown case TILTED: 647c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 648c4762a1bSJed Brown else user->initialGuess[1] = tilted_phi_coupled_2d; 649c4762a1bSJed Brown break; 650c4762a1bSJed Brown } 651348a1646SMatthew G. Knepley break; 65263a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 653348a1646SMatthew G. Knepley } 654348a1646SMatthew G. Knepley exactFuncs[0] = user->initialGuess[0]; 655348a1646SMatthew G. Knepley exactFuncs[1] = user->initialGuess[1]; 65630602db0SMatthew G. Knepley switch (dim) { 657348a1646SMatthew G. Knepley case 2: 658c4762a1bSJed Brown switch (user->velocityDist) { 659c4762a1bSJed Brown case VEL_ZERO: 660348a1646SMatthew G. Knepley exactFuncs[0] = zero_u_2d; break; 661c4762a1bSJed Brown case VEL_CONSTANT: 662348a1646SMatthew G. Knepley exactFuncs[0] = constant_u_2d; break; 663c4762a1bSJed Brown case VEL_HARMONIC: 66430602db0SMatthew G. Knepley switch (bdt[0]) { 665c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 66630602db0SMatthew G. Knepley switch (bdt[1]) { 667c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 668348a1646SMatthew G. Knepley exactFuncs[0] = doubly_periodic_u_2d; break; 669c4762a1bSJed Brown default: 670348a1646SMatthew G. Knepley exactFuncs[0] = periodic_u_2d; break; 671c4762a1bSJed Brown } 672c4762a1bSJed Brown break; 673c4762a1bSJed Brown default: 674348a1646SMatthew G. Knepley exactFuncs[0] = quadratic_u_2d; break; 675c4762a1bSJed Brown } 676c4762a1bSJed Brown break; 677c4762a1bSJed Brown case VEL_SHEAR: 678348a1646SMatthew G. Knepley exactFuncs[0] = shear_bc; break; 67963a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 680c4762a1bSJed Brown } 681348a1646SMatthew G. Knepley break; 68263a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 683c4762a1bSJed Brown } 684c4762a1bSJed Brown { 685c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 686c4762a1bSJed Brown 6879566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 688348a1646SMatthew G. Knepley if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 689c4762a1bSJed Brown } 6909566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-dmts_check", &check)); 691c4762a1bSJed Brown if (check) { 692348a1646SMatthew G. Knepley user->initialGuess[0] = exactFuncs[0]; 693348a1646SMatthew G. Knepley user->initialGuess[1] = exactFuncs[1]; 694c4762a1bSJed Brown } 695c4762a1bSJed Brown /* Set BC */ 6969566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 6979566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 6989566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user)); 6999566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user)); 700c4762a1bSJed Brown if (label) { 701c4762a1bSJed Brown const PetscInt id = 1; 702c4762a1bSJed Brown 7039566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL)); 704c4762a1bSJed Brown } 7059566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &label)); 706c4762a1bSJed Brown if (label && user->useFV) { 707c4762a1bSJed Brown const PetscInt inflowids[] = {100,200,300}, outflowids[] = {101}; 708c4762a1bSJed Brown 709dd39110bSPierre 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)); 710dd39110bSPierre 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)); 711c4762a1bSJed Brown } 712c4762a1bSJed Brown PetscFunctionReturn(0); 713c4762a1bSJed Brown } 714c4762a1bSJed Brown 715c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 716c4762a1bSJed Brown { 71730602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 718c4762a1bSJed Brown PetscDS prob; 71930602db0SMatthew G. Knepley PetscInt n = 3; 72030602db0SMatthew G. Knepley const char *prefix; 721c4762a1bSJed Brown 722c4762a1bSJed Brown PetscFunctionBeginUser; 7239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 7249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *) bdt, &n, NULL)); 7259566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 726c4762a1bSJed Brown switch (user->velocityDist) { 727c4762a1bSJed Brown case VEL_ZERO: 7289566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u)); 729c4762a1bSJed Brown break; 730c4762a1bSJed Brown case VEL_CONSTANT: 7319566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u)); 7329566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL)); 7339566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL)); 734c4762a1bSJed Brown break; 735c4762a1bSJed Brown case VEL_HARMONIC: 73630602db0SMatthew G. Knepley switch (bdt[0]) { 737c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 73830602db0SMatthew G. Knepley switch (bdt[1]) { 739c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 7409566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u)); 741c4762a1bSJed Brown break; 742c4762a1bSJed Brown default: 7439566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u)); 744c4762a1bSJed Brown break; 745c4762a1bSJed Brown } 746c4762a1bSJed Brown break; 747c4762a1bSJed Brown default: 7489566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u)); 749c4762a1bSJed Brown break; 750c4762a1bSJed Brown } 7519566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 752c4762a1bSJed Brown break; 753c4762a1bSJed Brown case VEL_SHEAR: 7549566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u)); 7559566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 756c4762a1bSJed Brown break; 757c4762a1bSJed Brown } 7589566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_advection, f1_advection)); 7599566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL)); 7609566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL)); 7619566063dSJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_advection)); 7629566063dSJacob Faibussowitsch else PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection)); 763c4762a1bSJed Brown 7649566063dSJacob Faibussowitsch PetscCall(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user)); 765c4762a1bSJed Brown PetscFunctionReturn(0); 766c4762a1bSJed Brown } 767c4762a1bSJed Brown 768c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 769c4762a1bSJed Brown { 770c4762a1bSJed Brown DM cdm = dm; 771c4762a1bSJed Brown PetscQuadrature q; 772c4762a1bSJed Brown PetscFE fe[2]; 773c4762a1bSJed Brown PetscFV fv; 774c4762a1bSJed Brown MPI_Comm comm; 77530602db0SMatthew G. Knepley PetscInt dim; 776c4762a1bSJed Brown 777c4762a1bSJed Brown PetscFunctionBeginUser; 778c4762a1bSJed Brown /* Create finite element */ 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 7809566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7819566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0])); 7829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity")); 7839566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1])); 7849566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 7859566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[1], "porosity")); 786c4762a1bSJed Brown 7879566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(PetscObjectComm((PetscObject) dm), &fv)); 7889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fv, "porosity")); 7899566063dSJacob Faibussowitsch PetscCall(PetscFVSetFromOptions(fv)); 7909566063dSJacob Faibussowitsch PetscCall(PetscFVSetNumComponents(fv, 1)); 7919566063dSJacob Faibussowitsch PetscCall(PetscFVSetSpatialDimension(fv, dim)); 7929566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe[0], &q)); 7939566063dSJacob Faibussowitsch PetscCall(PetscFVSetQuadrature(fv, q)); 794c4762a1bSJed Brown 7959566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 7969566063dSJacob Faibussowitsch if (user->useFV) PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fv)); 7979566063dSJacob Faibussowitsch else PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 7989566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 7999566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 800c4762a1bSJed Brown 801c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 802c4762a1bSJed Brown while (cdm) { 8039566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 8049566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 805c4762a1bSJed Brown /* Coordinates were never localized for coarse meshes */ 8069566063dSJacob Faibussowitsch if (cdm) PetscCall(DMLocalizeCoordinates(cdm)); 807c4762a1bSJed Brown } 8089566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 8099566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 8109566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fv)); 811c4762a1bSJed Brown PetscFunctionReturn(0); 812c4762a1bSJed Brown } 813c4762a1bSJed Brown 814c4762a1bSJed Brown static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 815c4762a1bSJed Brown { 816c4762a1bSJed Brown PetscFunctionBeginUser; 8179566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, user, dm)); 818c4762a1bSJed Brown /* Handle refinement, etc. */ 8199566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 820c4762a1bSJed Brown /* Construct ghost cells */ 821c4762a1bSJed Brown if (user->useFV) { 822c4762a1bSJed Brown DM gdm; 823c4762a1bSJed Brown 8249566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm)); 8259566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 826c4762a1bSJed Brown *dm = gdm; 827c4762a1bSJed Brown } 828c4762a1bSJed Brown /* Localize coordinates */ 8299566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 8309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)(*dm),"Mesh")); 8319566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 832c4762a1bSJed Brown /* Setup problem */ 8339566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*dm, user)); 834c4762a1bSJed Brown /* Setup BC */ 8359566063dSJacob Faibussowitsch PetscCall(SetupBC(*dm, user)); 836c4762a1bSJed Brown PetscFunctionReturn(0); 837c4762a1bSJed Brown } 838c4762a1bSJed Brown 839c4762a1bSJed Brown static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx) 840c4762a1bSJed Brown { 841c4762a1bSJed Brown PetscDS prob; 842c4762a1bSJed Brown DM dmCell; 843c4762a1bSJed Brown Vec cellgeom; 844c4762a1bSJed Brown const PetscScalar *cgeom; 845c4762a1bSJed Brown PetscScalar *x; 846c4762a1bSJed Brown PetscInt dim, Nf, cStart, cEnd, c; 847c4762a1bSJed Brown 848c4762a1bSJed Brown PetscFunctionBeginUser; 8499566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 8509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8519566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 8529566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8539566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 8549566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 8559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 8569566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 857c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 858c4762a1bSJed Brown PetscFVCellGeom *cg; 859c4762a1bSJed Brown PetscScalar *xc; 860c4762a1bSJed Brown 8619566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 8629566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc)); 8639566063dSJacob Faibussowitsch if (xc) PetscCall((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx)); 864c4762a1bSJed Brown } 8659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 8669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 867c4762a1bSJed Brown PetscFunctionReturn(0); 868c4762a1bSJed Brown } 869c4762a1bSJed Brown 870c4762a1bSJed Brown static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 871c4762a1bSJed Brown { 872c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 873c4762a1bSJed Brown char *ftable = NULL; 874c4762a1bSJed Brown DM dm; 875c4762a1bSJed Brown PetscSection s; 876c4762a1bSJed Brown Vec cellgeom; 877c4762a1bSJed Brown const PetscScalar *x; 878c4762a1bSJed Brown PetscScalar *a; 879c4762a1bSJed Brown PetscReal *xnorms; 880c4762a1bSJed Brown PetscInt pStart, pEnd, p, Nf, f; 881c4762a1bSJed Brown 882c4762a1bSJed Brown PetscFunctionBeginUser; 8839566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X, (PetscObject) ts, "-view_solution")); 8849566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 8859566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 8869566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 8879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &Nf)); 8889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 8899566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nf*2, &xnorms)); 8909566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 891c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 892c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 893c4762a1bSJed Brown PetscInt dof, cdof, d; 894c4762a1bSJed Brown 8959566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 8969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 8979566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, x, &a)); 898c4762a1bSJed Brown /* TODO Use constrained indices here */ 899c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+0] = PetscMax(xnorms[f*2+0], PetscAbsScalar(a[d])); 900c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+1] += PetscAbsScalar(a[d]); 901c4762a1bSJed Brown } 902c4762a1bSJed Brown } 9039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 904c4762a1bSJed Brown if (stepnum >= 0) { /* No summary for final time */ 905c4762a1bSJed Brown DM dmCell, *fdm; 906c4762a1bSJed Brown Vec *fv; 907c4762a1bSJed Brown const PetscScalar *cgeom; 908c4762a1bSJed Brown PetscScalar **fx; 909c4762a1bSJed Brown PetscReal *fmin, *fmax, *fint, *ftmp, t; 910c4762a1bSJed Brown PetscInt cStart, cEnd, c, fcount, f, num; 911c4762a1bSJed Brown 912c4762a1bSJed Brown size_t ftableused,ftablealloc; 913c4762a1bSJed Brown 914c4762a1bSJed Brown /* Functionals have indices after registering, this is an upper bound */ 915c4762a1bSJed Brown fcount = user->numMonitorFuncs; 9169566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fint,fcount,&ftmp)); 9179566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fcount,&fdm,fcount,&fv,fcount,&fx)); 918c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 919c4762a1bSJed Brown PetscSection fs; 920c4762a1bSJed Brown const char *name = user->monitorFuncs[f]->name; 921c4762a1bSJed Brown 922c4762a1bSJed Brown fmin[f] = PETSC_MAX_REAL; 923c4762a1bSJed Brown fmax[f] = PETSC_MIN_REAL; 924c4762a1bSJed Brown fint[f] = 0; 925c4762a1bSJed Brown /* Make monitor vecs */ 9269566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &fdm[f])); 9279566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &num, &t)); 9289566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(fdm[f], num, t)); 9299566063dSJacob Faibussowitsch PetscCall(PetscSectionClone(s, &fs)); 9309566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 0, NULL)); 9319566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 1, name)); 9329566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(fdm[f], fs)); 9339566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&fs)); 9349566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm[f], &fv[f])); 9359566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fv[f], name)); 9369566063dSJacob Faibussowitsch PetscCall(VecGetArray(fv[f], &fx[f])); 937c4762a1bSJed Brown } 9389566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 9399566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 9409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 9419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 942c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 943c4762a1bSJed Brown PetscFVCellGeom *cg; 944c4762a1bSJed Brown PetscScalar *cx; 945c4762a1bSJed Brown 9469566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 9479566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx)); 948c4762a1bSJed Brown if (!cx) continue; /* not a global cell */ 949c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 950c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 951c4762a1bSJed Brown PetscScalar *fxc; 952c4762a1bSJed Brown 9539566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc)); 954c4762a1bSJed Brown /* I need to make it easier to get interpolated values here */ 9559566063dSJacob Faibussowitsch PetscCall((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx)); 956c4762a1bSJed Brown fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 957c4762a1bSJed Brown } 958c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 959c4762a1bSJed Brown fmin[f] = PetscMin(fmin[f], ftmp[f]); 960c4762a1bSJed Brown fmax[f] = PetscMax(fmax[f], ftmp[f]); 961c4762a1bSJed Brown fint[f] += cg->volume * ftmp[f]; 962c4762a1bSJed Brown } 963c4762a1bSJed Brown } 9649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 9659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 9669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 9679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 9689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 969c4762a1bSJed Brown /* Output functional data */ 970c4762a1bSJed Brown ftablealloc = fcount * 100; 971c4762a1bSJed Brown ftableused = 0; 9729566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ftablealloc, &ftable)); 973c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 974c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 975c4762a1bSJed Brown PetscInt id = func->offset; 976c4762a1bSJed Brown char newline[] = "\n"; 977c4762a1bSJed Brown char buffer[256], *p, *prefix; 978c4762a1bSJed Brown size_t countused, len; 979c4762a1bSJed Brown 980c4762a1bSJed Brown /* Create string with functional outputs */ 981c4762a1bSJed Brown if (f % 3) { 9829566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, " ", 2)); 983c4762a1bSJed Brown p = buffer + 2; 984c4762a1bSJed Brown } else if (f) { 9859566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, newline, sizeof(newline)-1)); 986c4762a1bSJed Brown p = buffer + sizeof(newline) - 1; 987c4762a1bSJed Brown } else { 988c4762a1bSJed Brown p = buffer; 989c4762a1bSJed Brown } 9909566063dSJacob 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])); 991c4762a1bSJed Brown countused += p - buffer; 992c4762a1bSJed Brown /* reallocate */ 993c4762a1bSJed Brown if (countused > ftablealloc-ftableused-1) { 994c4762a1bSJed Brown char *ftablenew; 995c4762a1bSJed Brown 996c4762a1bSJed Brown ftablealloc = 2*ftablealloc + countused; 9979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ftablealloc, &ftablenew)); 9989566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); 9999566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 1000c4762a1bSJed Brown ftable = ftablenew; 1001c4762a1bSJed Brown } 10029566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftable+ftableused, buffer, countused)); 1003c4762a1bSJed Brown ftableused += countused; 1004c4762a1bSJed Brown ftable[ftableused] = 0; 1005c4762a1bSJed Brown /* Output vecs */ 10069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fv[f], &fx[f])); 10079566063dSJacob Faibussowitsch PetscCall(PetscStrlen(func->name, &len)); 10089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len+2,&prefix)); 10099566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(prefix, func->name)); 10109566063dSJacob Faibussowitsch PetscCall(PetscStrcat(prefix, "_")); 10119566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix)); 10129566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fv[f], NULL, "-vec_view")); 10139566063dSJacob Faibussowitsch PetscCall(PetscFree(prefix)); 10149566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm[f], &fv[f])); 10159566063dSJacob Faibussowitsch PetscCall(DMDestroy(&fdm[f])); 1016c4762a1bSJed Brown } 10179566063dSJacob Faibussowitsch PetscCall(PetscFree4(fmin, fmax, fint, ftmp)); 10189566063dSJacob Faibussowitsch PetscCall(PetscFree3(fdm, fv, fx)); 101963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "% 3" PetscInt_FMT " time %8.4g |x| (", stepnum, (double) time)); 1020c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10219566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 10229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+0])); 1023c4762a1bSJed Brown } 10249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ") |x|_1 (")); 1025c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 10269566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ", ")); 10279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+1])); 1028c4762a1bSJed Brown } 10299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), ") %s\n", ftable ? ftable : "")); 10309566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 1031c4762a1bSJed Brown } 10329566063dSJacob Faibussowitsch PetscCall(PetscFree(xnorms)); 1033c4762a1bSJed Brown PetscFunctionReturn(0); 1034c4762a1bSJed Brown } 1035c4762a1bSJed Brown 1036c4762a1bSJed Brown int main(int argc, char **argv) 1037c4762a1bSJed Brown { 1038c4762a1bSJed Brown MPI_Comm comm; 1039c4762a1bSJed Brown TS ts; 1040c4762a1bSJed Brown DM dm; 1041c4762a1bSJed Brown Vec u; 1042c4762a1bSJed Brown AppCtx user; 1043c4762a1bSJed Brown PetscReal t0, t = 0.0; 1044c4762a1bSJed Brown void *ctxs[2]; 1045c4762a1bSJed Brown 1046c4762a1bSJed Brown ctxs[0] = &t; 1047c4762a1bSJed Brown ctxs[1] = &t; 1048*327415f7SBarry Smith PetscFunctionBeginUser; 10499566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char*) 0, help)); 1050c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1051c4762a1bSJed Brown user.functionalRegistry = NULL; 1052c4762a1bSJed Brown globalUser = &user; 10539566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 10549566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 10559566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 10569566063dSJacob Faibussowitsch PetscCall(CreateDM(comm, &user, &dm)); 10579566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 10589566063dSJacob Faibussowitsch PetscCall(ProcessMonitorOptions(comm, &user)); 1059c4762a1bSJed Brown 10609566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 10619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "solution")); 1062c4762a1bSJed Brown if (user.useFV) { 1063c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 1064c4762a1bSJed Brown 10659566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit)); 1066c4762a1bSJed Brown if (isImplicit) { 10679566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10689566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1069c4762a1bSJed Brown } 10709566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10719566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user)); 1072c4762a1bSJed Brown } else { 10739566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 10749566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 10759566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 1076c4762a1bSJed Brown } 10779566063dSJacob Faibussowitsch if (user.useFV) PetscCall(TSMonitorSet(ts, MonitorFunctionals, &user, NULL)); 10789566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1)); 10799566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2.0)); 10809566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,0.01)); 10819566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 10829566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1083c4762a1bSJed Brown 10849566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u)); 10859566063dSJacob Faibussowitsch if (user.useFV) PetscCall(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1])); 10869566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-init_vec_view")); 10879566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1088c4762a1bSJed Brown t0 = t; 10899566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 10909566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 10919566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 10929566063dSJacob Faibussowitsch if (t > t0) PetscCall(DMTSCheckFromOptions(ts, u)); 10939566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 1094c4762a1bSJed Brown { 1095c4762a1bSJed Brown PetscReal ftime; 1096c4762a1bSJed Brown PetscInt nsteps; 1097c4762a1bSJed Brown TSConvergedReason reason; 1098c4762a1bSJed Brown 10999566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 11009566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &nsteps)); 11019566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 110263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double) ftime, nsteps)); 1103c4762a1bSJed Brown } 1104c4762a1bSJed Brown 11059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 11069566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 11079566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 11089566063dSJacob Faibussowitsch PetscCall(PetscFree(user.monitorFuncs)); 11099566063dSJacob Faibussowitsch PetscCall(FunctionalDestroy(&user.functionalRegistry)); 11109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1111b122ec5aSJacob Faibussowitsch return 0; 1112c4762a1bSJed Brown } 1113c4762a1bSJed Brown 1114c4762a1bSJed Brown /*TEST 1115c4762a1bSJed Brown 111630602db0SMatthew G. Knepley testset: 111730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 111830602db0SMatthew G. Knepley 1119c4762a1bSJed Brown # 2D harmonic velocity, no porosity 1120c4762a1bSJed Brown test: 1121c4762a1bSJed Brown suffix: p1p1 1122c4762a1bSJed Brown requires: !complex !single 112330602db0SMatthew 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 1124c4762a1bSJed Brown test: 1125c4762a1bSJed Brown suffix: p1p1_xper 1126c4762a1bSJed Brown requires: !complex !single 112730602db0SMatthew 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 1128c4762a1bSJed Brown test: 1129c4762a1bSJed Brown suffix: p1p1_xper_ref 1130c4762a1bSJed Brown requires: !complex !single 113130602db0SMatthew 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 1132c4762a1bSJed Brown test: 1133c4762a1bSJed Brown suffix: p1p1_xyper 1134c4762a1bSJed Brown requires: !complex !single 113530602db0SMatthew 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 1136c4762a1bSJed Brown test: 1137c4762a1bSJed Brown suffix: p1p1_xyper_ref 1138c4762a1bSJed Brown requires: !complex !single 113930602db0SMatthew 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 1140c4762a1bSJed Brown test: 1141c4762a1bSJed Brown suffix: p2p1 1142c4762a1bSJed Brown requires: !complex !single 114330602db0SMatthew 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 1144c4762a1bSJed Brown test: 1145c4762a1bSJed Brown suffix: p2p1_xyper 1146c4762a1bSJed Brown requires: !complex !single 114730602db0SMatthew 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 114830602db0SMatthew G. Knepley 114930602db0SMatthew G. Knepley test: 115030602db0SMatthew G. Knepley suffix: adv_1 115130602db0SMatthew G. Knepley requires: !complex !single 115230602db0SMatthew 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 115330602db0SMatthew G. Knepley 115430602db0SMatthew G. Knepley test: 115530602db0SMatthew G. Knepley suffix: adv_2 115630602db0SMatthew G. Knepley requires: !complex 115730602db0SMatthew G. Knepley TODO: broken memory corruption 115830602db0SMatthew 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 115930602db0SMatthew G. Knepley 116030602db0SMatthew G. Knepley test: 116130602db0SMatthew G. Knepley suffix: adv_3 116230602db0SMatthew G. Knepley requires: !complex 116330602db0SMatthew G. Knepley TODO: broken memory corruption 116430602db0SMatthew 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 116530602db0SMatthew G. Knepley 116630602db0SMatthew G. Knepley test: 116730602db0SMatthew G. Knepley suffix: adv_3_ex 116830602db0SMatthew G. Knepley requires: !complex 116930602db0SMatthew 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 117030602db0SMatthew G. Knepley 117130602db0SMatthew G. Knepley test: 117230602db0SMatthew G. Knepley suffix: adv_4 117330602db0SMatthew G. Knepley requires: !complex 117430602db0SMatthew G. Knepley TODO: broken memory corruption 117530602db0SMatthew 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 117630602db0SMatthew G. Knepley 117730602db0SMatthew G. Knepley # 2D Advection, box, delta 117830602db0SMatthew G. Knepley test: 117930602db0SMatthew G. Knepley suffix: adv_delta_yper_0 118030602db0SMatthew G. Knepley requires: !complex 118130602db0SMatthew G. Knepley TODO: broken 118230602db0SMatthew 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 118330602db0SMatthew G. Knepley 118430602db0SMatthew G. Knepley test: 118530602db0SMatthew G. Knepley suffix: adv_delta_yper_1 118630602db0SMatthew G. Knepley requires: !complex 118730602db0SMatthew G. Knepley TODO: broken 118830602db0SMatthew 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 118930602db0SMatthew G. Knepley 119030602db0SMatthew G. Knepley test: 119130602db0SMatthew G. Knepley suffix: adv_delta_yper_2 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 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 119530602db0SMatthew G. Knepley 119630602db0SMatthew G. Knepley test: 119730602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_0 119830602db0SMatthew G. Knepley requires: !complex 119930602db0SMatthew G. Knepley TODO: broken 120030602db0SMatthew G. Knepley args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 0 -velocity_dist constant -porosity_dist delta -inflow_state 0.0 -ts_type mimex -ts_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 120130602db0SMatthew G. Knepley 120230602db0SMatthew G. Knepley test: 120330602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_1 120430602db0SMatthew G. Knepley requires: !complex 120530602db0SMatthew G. Knepley TODO: broken 120630602db0SMatthew G. Knepley args: -dm_plex_box_bd none,periodic -use_fv -use_implicit -velocity_petscspace_degree 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 120730602db0SMatthew G. Knepley 120830602db0SMatthew G. Knepley test: 120930602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_2 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 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 121330602db0SMatthew G. Knepley 121430602db0SMatthew G. Knepley test: 121530602db0SMatthew G. Knepley suffix: adv_delta_yper_im_0 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 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 121930602db0SMatthew G. Knepley 122030602db0SMatthew G. Knepley test: 122130602db0SMatthew G. Knepley suffix: adv_delta_yper_im_1 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 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 122530602db0SMatthew G. Knepley 122630602db0SMatthew G. Knepley test: 122730602db0SMatthew G. Knepley suffix: adv_delta_yper_im_2 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 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 123130602db0SMatthew G. Knepley 123230602db0SMatthew G. Knepley test: 123330602db0SMatthew G. Knepley suffix: adv_delta_yper_im_3 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 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 123730602db0SMatthew G. Knepley 123830602db0SMatthew G. Knepley # I believe the nullspace is sin(pi y) 123930602db0SMatthew G. Knepley test: 124030602db0SMatthew G. Knepley suffix: adv_delta_yper_im_4 124130602db0SMatthew G. Knepley requires: !complex 124230602db0SMatthew G. Knepley TODO: broken 124330602db0SMatthew 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 124430602db0SMatthew G. Knepley 124530602db0SMatthew G. Knepley test: 124630602db0SMatthew G. Knepley suffix: adv_delta_yper_im_5 124730602db0SMatthew G. Knepley requires: !complex 124830602db0SMatthew G. Knepley TODO: broken 124930602db0SMatthew 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 125030602db0SMatthew G. Knepley 125130602db0SMatthew G. Knepley test: 125230602db0SMatthew G. Knepley suffix: adv_delta_yper_im_6 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 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 125630602db0SMatthew G. Knepley # 2D Advection, magma benchmark 1 125730602db0SMatthew G. Knepley 125830602db0SMatthew G. Knepley test: 125930602db0SMatthew G. Knepley suffix: adv_delta_shear_im_0 126030602db0SMatthew G. Knepley requires: !complex 126130602db0SMatthew G. Knepley TODO: broken 126230602db0SMatthew 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 126330602db0SMatthew G. Knepley # 2D Advection, box, gaussian 126430602db0SMatthew G. Knepley 126530602db0SMatthew G. Knepley test: 126630602db0SMatthew G. Knepley suffix: adv_gauss 126730602db0SMatthew G. Knepley requires: !complex 126830602db0SMatthew G. Knepley TODO: broken 126930602db0SMatthew G. Knepley args: -use_fv -velocity_dist constant -porosity_dist gaussian -inflow_state 0.0 -ts_type ssp -ts_max_time 2.0 -ts_max_steps 100 -ts_dt 0.01 -bc_inflow 1 -bc_outflow 3 -ts_view -dm_view 127030602db0SMatthew G. Knepley 127130602db0SMatthew G. Knepley test: 127230602db0SMatthew G. Knepley suffix: adv_gauss_im 127330602db0SMatthew G. Knepley requires: !complex 127430602db0SMatthew G. Knepley TODO: broken 127530602db0SMatthew 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 127630602db0SMatthew G. Knepley 127730602db0SMatthew G. Knepley test: 127830602db0SMatthew G. Knepley suffix: adv_gauss_im_1 127930602db0SMatthew G. Knepley requires: !complex 128030602db0SMatthew G. Knepley TODO: broken 128130602db0SMatthew 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 128230602db0SMatthew G. Knepley 128330602db0SMatthew G. Knepley test: 128430602db0SMatthew G. Knepley suffix: adv_gauss_im_2 128530602db0SMatthew G. Knepley requires: !complex 128630602db0SMatthew G. Knepley TODO: broken 128730602db0SMatthew 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 128830602db0SMatthew G. Knepley 128930602db0SMatthew G. Knepley test: 129030602db0SMatthew G. Knepley suffix: adv_gauss_corner 129130602db0SMatthew G. Knepley requires: !complex 129230602db0SMatthew G. Knepley TODO: broken 129330602db0SMatthew 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 129430602db0SMatthew G. Knepley 129530602db0SMatthew G. Knepley # 2D Advection+Harmonic 12- 129630602db0SMatthew G. Knepley test: 129730602db0SMatthew G. Knepley suffix: adv_harm_0 129830602db0SMatthew G. Knepley requires: !complex 129930602db0SMatthew G. Knepley TODO: broken memory corruption 130030602db0SMatthew 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 130130602db0SMatthew G. Knepley 1302c4762a1bSJed Brown # Must check that FV BCs propagate to coarse meshes 1303c4762a1bSJed Brown # Must check that FV BC ids propagate to coarse meshes 1304c4762a1bSJed Brown # Must check that FE+FV BCs work at the same time 1305c4762a1bSJed Brown # 2D Advection, matching wind in ex11 8-11 1306c4762a1bSJed Brown # NOTE implicit solves are limited by accuracy of FD Jacobian 1307c4762a1bSJed Brown test: 1308c4762a1bSJed Brown suffix: adv_0 1309c4762a1bSJed Brown requires: !complex !single exodusii 131030602db0SMatthew 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 1311c4762a1bSJed Brown 1312c4762a1bSJed Brown test: 1313c4762a1bSJed Brown suffix: adv_0_im 1314c4762a1bSJed Brown requires: !complex exodusii 1315c4762a1bSJed Brown TODO: broken memory corruption 131630602db0SMatthew 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 1317c4762a1bSJed Brown 1318c4762a1bSJed Brown test: 1319c4762a1bSJed Brown suffix: adv_0_im_2 1320c4762a1bSJed Brown requires: !complex exodusii 1321c4762a1bSJed Brown TODO: broken 132230602db0SMatthew 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 1323c4762a1bSJed Brown 1324c4762a1bSJed Brown test: 1325c4762a1bSJed Brown suffix: adv_0_im_3 1326c4762a1bSJed Brown requires: !complex exodusii 1327c4762a1bSJed Brown TODO: broken 132830602db0SMatthew 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 1329c4762a1bSJed Brown 1330c4762a1bSJed Brown test: 1331c4762a1bSJed Brown suffix: adv_0_im_4 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_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 1335c4762a1bSJed Brown # 2D Advection, misc 1336c4762a1bSJed Brown 1337c4762a1bSJed Brown TEST*/ 1338