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 17*9371c9d4SSatish Balay typedef enum { 18*9371c9d4SSatish Balay VEL_ZERO, 19*9371c9d4SSatish Balay VEL_CONSTANT, 20*9371c9d4SSatish Balay VEL_HARMONIC, 21*9371c9d4SSatish Balay VEL_SHEAR 22*9371c9d4SSatish Balay } VelocityDistribution; 23c4762a1bSJed Brown 24*9371c9d4SSatish Balay typedef enum { 25*9371c9d4SSatish Balay ZERO, 26*9371c9d4SSatish Balay CONSTANT, 27*9371c9d4SSatish Balay GAUSSIAN, 28*9371c9d4SSatish Balay TILTED, 29*9371c9d4SSatish Balay DELTA 30*9371c9d4SSatish Balay } PorosityDistribution; 31c4762a1bSJed Brown 32c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown FunctionalFunc - Calculates the value of a functional of the solution at a point 36c4762a1bSJed Brown 37c4762a1bSJed Brown Input Parameters: 38c4762a1bSJed Brown + dm - The DM 39c4762a1bSJed Brown . time - The TS time 40c4762a1bSJed Brown . x - The coordinates of the evaluation point 41c4762a1bSJed Brown . u - The field values at point x 42c4762a1bSJed Brown - ctx - A user context, or NULL 43c4762a1bSJed Brown 44c4762a1bSJed Brown Output Parameter: 45c4762a1bSJed Brown . f - The value of the functional at point x 46c4762a1bSJed Brown 47c4762a1bSJed Brown */ 48c4762a1bSJed Brown typedef PetscErrorCode (*FunctionalFunc)(DM, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); 49c4762a1bSJed Brown 50c4762a1bSJed Brown typedef struct _n_Functional *Functional; 51c4762a1bSJed Brown struct _n_Functional { 52c4762a1bSJed Brown char *name; 53c4762a1bSJed Brown FunctionalFunc func; 54c4762a1bSJed Brown void *ctx; 55c4762a1bSJed Brown PetscInt offset; 56c4762a1bSJed Brown Functional next; 57c4762a1bSJed Brown }; 58c4762a1bSJed Brown 59c4762a1bSJed Brown typedef struct { 60c4762a1bSJed Brown /* Problem definition */ 61c4762a1bSJed Brown PetscBool useFV; /* Use a finite volume scheme for advection */ 62c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 63c4762a1bSJed Brown VelocityDistribution velocityDist; 64c4762a1bSJed Brown PorosityDistribution porosityDist; 65c4762a1bSJed Brown PetscReal inflowState; 66c4762a1bSJed Brown PetscReal source[3]; 67c4762a1bSJed Brown /* Monitoring */ 68c4762a1bSJed Brown PetscInt numMonitorFuncs, maxMonitorFunc; 69c4762a1bSJed Brown Functional *monitorFuncs; 70c4762a1bSJed Brown PetscInt errorFunctional; 71c4762a1bSJed Brown Functional functionalRegistry; 72c4762a1bSJed Brown } AppCtx; 73c4762a1bSJed Brown 74c4762a1bSJed Brown static AppCtx *globalUser; 75c4762a1bSJed Brown 76*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 77c4762a1bSJed Brown const char *velocityDist[4] = {"zero", "constant", "harmonic", "shear"}; 78c4762a1bSJed Brown const char *porosityDist[5] = {"zero", "constant", "gaussian", "tilted", "delta"}; 7930602db0SMatthew G. Knepley PetscInt vd, pd, d; 80c4762a1bSJed Brown PetscBool flg; 81c4762a1bSJed Brown 82c4762a1bSJed Brown PetscFunctionBeginUser; 83c4762a1bSJed Brown options->useFV = PETSC_FALSE; 84c4762a1bSJed Brown options->velocityDist = VEL_HARMONIC; 85c4762a1bSJed Brown options->porosityDist = ZERO; 86c4762a1bSJed Brown options->inflowState = -2.0; 87c4762a1bSJed Brown options->numMonitorFuncs = 0; 88c4762a1bSJed Brown options->source[0] = 0.5; 89c4762a1bSJed Brown options->source[1] = 0.5; 90c4762a1bSJed Brown options->source[2] = 0.5; 91c4762a1bSJed Brown 92d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX"); 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL)); 94c4762a1bSJed Brown vd = options->velocityDist; 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-velocity_dist", "Velocity distribution type", "ex18.c", velocityDist, 4, velocityDist[options->velocityDist], &vd, NULL)); 96c4762a1bSJed Brown options->velocityDist = (VelocityDistribution)vd; 97c4762a1bSJed Brown pd = options->porosityDist; 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-porosity_dist", "Initial porosity distribution type", "ex18.c", porosityDist, 5, porosityDist[options->porosityDist], &pd, NULL)); 99c4762a1bSJed Brown options->porosityDist = (PorosityDistribution)pd; 1009566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL)); 10130602db0SMatthew G. Knepley d = 2; 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg)); 10363a3b9bcSJacob Faibussowitsch PetscCheck(!flg || d == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %" PetscInt_FMT, d); 104d0609cedSBarry Smith PetscOptionsEnd(); 105c4762a1bSJed Brown 106c4762a1bSJed Brown PetscFunctionReturn(0); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109*9371c9d4SSatish Balay static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) { 110c4762a1bSJed Brown Functional func; 111c4762a1bSJed Brown char *names[256]; 112c4762a1bSJed Brown PetscInt f; 113c4762a1bSJed Brown 114c4762a1bSJed Brown PetscFunctionBeginUser; 115d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX"); 116dd39110bSPierre Jolivet options->numMonitorFuncs = PETSC_STATIC_ARRAY_LENGTH(names); 1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL)); 1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs)); 119c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 120c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 121c4762a1bSJed Brown PetscBool match; 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(names[f], func->name, &match)); 124c4762a1bSJed Brown if (match) break; 125c4762a1bSJed Brown } 1263c633725SBarry Smith PetscCheck(func, comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 127c4762a1bSJed Brown options->monitorFuncs[f] = func; 128c4762a1bSJed Brown /* Jed inserts a de-duplication of functionals here */ 1299566063dSJacob Faibussowitsch PetscCall(PetscFree(names[f])); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 132c4762a1bSJed Brown options->maxMonitorFunc = -1; 133c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 134c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 135c4762a1bSJed Brown Functional call = options->monitorFuncs[f]; 136c4762a1bSJed Brown 137c4762a1bSJed Brown if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown } 140d0609cedSBarry Smith PetscOptionsEnd(); 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144*9371c9d4SSatish Balay static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) { 145c4762a1bSJed Brown Functional *ptr, f; 146c4762a1bSJed Brown PetscInt lastoffset = -1; 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscFunctionBeginUser; 149c4762a1bSJed Brown for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 1509566063dSJacob Faibussowitsch PetscCall(PetscNew(&f)); 1519566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &f->name)); 152c4762a1bSJed Brown f->offset = lastoffset + 1; 153c4762a1bSJed Brown f->func = func; 154c4762a1bSJed Brown f->ctx = ctx; 155c4762a1bSJed Brown f->next = NULL; 156c4762a1bSJed Brown *ptr = f; 157c4762a1bSJed Brown *offset = f->offset; 158c4762a1bSJed Brown PetscFunctionReturn(0); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161*9371c9d4SSatish Balay static PetscErrorCode FunctionalDestroy(Functional *link) { 162c4762a1bSJed Brown Functional next, l; 163c4762a1bSJed Brown 164c4762a1bSJed Brown PetscFunctionBeginUser; 165c4762a1bSJed Brown if (!link) PetscFunctionReturn(0); 166c4762a1bSJed Brown l = *link; 167c4762a1bSJed Brown *link = NULL; 168c4762a1bSJed Brown for (; l; l = next) { 169c4762a1bSJed Brown next = l->next; 1709566063dSJacob Faibussowitsch PetscCall(PetscFree(l->name)); 1719566063dSJacob Faibussowitsch PetscCall(PetscFree(l)); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown PetscFunctionReturn(0); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176*9371c9d4SSatish Balay static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 177c4762a1bSJed Brown PetscInt comp; 178c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181*9371c9d4SSatish Balay static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 182c4762a1bSJed Brown PetscScalar wind[3] = {0.0, 0.0, 0.0}; 183c4762a1bSJed Brown PetscInt comp; 184c4762a1bSJed Brown 185c4762a1bSJed Brown constant_u_2d(dim, t, x, Nf, wind, NULL); 186c4762a1bSJed Brown for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp]; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189*9371c9d4SSatish Balay static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) { 190c4762a1bSJed Brown PetscInt comp; 191c4762a1bSJed Brown for (comp = 0; comp < dim * dim; ++comp) f1[comp] = 0.0; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194*9371c9d4SSatish Balay static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 195c4762a1bSJed Brown PetscInt d; 196c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d * dim + d] = 1.0; 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199*9371c9d4SSatish Balay static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 200c4762a1bSJed Brown g0[0] = 1.0; 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203*9371c9d4SSatish Balay static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 204c4762a1bSJed Brown PetscInt comp; 205c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208*9371c9d4SSatish Balay static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) { 209c4762a1bSJed Brown PetscInt comp, d; 210c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) { 211*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { f1[comp * dim + d] = u_x[comp * dim + d]; } 212c4762a1bSJed Brown } 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215*9371c9d4SSatish Balay static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 216c4762a1bSJed Brown f0[0] = -PetscSinReal(2.0 * PETSC_PI * x[0]); 217c4762a1bSJed Brown f0[1] = 2.0 * PETSC_PI * x[1] * PetscCosReal(2.0 * PETSC_PI * x[0]); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220*9371c9d4SSatish Balay static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 221c4762a1bSJed Brown f0[0] = -2.0 * PetscSinReal(2.0 * PETSC_PI * x[0]) * PetscCosReal(2.0 * PETSC_PI * x[1]); 222c4762a1bSJed Brown f0[1] = 2.0 * PetscSinReal(2.0 * PETSC_PI * x[1]) * PetscCosReal(2.0 * PETSC_PI * x[0]); 223c4762a1bSJed Brown } 224c4762a1bSJed Brown 225*9371c9d4SSatish Balay void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) { 226c4762a1bSJed Brown const PetscInt Ncomp = dim; 227c4762a1bSJed Brown PetscInt compI, d; 228c4762a1bSJed Brown 229c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 230*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { g3[((compI * Ncomp + compI) * dim + d) * dim + d] = 1.0; } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 235*9371c9d4SSatish Balay static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 236c4762a1bSJed Brown PetscInt d; 237c4762a1bSJed Brown f0[0] = u_t[dim]; 238c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += u[dim] * u_x[d * dim + d] + u_x[dim * dim + d] * u[d]; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241*9371c9d4SSatish Balay static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) { 242c4762a1bSJed Brown PetscInt d; 243c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[0] = 0.0; 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246*9371c9d4SSatish Balay void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 247c4762a1bSJed Brown PetscInt d; 248c4762a1bSJed Brown g0[0] = u_tShift; 249c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[d * dim + d]; 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252*9371c9d4SSatish Balay void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) { 253c4762a1bSJed Brown PetscInt d; 254c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = u[d]; 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257*9371c9d4SSatish Balay void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 258c4762a1bSJed Brown PetscInt d; 259c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[dim * dim + d]; 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262*9371c9d4SSatish Balay void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) { 263c4762a1bSJed Brown PetscInt d; 264c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = u[dim]; 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267*9371c9d4SSatish Balay 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) { 268c4762a1bSJed Brown PetscReal wind[3] = {0.0, 1.0, 0.0}; 269c4762a1bSJed Brown PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim, 3), wind, n); 270c4762a1bSJed Brown 271c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274*9371c9d4SSatish Balay 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) { 275c4762a1bSJed Brown PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 276c4762a1bSJed Brown 277c4762a1bSJed Brown #if 1 278c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 279c4762a1bSJed Brown #else 280c4762a1bSJed 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]); */ 281c4762a1bSJed Brown /* Smear it out */ 282c4762a1bSJed Brown flux[0] = 0.5 * ((uL[dim] + uR[dim]) + (uL[dim] - uR[dim]) * tanh(1.0e5 * wn)) * wn; 283c4762a1bSJed Brown #endif 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286*9371c9d4SSatish Balay static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 287c4762a1bSJed Brown u[0] = 0.0; 288c4762a1bSJed Brown u[1] = 0.0; 289c4762a1bSJed Brown return 0; 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292*9371c9d4SSatish Balay static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 293c4762a1bSJed Brown u[0] = 0.0; 294c4762a1bSJed Brown u[1] = 1.0; 295c4762a1bSJed Brown return 0; 296c4762a1bSJed Brown } 297c4762a1bSJed Brown 298c4762a1bSJed Brown /* Coordinates of the point which was at x at t = 0 */ 299*9371c9d4SSatish Balay static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 300c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 301c4762a1bSJed Brown u[0] = x[0]; 302c4762a1bSJed Brown u[1] = x[1] + t; 30330602db0SMatthew G. Knepley #if 0 3049566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u)); 30530602db0SMatthew G. Knepley #else 30630602db0SMatthew G. Knepley u[1] = u[1] - (int)PetscRealPart(u[1]); 30730602db0SMatthew G. Knepley #endif 308c4762a1bSJed Brown return 0; 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown /* 312c4762a1bSJed Brown In 2D we use the exact solution: 313c4762a1bSJed Brown 314c4762a1bSJed Brown u = x^2 + y^2 315c4762a1bSJed Brown v = 2 x^2 - 2xy 316c4762a1bSJed Brown phi = h(x + y + (u + v) t) 317c4762a1bSJed Brown f_x = f_y = 4 318c4762a1bSJed Brown 319c4762a1bSJed Brown so that 320c4762a1bSJed Brown 321c4762a1bSJed Brown -\Delta u + f = <-4, -4> + <4, 4> = 0 322c4762a1bSJed Brown {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 323c4762a1bSJed Brown h_t(x + y + (u + v) t) - u . grad phi - phi div u 324c4762a1bSJed Brown = u h' + v h' - u h_x - v h_y 325c4762a1bSJed Brown = 0 326c4762a1bSJed Brown 327c4762a1bSJed Brown We will conserve phi since 328c4762a1bSJed Brown 329c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 330c4762a1bSJed Brown 331c4762a1bSJed Brown Also try h((x + ut)^2 + (y + vt)^2), so that 332c4762a1bSJed Brown 333c4762a1bSJed Brown h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 334c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 335c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 336c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 337c4762a1bSJed Brown = 0 338c4762a1bSJed Brown 339c4762a1bSJed Brown */ 340*9371c9d4SSatish Balay static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 341c4762a1bSJed Brown u[0] = x[0] * x[0] + x[1] * x[1]; 342c4762a1bSJed Brown u[1] = 2.0 * x[0] * x[0] - 2.0 * x[0] * x[1]; 343c4762a1bSJed Brown return 0; 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 346c4762a1bSJed Brown /* 347c4762a1bSJed Brown In 2D we use the exact, periodic solution: 348c4762a1bSJed Brown 349c4762a1bSJed Brown u = sin(2 pi x)/4 pi^2 350c4762a1bSJed Brown v = -y cos(2 pi x)/2 pi 351c4762a1bSJed Brown phi = h(x + y + (u + v) t) 352c4762a1bSJed Brown f_x = -sin(2 pi x) 353c4762a1bSJed Brown f_y = 2 pi y cos(2 pi x) 354c4762a1bSJed Brown 355c4762a1bSJed Brown so that 356c4762a1bSJed Brown 357c4762a1bSJed Brown -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 358c4762a1bSJed Brown 359c4762a1bSJed Brown We will conserve phi since 360c4762a1bSJed Brown 361c4762a1bSJed Brown \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 362c4762a1bSJed Brown */ 363*9371c9d4SSatish Balay static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 364c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]) / PetscSqr(2.0 * PETSC_PI); 365c4762a1bSJed Brown u[1] = -x[1] * PetscCosReal(2.0 * PETSC_PI * x[0]) / (2.0 * PETSC_PI); 366c4762a1bSJed Brown return 0; 367c4762a1bSJed Brown } 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* 370c4762a1bSJed Brown In 2D we use the exact, doubly periodic solution: 371c4762a1bSJed Brown 372c4762a1bSJed Brown u = sin(2 pi x) cos(2 pi y)/4 pi^2 373c4762a1bSJed Brown v = -sin(2 pi y) cos(2 pi x)/4 pi^2 374c4762a1bSJed Brown phi = h(x + y + (u + v) t) 375c4762a1bSJed Brown f_x = -2sin(2 pi x) cos(2 pi y) 376c4762a1bSJed Brown f_y = 2sin(2 pi y) cos(2 pi x) 377c4762a1bSJed Brown 378c4762a1bSJed Brown so that 379c4762a1bSJed Brown 380c4762a1bSJed 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 381c4762a1bSJed Brown 382c4762a1bSJed Brown We will conserve phi since 383c4762a1bSJed Brown 384c4762a1bSJed Brown \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 385c4762a1bSJed Brown */ 386*9371c9d4SSatish Balay static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 387c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]) * PetscCosReal(2.0 * PETSC_PI * x[1]) / PetscSqr(2.0 * PETSC_PI); 388c4762a1bSJed Brown u[1] = -PetscSinReal(2.0 * PETSC_PI * x[1]) * PetscCosReal(2.0 * PETSC_PI * x[0]) / PetscSqr(2.0 * PETSC_PI); 389c4762a1bSJed Brown return 0; 390c4762a1bSJed Brown } 391c4762a1bSJed Brown 392*9371c9d4SSatish Balay static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 393c4762a1bSJed Brown u[0] = x[1] - 0.5; 394c4762a1bSJed Brown u[1] = 0.0; 395c4762a1bSJed Brown return 0; 396c4762a1bSJed Brown } 397c4762a1bSJed Brown 398*9371c9d4SSatish Balay static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 399c4762a1bSJed Brown PetscInt d; 400c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 401c4762a1bSJed Brown return 0; 402c4762a1bSJed Brown } 403c4762a1bSJed Brown 404*9371c9d4SSatish Balay static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 405c4762a1bSJed Brown u[0] = 0.0; 406c4762a1bSJed Brown return 0; 407c4762a1bSJed Brown } 408c4762a1bSJed Brown 409*9371c9d4SSatish Balay static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 410c4762a1bSJed Brown u[0] = 1.0; 411c4762a1bSJed Brown return 0; 412c4762a1bSJed Brown } 413c4762a1bSJed Brown 414*9371c9d4SSatish Balay static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 415c4762a1bSJed Brown PetscReal x0[2]; 416c4762a1bSJed Brown PetscScalar xn[2]; 417c4762a1bSJed Brown 418c4762a1bSJed Brown x0[0] = globalUser->source[0]; 419c4762a1bSJed Brown x0[1] = globalUser->source[1]; 420c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 421c4762a1bSJed Brown { 422c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 423c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 424c4762a1bSJed Brown const PetscReal r2 = xi * xi + eta * eta; 425c4762a1bSJed Brown 426c4762a1bSJed Brown u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 427c4762a1bSJed Brown } 428c4762a1bSJed Brown return 0; 429c4762a1bSJed Brown } 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* 432c4762a1bSJed Brown Gaussian blob, initially centered on (0.5, 0.5) 433c4762a1bSJed Brown 434c4762a1bSJed Brown xi = x(t) - x0, eta = y(t) - y0 435c4762a1bSJed Brown 436c4762a1bSJed Brown where x(t), y(t) are the integral curves of v(t), 437c4762a1bSJed Brown 438c4762a1bSJed Brown dx/dt . grad f = v . f 439c4762a1bSJed Brown 440c4762a1bSJed Brown Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 441c4762a1bSJed Brown 442c4762a1bSJed Brown v0 f_x + w0 f_y = v . f 443c4762a1bSJed Brown */ 444*9371c9d4SSatish Balay static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 445c4762a1bSJed Brown const PetscReal x0[2] = {0.5, 0.5}; 446c4762a1bSJed Brown const PetscReal sigma = 1.0 / 6.0; 447c4762a1bSJed Brown PetscScalar xn[2]; 448c4762a1bSJed Brown 449c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 450c4762a1bSJed Brown { 451c4762a1bSJed Brown /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 452c4762a1bSJed Brown /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 453c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 454c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 455c4762a1bSJed Brown const PetscReal r2 = xi * xi + eta * eta; 456c4762a1bSJed Brown 457c4762a1bSJed Brown u[0] = PetscExpReal(-r2 / (2.0 * sigma * sigma)) / (sigma * PetscSqrtReal(2.0 * PETSC_PI)); 458c4762a1bSJed Brown } 459c4762a1bSJed Brown return 0; 460c4762a1bSJed Brown } 461c4762a1bSJed Brown 462*9371c9d4SSatish Balay static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 463c4762a1bSJed Brown PetscReal x0[3]; 464c4762a1bSJed Brown const PetscReal wind[3] = {0.0, 1.0, 0.0}; 465c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 466c4762a1bSJed Brown 467c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 468c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0 * x[0] + 3.0 * x[1]; 469c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 470c4762a1bSJed Brown return 0; 471c4762a1bSJed Brown } 472c4762a1bSJed Brown 473*9371c9d4SSatish Balay static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 474c4762a1bSJed Brown PetscReal ur[3]; 475c4762a1bSJed Brown PetscReal x0[3]; 476c4762a1bSJed Brown const PetscReal t = *((PetscReal *)ctx); 477c4762a1bSJed Brown 478*9371c9d4SSatish Balay ur[0] = PetscRealPart(u[0]); 479*9371c9d4SSatish Balay ur[1] = PetscRealPart(u[1]); 480*9371c9d4SSatish Balay ur[2] = PetscRealPart(u[2]); 481c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 482c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0 * x[0] + 3.0 * x[1]; 483c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 484c4762a1bSJed Brown return 0; 485c4762a1bSJed Brown } 486c4762a1bSJed Brown 487*9371c9d4SSatish Balay static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) { 488c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 489c4762a1bSJed Brown 490c4762a1bSJed Brown PetscFunctionBeginUser; 491c4762a1bSJed Brown xG[0] = user->inflowState; 492c4762a1bSJed Brown PetscFunctionReturn(0); 493c4762a1bSJed Brown } 494c4762a1bSJed Brown 495*9371c9d4SSatish Balay static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) { 496c4762a1bSJed Brown PetscFunctionBeginUser; 49730602db0SMatthew G. Knepley //xG[0] = xI[dim]; 49830602db0SMatthew G. Knepley xG[0] = xI[2]; 499c4762a1bSJed Brown PetscFunctionReturn(0); 500c4762a1bSJed Brown } 501c4762a1bSJed Brown 502*9371c9d4SSatish Balay static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) { 503c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 504c4762a1bSJed Brown PetscInt dim; 505c4762a1bSJed Brown 506c4762a1bSJed Brown PetscFunctionBeginUser; 5079566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 508c4762a1bSJed Brown switch (user->porosityDist) { 509c4762a1bSJed Brown case TILTED: 510c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) tilted_phi_2d(dim, time, x, 2, u, (void *)&time); 511c4762a1bSJed Brown else tilted_phi_coupled_2d(dim, time, x, 2, u, (void *)&time); 512c4762a1bSJed Brown break; 513*9371c9d4SSatish Balay case GAUSSIAN: gaussian_phi_2d(dim, time, x, 2, u, (void *)&time); break; 514*9371c9d4SSatish Balay case DELTA: delta_phi_2d(dim, time, x, 2, u, (void *)&time); break; 515c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 516c4762a1bSJed Brown } 517c4762a1bSJed Brown PetscFunctionReturn(0); 518c4762a1bSJed Brown } 519c4762a1bSJed Brown 520*9371c9d4SSatish Balay static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) { 521c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 522c4762a1bSJed Brown PetscScalar yexact[3] = {0, 0, 0}; 523c4762a1bSJed Brown 524c4762a1bSJed Brown PetscFunctionBeginUser; 5259566063dSJacob Faibussowitsch PetscCall(ExactSolution(dm, time, x, yexact, ctx)); 526c4762a1bSJed Brown f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 527c4762a1bSJed Brown PetscFunctionReturn(0); 528c4762a1bSJed Brown } 529c4762a1bSJed Brown 530*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 531c4762a1bSJed Brown PetscFunctionBeginUser; 5329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 5339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 5349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 5359566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 536c4762a1bSJed Brown PetscFunctionReturn(0); 537c4762a1bSJed Brown } 538c4762a1bSJed Brown 539*9371c9d4SSatish Balay static PetscErrorCode SetupBC(DM dm, AppCtx *user) { 540348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 54130602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 542c4762a1bSJed Brown PetscDS prob; 543c4762a1bSJed Brown DMLabel label; 544c4762a1bSJed Brown PetscBool check; 54530602db0SMatthew G. Knepley PetscInt dim, n = 3; 54630602db0SMatthew G. Knepley const char *prefix; 547c4762a1bSJed Brown 548c4762a1bSJed Brown PetscFunctionBeginUser; 5499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 5509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *)bdt, &n, NULL)); 5519566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 552c4762a1bSJed Brown /* Set initial guesses and exact solutions */ 55330602db0SMatthew G. Knepley switch (dim) { 554c4762a1bSJed Brown case 2: 555c4762a1bSJed Brown user->initialGuess[0] = initialVelocity; 556c4762a1bSJed Brown switch (user->porosityDist) { 557c4762a1bSJed Brown case ZERO: user->initialGuess[1] = zero_phi; break; 558c4762a1bSJed Brown case CONSTANT: user->initialGuess[1] = constant_phi; break; 559c4762a1bSJed Brown case GAUSSIAN: user->initialGuess[1] = gaussian_phi_2d; break; 560c4762a1bSJed Brown case DELTA: user->initialGuess[1] = delta_phi_2d; break; 561c4762a1bSJed Brown case TILTED: 562c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 563c4762a1bSJed Brown else user->initialGuess[1] = tilted_phi_coupled_2d; 564c4762a1bSJed Brown break; 565c4762a1bSJed Brown } 566348a1646SMatthew G. Knepley break; 56763a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 568348a1646SMatthew G. Knepley } 569348a1646SMatthew G. Knepley exactFuncs[0] = user->initialGuess[0]; 570348a1646SMatthew G. Knepley exactFuncs[1] = user->initialGuess[1]; 57130602db0SMatthew G. Knepley switch (dim) { 572348a1646SMatthew G. Knepley case 2: 573c4762a1bSJed Brown switch (user->velocityDist) { 574*9371c9d4SSatish Balay case VEL_ZERO: exactFuncs[0] = zero_u_2d; break; 575*9371c9d4SSatish Balay case VEL_CONSTANT: exactFuncs[0] = constant_u_2d; break; 576c4762a1bSJed Brown case VEL_HARMONIC: 57730602db0SMatthew G. Knepley switch (bdt[0]) { 578c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 57930602db0SMatthew G. Knepley switch (bdt[1]) { 580*9371c9d4SSatish Balay case DM_BOUNDARY_PERIODIC: exactFuncs[0] = doubly_periodic_u_2d; break; 581*9371c9d4SSatish Balay default: exactFuncs[0] = periodic_u_2d; break; 582c4762a1bSJed Brown } 583c4762a1bSJed Brown break; 584*9371c9d4SSatish Balay default: exactFuncs[0] = quadratic_u_2d; break; 585c4762a1bSJed Brown } 586c4762a1bSJed Brown break; 587*9371c9d4SSatish Balay case VEL_SHEAR: exactFuncs[0] = shear_bc; break; 58863a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 589c4762a1bSJed Brown } 590348a1646SMatthew G. Knepley break; 59163a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", dim); 592c4762a1bSJed Brown } 593c4762a1bSJed Brown { 594c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 595c4762a1bSJed Brown 5969566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, "", "-use_implicit", &isImplicit)); 597348a1646SMatthew G. Knepley if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = exactFuncs[0]; 598c4762a1bSJed Brown } 5999566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-dmts_check", &check)); 600c4762a1bSJed Brown if (check) { 601348a1646SMatthew G. Knepley user->initialGuess[0] = exactFuncs[0]; 602348a1646SMatthew G. Knepley user->initialGuess[1] = exactFuncs[1]; 603c4762a1bSJed Brown } 604c4762a1bSJed Brown /* Set BC */ 6059566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 6069566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 6079566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], user)); 6089566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], user)); 609c4762a1bSJed Brown if (label) { 610c4762a1bSJed Brown const PetscInt id = 1; 611c4762a1bSJed Brown 6129566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, user, NULL)); 613c4762a1bSJed Brown } 6149566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &label)); 615c4762a1bSJed Brown if (label && user->useFV) { 616c4762a1bSJed Brown const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101}; 617c4762a1bSJed Brown 618dd39110bSPierre 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)); 619dd39110bSPierre 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)); 620c4762a1bSJed Brown } 621c4762a1bSJed Brown PetscFunctionReturn(0); 622c4762a1bSJed Brown } 623c4762a1bSJed Brown 624*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *user) { 62530602db0SMatthew G. Knepley DMBoundaryType bdt[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 626c4762a1bSJed Brown PetscDS prob; 62730602db0SMatthew G. Knepley PetscInt n = 3; 62830602db0SMatthew G. Knepley const char *prefix; 629c4762a1bSJed Brown 630c4762a1bSJed Brown PetscFunctionBeginUser; 6319566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 6329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnumArray(NULL, prefix, "-dm_plex_box_bd", DMBoundaryTypes, (PetscEnum *)bdt, &n, NULL)); 6339566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 634c4762a1bSJed Brown switch (user->velocityDist) { 635*9371c9d4SSatish Balay case VEL_ZERO: PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u)); break; 636c4762a1bSJed Brown case VEL_CONSTANT: 6379566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u)); 6389566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL)); 6399566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL)); 640c4762a1bSJed Brown break; 641c4762a1bSJed Brown case VEL_HARMONIC: 64230602db0SMatthew G. Knepley switch (bdt[0]) { 643c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 64430602db0SMatthew G. Knepley switch (bdt[1]) { 645*9371c9d4SSatish Balay case DM_BOUNDARY_PERIODIC: PetscCall(PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u)); break; 646*9371c9d4SSatish Balay default: PetscCall(PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u)); break; 647c4762a1bSJed Brown } 648c4762a1bSJed Brown break; 649*9371c9d4SSatish Balay default: PetscCall(PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u)); break; 650c4762a1bSJed Brown } 6519566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 652c4762a1bSJed Brown break; 653c4762a1bSJed Brown case VEL_SHEAR: 6549566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u)); 6559566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu)); 656c4762a1bSJed Brown break; 657c4762a1bSJed Brown } 6589566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_advection, f1_advection)); 6599566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL)); 6609566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL)); 6619566063dSJacob Faibussowitsch if (user->velocityDist == VEL_ZERO) PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_advection)); 6629566063dSJacob Faibussowitsch else PetscCall(PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection)); 663c4762a1bSJed Brown 6649566063dSJacob Faibussowitsch PetscCall(FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user)); 665c4762a1bSJed Brown PetscFunctionReturn(0); 666c4762a1bSJed Brown } 667c4762a1bSJed Brown 668*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) { 669c4762a1bSJed Brown DM cdm = dm; 670c4762a1bSJed Brown PetscQuadrature q; 671c4762a1bSJed Brown PetscFE fe[2]; 672c4762a1bSJed Brown PetscFV fv; 673c4762a1bSJed Brown MPI_Comm comm; 67430602db0SMatthew G. Knepley PetscInt dim; 675c4762a1bSJed Brown 676c4762a1bSJed Brown PetscFunctionBeginUser; 677c4762a1bSJed Brown /* Create finite element */ 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 6799566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 6809566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0])); 6819566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 6829566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1])); 6839566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 6849566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "porosity")); 685c4762a1bSJed Brown 6869566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)dm), &fv)); 6879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fv, "porosity")); 6889566063dSJacob Faibussowitsch PetscCall(PetscFVSetFromOptions(fv)); 6899566063dSJacob Faibussowitsch PetscCall(PetscFVSetNumComponents(fv, 1)); 6909566063dSJacob Faibussowitsch PetscCall(PetscFVSetSpatialDimension(fv, dim)); 6919566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe[0], &q)); 6929566063dSJacob Faibussowitsch PetscCall(PetscFVSetQuadrature(fv, q)); 693c4762a1bSJed Brown 6949566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 6959566063dSJacob Faibussowitsch if (user->useFV) PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fv)); 6969566063dSJacob Faibussowitsch else PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 6979566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 6989566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 699c4762a1bSJed Brown 700c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 701c4762a1bSJed Brown while (cdm) { 7029566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 7039566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 704c4762a1bSJed Brown /* Coordinates were never localized for coarse meshes */ 7059566063dSJacob Faibussowitsch if (cdm) PetscCall(DMLocalizeCoordinates(cdm)); 706c4762a1bSJed Brown } 7079566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 7089566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 7099566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fv)); 710c4762a1bSJed Brown PetscFunctionReturn(0); 711c4762a1bSJed Brown } 712c4762a1bSJed Brown 713*9371c9d4SSatish Balay static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) { 714c4762a1bSJed Brown PetscFunctionBeginUser; 7159566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, user, dm)); 716c4762a1bSJed Brown /* Handle refinement, etc. */ 7179566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 718c4762a1bSJed Brown /* Construct ghost cells */ 719c4762a1bSJed Brown if (user->useFV) { 720c4762a1bSJed Brown DM gdm; 721c4762a1bSJed Brown 7229566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm)); 7239566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 724c4762a1bSJed Brown *dm = gdm; 725c4762a1bSJed Brown } 726c4762a1bSJed Brown /* Localize coordinates */ 7279566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 7289566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)(*dm), "Mesh")); 7299566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 730c4762a1bSJed Brown /* Setup problem */ 7319566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*dm, user)); 732c4762a1bSJed Brown /* Setup BC */ 7339566063dSJacob Faibussowitsch PetscCall(SetupBC(*dm, user)); 734c4762a1bSJed Brown PetscFunctionReturn(0); 735c4762a1bSJed Brown } 736c4762a1bSJed Brown 737*9371c9d4SSatish Balay static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx) { 738c4762a1bSJed Brown PetscDS prob; 739c4762a1bSJed Brown DM dmCell; 740c4762a1bSJed Brown Vec cellgeom; 741c4762a1bSJed Brown const PetscScalar *cgeom; 742c4762a1bSJed Brown PetscScalar *x; 743c4762a1bSJed Brown PetscInt dim, Nf, cStart, cEnd, c; 744c4762a1bSJed Brown 745c4762a1bSJed Brown PetscFunctionBeginUser; 7469566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 7479566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7489566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 7499566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 7509566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 7519566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 7529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 7539566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 754c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 755c4762a1bSJed Brown PetscFVCellGeom *cg; 756c4762a1bSJed Brown PetscScalar *xc; 757c4762a1bSJed Brown 7589566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 7599566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, field, x, &xc)); 7609566063dSJacob Faibussowitsch if (xc) PetscCall((*func)(dim, 0.0, cg->centroid, Nf, xc, ctx)); 761c4762a1bSJed Brown } 7629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 7639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 764c4762a1bSJed Brown PetscFunctionReturn(0); 765c4762a1bSJed Brown } 766c4762a1bSJed Brown 767*9371c9d4SSatish Balay static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) { 768c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 769c4762a1bSJed Brown char *ftable = NULL; 770c4762a1bSJed Brown DM dm; 771c4762a1bSJed Brown PetscSection s; 772c4762a1bSJed Brown Vec cellgeom; 773c4762a1bSJed Brown const PetscScalar *x; 774c4762a1bSJed Brown PetscScalar *a; 775c4762a1bSJed Brown PetscReal *xnorms; 776c4762a1bSJed Brown PetscInt pStart, pEnd, p, Nf, f; 777c4762a1bSJed Brown 778c4762a1bSJed Brown PetscFunctionBeginUser; 7799566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X, (PetscObject)ts, "-view_solution")); 7809566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 7819566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dm, NULL, &cellgeom, NULL)); 7829566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 7839566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &Nf)); 7849566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 7859566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nf * 2, &xnorms)); 7869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 787c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 788c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 789c4762a1bSJed Brown PetscInt dof, cdof, d; 790c4762a1bSJed Brown 7919566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 7929566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 7939566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, p, f, x, &a)); 794c4762a1bSJed Brown /* TODO Use constrained indices here */ 795c4762a1bSJed Brown for (d = 0; d < dof - cdof; ++d) xnorms[f * 2 + 0] = PetscMax(xnorms[f * 2 + 0], PetscAbsScalar(a[d])); 796c4762a1bSJed Brown for (d = 0; d < dof - cdof; ++d) xnorms[f * 2 + 1] += PetscAbsScalar(a[d]); 797c4762a1bSJed Brown } 798c4762a1bSJed Brown } 7999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 800c4762a1bSJed Brown if (stepnum >= 0) { /* No summary for final time */ 801c4762a1bSJed Brown DM dmCell, *fdm; 802c4762a1bSJed Brown Vec *fv; 803c4762a1bSJed Brown const PetscScalar *cgeom; 804c4762a1bSJed Brown PetscScalar **fx; 805c4762a1bSJed Brown PetscReal *fmin, *fmax, *fint, *ftmp, t; 806c4762a1bSJed Brown PetscInt cStart, cEnd, c, fcount, f, num; 807c4762a1bSJed Brown 808c4762a1bSJed Brown size_t ftableused, ftablealloc; 809c4762a1bSJed Brown 810c4762a1bSJed Brown /* Functionals have indices after registering, this is an upper bound */ 811c4762a1bSJed Brown fcount = user->numMonitorFuncs; 8129566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fint, fcount, &ftmp)); 8139566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fcount, &fdm, fcount, &fv, fcount, &fx)); 814c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 815c4762a1bSJed Brown PetscSection fs; 816c4762a1bSJed Brown const char *name = user->monitorFuncs[f]->name; 817c4762a1bSJed Brown 818c4762a1bSJed Brown fmin[f] = PETSC_MAX_REAL; 819c4762a1bSJed Brown fmax[f] = PETSC_MIN_REAL; 820c4762a1bSJed Brown fint[f] = 0; 821c4762a1bSJed Brown /* Make monitor vecs */ 8229566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &fdm[f])); 8239566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &num, &t)); 8249566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(fdm[f], num, t)); 8259566063dSJacob Faibussowitsch PetscCall(PetscSectionClone(s, &fs)); 8269566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 0, NULL)); 8279566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(fs, 1, name)); 8289566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(fdm[f], fs)); 8299566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&fs)); 8309566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm[f], &fv[f])); 8319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fv[f], name)); 8329566063dSJacob Faibussowitsch PetscCall(VecGetArray(fv[f], &fx[f])); 833c4762a1bSJed Brown } 8349566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 8359566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 8369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 8379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 838c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 839c4762a1bSJed Brown PetscFVCellGeom *cg; 840c4762a1bSJed Brown PetscScalar *cx; 841c4762a1bSJed Brown 8429566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 8439566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx)); 844c4762a1bSJed Brown if (!cx) continue; /* not a global cell */ 845c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 846c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 847c4762a1bSJed Brown PetscScalar *fxc; 848c4762a1bSJed Brown 8499566063dSJacob Faibussowitsch PetscCall(DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc)); 850c4762a1bSJed Brown /* I need to make it easier to get interpolated values here */ 8519566063dSJacob Faibussowitsch PetscCall((*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx)); 852c4762a1bSJed Brown fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 853c4762a1bSJed Brown } 854c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 855c4762a1bSJed Brown fmin[f] = PetscMin(fmin[f], ftmp[f]); 856c4762a1bSJed Brown fmax[f] = PetscMax(fmax[f], ftmp[f]); 857c4762a1bSJed Brown fint[f] += cg->volume * ftmp[f]; 858c4762a1bSJed Brown } 859c4762a1bSJed Brown } 8609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 8619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 8629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 8639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); 8649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 865c4762a1bSJed Brown /* Output functional data */ 866c4762a1bSJed Brown ftablealloc = fcount * 100; 867c4762a1bSJed Brown ftableused = 0; 8689566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ftablealloc, &ftable)); 869c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 870c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 871c4762a1bSJed Brown PetscInt id = func->offset; 872c4762a1bSJed Brown char newline[] = "\n"; 873c4762a1bSJed Brown char buffer[256], *p, *prefix; 874c4762a1bSJed Brown size_t countused, len; 875c4762a1bSJed Brown 876c4762a1bSJed Brown /* Create string with functional outputs */ 877c4762a1bSJed Brown if (f % 3) { 8789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, " ", 2)); 879c4762a1bSJed Brown p = buffer + 2; 880c4762a1bSJed Brown } else if (f) { 8819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(buffer, newline, sizeof(newline) - 1)); 882c4762a1bSJed Brown p = buffer + sizeof(newline) - 1; 883c4762a1bSJed Brown } else { 884c4762a1bSJed Brown p = buffer; 885c4762a1bSJed Brown } 8869566063dSJacob 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])); 887c4762a1bSJed Brown countused += p - buffer; 888c4762a1bSJed Brown /* reallocate */ 889c4762a1bSJed Brown if (countused > ftablealloc - ftableused - 1) { 890c4762a1bSJed Brown char *ftablenew; 891c4762a1bSJed Brown 892c4762a1bSJed Brown ftablealloc = 2 * ftablealloc + countused; 8939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ftablealloc, &ftablenew)); 8949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); 8959566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 896c4762a1bSJed Brown ftable = ftablenew; 897c4762a1bSJed Brown } 8989566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused)); 899c4762a1bSJed Brown ftableused += countused; 900c4762a1bSJed Brown ftable[ftableused] = 0; 901c4762a1bSJed Brown /* Output vecs */ 9029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(fv[f], &fx[f])); 9039566063dSJacob Faibussowitsch PetscCall(PetscStrlen(func->name, &len)); 9049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len + 2, &prefix)); 9059566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(prefix, func->name)); 9069566063dSJacob Faibussowitsch PetscCall(PetscStrcat(prefix, "_")); 9079566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix)); 9089566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fv[f], NULL, "-vec_view")); 9099566063dSJacob Faibussowitsch PetscCall(PetscFree(prefix)); 9109566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm[f], &fv[f])); 9119566063dSJacob Faibussowitsch PetscCall(DMDestroy(&fdm[f])); 912c4762a1bSJed Brown } 9139566063dSJacob Faibussowitsch PetscCall(PetscFree4(fmin, fmax, fint, ftmp)); 9149566063dSJacob Faibussowitsch PetscCall(PetscFree3(fdm, fv, fx)); 91563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT " time %8.4g |x| (", stepnum, (double)time)); 916c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 9179566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ", ")); 9189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%8.4g", (double)xnorms[f * 2 + 0])); 919c4762a1bSJed Brown } 9209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ") |x|_1 (")); 921c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 9229566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ", ")); 9239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%8.4g", (double)xnorms[f * 2 + 1])); 924c4762a1bSJed Brown } 9259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), ") %s\n", ftable ? ftable : "")); 9269566063dSJacob Faibussowitsch PetscCall(PetscFree(ftable)); 927c4762a1bSJed Brown } 9289566063dSJacob Faibussowitsch PetscCall(PetscFree(xnorms)); 929c4762a1bSJed Brown PetscFunctionReturn(0); 930c4762a1bSJed Brown } 931c4762a1bSJed Brown 932*9371c9d4SSatish Balay int main(int argc, char **argv) { 933c4762a1bSJed Brown MPI_Comm comm; 934c4762a1bSJed Brown TS ts; 935c4762a1bSJed Brown DM dm; 936c4762a1bSJed Brown Vec u; 937c4762a1bSJed Brown AppCtx user; 938c4762a1bSJed Brown PetscReal t0, t = 0.0; 939c4762a1bSJed Brown void *ctxs[2]; 940c4762a1bSJed Brown 941c4762a1bSJed Brown ctxs[0] = &t; 942c4762a1bSJed Brown ctxs[1] = &t; 943327415f7SBarry Smith PetscFunctionBeginUser; 9449566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 945c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 946c4762a1bSJed Brown user.functionalRegistry = NULL; 947c4762a1bSJed Brown globalUser = &user; 9489566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 9499566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 9509566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 9519566063dSJacob Faibussowitsch PetscCall(CreateDM(comm, &user, &dm)); 9529566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 9539566063dSJacob Faibussowitsch PetscCall(ProcessMonitorOptions(comm, &user)); 954c4762a1bSJed Brown 9559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 9569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 957c4762a1bSJed Brown if (user.useFV) { 958c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 959c4762a1bSJed Brown 9609566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, "", "-use_implicit", &isImplicit)); 961c4762a1bSJed Brown if (isImplicit) { 9629566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 9639566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 964c4762a1bSJed Brown } 9659566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 9669566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user)); 967c4762a1bSJed Brown } else { 9689566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 9699566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 9709566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 971c4762a1bSJed Brown } 9729566063dSJacob Faibussowitsch if (user.useFV) PetscCall(TSMonitorSet(ts, MonitorFunctionals, &user, NULL)); 9739566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1)); 9749566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2.0)); 9759566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01)); 9769566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 9779566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 978c4762a1bSJed Brown 9799566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u)); 9809566063dSJacob Faibussowitsch if (user.useFV) PetscCall(SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1])); 9819566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-init_vec_view")); 9829566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 983c4762a1bSJed Brown t0 = t; 9849566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 9859566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 9869566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 9879566063dSJacob Faibussowitsch if (t > t0) PetscCall(DMTSCheckFromOptions(ts, u)); 9889566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 989c4762a1bSJed Brown { 990c4762a1bSJed Brown PetscReal ftime; 991c4762a1bSJed Brown PetscInt nsteps; 992c4762a1bSJed Brown TSConvergedReason reason; 993c4762a1bSJed Brown 9949566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 9959566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &nsteps)); 9969566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 99763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps)); 998c4762a1bSJed Brown } 999c4762a1bSJed Brown 10009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 10019566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 10029566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 10039566063dSJacob Faibussowitsch PetscCall(PetscFree(user.monitorFuncs)); 10049566063dSJacob Faibussowitsch PetscCall(FunctionalDestroy(&user.functionalRegistry)); 10059566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1006b122ec5aSJacob Faibussowitsch return 0; 1007c4762a1bSJed Brown } 1008c4762a1bSJed Brown 1009c4762a1bSJed Brown /*TEST 1010c4762a1bSJed Brown 101130602db0SMatthew G. Knepley testset: 101230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3,3 101330602db0SMatthew G. Knepley 1014c4762a1bSJed Brown # 2D harmonic velocity, no porosity 1015c4762a1bSJed Brown test: 1016c4762a1bSJed Brown suffix: p1p1 1017c4762a1bSJed Brown requires: !complex !single 101830602db0SMatthew 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 1019c4762a1bSJed Brown test: 1020c4762a1bSJed Brown suffix: p1p1_xper 1021c4762a1bSJed Brown requires: !complex !single 102230602db0SMatthew 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 1023c4762a1bSJed Brown test: 1024c4762a1bSJed Brown suffix: p1p1_xper_ref 1025c4762a1bSJed Brown requires: !complex !single 102630602db0SMatthew 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 1027c4762a1bSJed Brown test: 1028c4762a1bSJed Brown suffix: p1p1_xyper 1029c4762a1bSJed Brown requires: !complex !single 103030602db0SMatthew 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 1031c4762a1bSJed Brown test: 1032c4762a1bSJed Brown suffix: p1p1_xyper_ref 1033c4762a1bSJed Brown requires: !complex !single 103430602db0SMatthew 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 1035c4762a1bSJed Brown test: 1036c4762a1bSJed Brown suffix: p2p1 1037c4762a1bSJed Brown requires: !complex !single 103830602db0SMatthew 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 1039c4762a1bSJed Brown test: 1040c4762a1bSJed Brown suffix: p2p1_xyper 1041c4762a1bSJed Brown requires: !complex !single 104230602db0SMatthew 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 104330602db0SMatthew G. Knepley 104430602db0SMatthew G. Knepley test: 104530602db0SMatthew G. Knepley suffix: adv_1 104630602db0SMatthew G. Knepley requires: !complex !single 104730602db0SMatthew 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 104830602db0SMatthew G. Knepley 104930602db0SMatthew G. Knepley test: 105030602db0SMatthew G. Knepley suffix: adv_2 105130602db0SMatthew G. Knepley requires: !complex 105230602db0SMatthew G. Knepley TODO: broken memory corruption 105330602db0SMatthew 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 105430602db0SMatthew G. Knepley 105530602db0SMatthew G. Knepley test: 105630602db0SMatthew G. Knepley suffix: adv_3 105730602db0SMatthew G. Knepley requires: !complex 105830602db0SMatthew G. Knepley TODO: broken memory corruption 105930602db0SMatthew 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 106030602db0SMatthew G. Knepley 106130602db0SMatthew G. Knepley test: 106230602db0SMatthew G. Knepley suffix: adv_3_ex 106330602db0SMatthew G. Knepley requires: !complex 106430602db0SMatthew 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 106530602db0SMatthew G. Knepley 106630602db0SMatthew G. Knepley test: 106730602db0SMatthew G. Knepley suffix: adv_4 106830602db0SMatthew G. Knepley requires: !complex 106930602db0SMatthew G. Knepley TODO: broken memory corruption 107030602db0SMatthew 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 107130602db0SMatthew G. Knepley 107230602db0SMatthew G. Knepley # 2D Advection, box, delta 107330602db0SMatthew G. Knepley test: 107430602db0SMatthew G. Knepley suffix: adv_delta_yper_0 107530602db0SMatthew G. Knepley requires: !complex 107630602db0SMatthew G. Knepley TODO: broken 107730602db0SMatthew 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 107830602db0SMatthew G. Knepley 107930602db0SMatthew G. Knepley test: 108030602db0SMatthew G. Knepley suffix: adv_delta_yper_1 108130602db0SMatthew G. Knepley requires: !complex 108230602db0SMatthew G. Knepley TODO: broken 108330602db0SMatthew 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 108430602db0SMatthew G. Knepley 108530602db0SMatthew G. Knepley test: 108630602db0SMatthew G. Knepley suffix: adv_delta_yper_2 108730602db0SMatthew G. Knepley requires: !complex 108830602db0SMatthew G. Knepley TODO: broken 108930602db0SMatthew 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 109030602db0SMatthew G. Knepley 109130602db0SMatthew G. Knepley test: 109230602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_0 109330602db0SMatthew G. Knepley requires: !complex 109430602db0SMatthew G. Knepley TODO: broken 109530602db0SMatthew 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 109630602db0SMatthew G. Knepley 109730602db0SMatthew G. Knepley test: 109830602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_1 109930602db0SMatthew G. Knepley requires: !complex 110030602db0SMatthew G. Knepley TODO: broken 110130602db0SMatthew 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 110230602db0SMatthew G. Knepley 110330602db0SMatthew G. Knepley test: 110430602db0SMatthew G. Knepley suffix: adv_delta_yper_fim_2 110530602db0SMatthew G. Knepley requires: !complex 110630602db0SMatthew G. Knepley TODO: broken 110730602db0SMatthew 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 110830602db0SMatthew G. Knepley 110930602db0SMatthew G. Knepley test: 111030602db0SMatthew G. Knepley suffix: adv_delta_yper_im_0 111130602db0SMatthew G. Knepley requires: !complex 111230602db0SMatthew G. Knepley TODO: broken 111330602db0SMatthew 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 111430602db0SMatthew G. Knepley 111530602db0SMatthew G. Knepley test: 111630602db0SMatthew G. Knepley suffix: adv_delta_yper_im_1 111730602db0SMatthew G. Knepley requires: !complex 111830602db0SMatthew G. Knepley TODO: broken 111930602db0SMatthew 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 112030602db0SMatthew G. Knepley 112130602db0SMatthew G. Knepley test: 112230602db0SMatthew G. Knepley suffix: adv_delta_yper_im_2 112330602db0SMatthew G. Knepley requires: !complex 112430602db0SMatthew G. Knepley TODO: broken 112530602db0SMatthew 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 112630602db0SMatthew G. Knepley 112730602db0SMatthew G. Knepley test: 112830602db0SMatthew G. Knepley suffix: adv_delta_yper_im_3 112930602db0SMatthew G. Knepley requires: !complex 113030602db0SMatthew G. Knepley TODO: broken 113130602db0SMatthew 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 113230602db0SMatthew G. Knepley 113330602db0SMatthew G. Knepley # I believe the nullspace is sin(pi y) 113430602db0SMatthew G. Knepley test: 113530602db0SMatthew G. Knepley suffix: adv_delta_yper_im_4 113630602db0SMatthew G. Knepley requires: !complex 113730602db0SMatthew G. Knepley TODO: broken 113830602db0SMatthew 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 113930602db0SMatthew G. Knepley 114030602db0SMatthew G. Knepley test: 114130602db0SMatthew G. Knepley suffix: adv_delta_yper_im_5 114230602db0SMatthew G. Knepley requires: !complex 114330602db0SMatthew G. Knepley TODO: broken 114430602db0SMatthew 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 114530602db0SMatthew G. Knepley 114630602db0SMatthew G. Knepley test: 114730602db0SMatthew G. Knepley suffix: adv_delta_yper_im_6 114830602db0SMatthew G. Knepley requires: !complex 114930602db0SMatthew G. Knepley TODO: broken 115030602db0SMatthew 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 115130602db0SMatthew G. Knepley # 2D Advection, magma benchmark 1 115230602db0SMatthew G. Knepley 115330602db0SMatthew G. Knepley test: 115430602db0SMatthew G. Knepley suffix: adv_delta_shear_im_0 115530602db0SMatthew G. Knepley requires: !complex 115630602db0SMatthew G. Knepley TODO: broken 115730602db0SMatthew 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 115830602db0SMatthew G. Knepley # 2D Advection, box, gaussian 115930602db0SMatthew G. Knepley 116030602db0SMatthew G. Knepley test: 116130602db0SMatthew G. Knepley suffix: adv_gauss 116230602db0SMatthew G. Knepley requires: !complex 116330602db0SMatthew G. Knepley TODO: broken 116430602db0SMatthew 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 116530602db0SMatthew G. Knepley 116630602db0SMatthew G. Knepley test: 116730602db0SMatthew G. Knepley suffix: adv_gauss_im 116830602db0SMatthew G. Knepley requires: !complex 116930602db0SMatthew G. Knepley TODO: broken 117030602db0SMatthew 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 117130602db0SMatthew G. Knepley 117230602db0SMatthew G. Knepley test: 117330602db0SMatthew G. Knepley suffix: adv_gauss_im_1 117430602db0SMatthew G. Knepley requires: !complex 117530602db0SMatthew G. Knepley TODO: broken 117630602db0SMatthew 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 117730602db0SMatthew G. Knepley 117830602db0SMatthew G. Knepley test: 117930602db0SMatthew G. Knepley suffix: adv_gauss_im_2 118030602db0SMatthew G. Knepley requires: !complex 118130602db0SMatthew G. Knepley TODO: broken 118230602db0SMatthew 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 118330602db0SMatthew G. Knepley 118430602db0SMatthew G. Knepley test: 118530602db0SMatthew G. Knepley suffix: adv_gauss_corner 118630602db0SMatthew G. Knepley requires: !complex 118730602db0SMatthew G. Knepley TODO: broken 118830602db0SMatthew 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 118930602db0SMatthew G. Knepley 119030602db0SMatthew G. Knepley # 2D Advection+Harmonic 12- 119130602db0SMatthew G. Knepley test: 119230602db0SMatthew G. Knepley suffix: adv_harm_0 119330602db0SMatthew G. Knepley requires: !complex 119430602db0SMatthew G. Knepley TODO: broken memory corruption 119530602db0SMatthew 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 119630602db0SMatthew G. Knepley 1197c4762a1bSJed Brown # Must check that FV BCs propagate to coarse meshes 1198c4762a1bSJed Brown # Must check that FV BC ids propagate to coarse meshes 1199c4762a1bSJed Brown # Must check that FE+FV BCs work at the same time 1200c4762a1bSJed Brown # 2D Advection, matching wind in ex11 8-11 1201c4762a1bSJed Brown # NOTE implicit solves are limited by accuracy of FD Jacobian 1202c4762a1bSJed Brown test: 1203c4762a1bSJed Brown suffix: adv_0 1204c4762a1bSJed Brown requires: !complex !single exodusii 120530602db0SMatthew 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 1206c4762a1bSJed Brown 1207c4762a1bSJed Brown test: 1208c4762a1bSJed Brown suffix: adv_0_im 1209c4762a1bSJed Brown requires: !complex exodusii 1210c4762a1bSJed Brown TODO: broken memory corruption 121130602db0SMatthew 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 1212c4762a1bSJed Brown 1213c4762a1bSJed Brown test: 1214c4762a1bSJed Brown suffix: adv_0_im_2 1215c4762a1bSJed Brown requires: !complex exodusii 1216c4762a1bSJed Brown TODO: broken 121730602db0SMatthew 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 1218c4762a1bSJed Brown 1219c4762a1bSJed Brown test: 1220c4762a1bSJed Brown suffix: adv_0_im_3 1221c4762a1bSJed Brown requires: !complex exodusii 1222c4762a1bSJed Brown TODO: broken 122330602db0SMatthew 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 1224c4762a1bSJed Brown 1225c4762a1bSJed Brown test: 1226c4762a1bSJed Brown suffix: adv_0_im_4 1227c4762a1bSJed Brown requires: !complex exodusii 1228c4762a1bSJed Brown TODO: broken 122930602db0SMatthew 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 1230c4762a1bSJed Brown # 2D Advection, misc 1231c4762a1bSJed Brown 1232c4762a1bSJed Brown TEST*/ 1233