1*c4762a1bSJed Brown static char help[] = "Hybrid Finite Element-Finite Volume Example.\n"; 2*c4762a1bSJed Brown /*F 3*c4762a1bSJed Brown Here we are advecting a passive tracer in a harmonic velocity field, defined by 4*c4762a1bSJed Brown a forcing function $f$: 5*c4762a1bSJed Brown \begin{align} 6*c4762a1bSJed Brown -\Delta \mathbf{u} + f &= 0 \\ 7*c4762a1bSJed Brown \frac{\partial\phi}{\partial t} + \nabla\cdot \phi \mathbf{u} &= 0 8*c4762a1bSJed Brown \end{align} 9*c4762a1bSJed Brown F*/ 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown #include <petscdmplex.h> 12*c4762a1bSJed Brown #include <petscds.h> 13*c4762a1bSJed Brown #include <petscts.h> 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h> /* For DotD */ 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown #define ALEN(a) (sizeof(a)/sizeof((a)[0])) 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown typedef enum {VEL_ZERO, VEL_CONSTANT, VEL_HARMONIC, VEL_SHEAR} VelocityDistribution; 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown typedef enum {ZERO, CONSTANT, GAUSSIAN, TILTED, DELTA} PorosityDistribution; 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown /* 26*c4762a1bSJed Brown FunctionalFunc - Calculates the value of a functional of the solution at a point 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown Input Parameters: 29*c4762a1bSJed Brown + dm - The DM 30*c4762a1bSJed Brown . time - The TS time 31*c4762a1bSJed Brown . x - The coordinates of the evaluation point 32*c4762a1bSJed Brown . u - The field values at point x 33*c4762a1bSJed Brown - ctx - A user context, or NULL 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown Output Parameter: 36*c4762a1bSJed Brown . f - The value of the functional at point x 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown */ 39*c4762a1bSJed Brown typedef PetscErrorCode (*FunctionalFunc)(DM, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown typedef struct _n_Functional *Functional; 42*c4762a1bSJed Brown struct _n_Functional { 43*c4762a1bSJed Brown char *name; 44*c4762a1bSJed Brown FunctionalFunc func; 45*c4762a1bSJed Brown void *ctx; 46*c4762a1bSJed Brown PetscInt offset; 47*c4762a1bSJed Brown Functional next; 48*c4762a1bSJed Brown }; 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown typedef struct { 51*c4762a1bSJed Brown /* Domain and mesh definition */ 52*c4762a1bSJed Brown DM dm; 53*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 54*c4762a1bSJed Brown DMBoundaryType bd[2]; /* The boundary type for the x- and y-boundary */ 55*c4762a1bSJed Brown char filename[2048]; /* The optional ExodusII file */ 56*c4762a1bSJed Brown /* Problem definition */ 57*c4762a1bSJed Brown PetscBool useFV; /* Use a finite volume scheme for advection */ 58*c4762a1bSJed Brown PetscErrorCode (*exactFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 59*c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 60*c4762a1bSJed Brown VelocityDistribution velocityDist; 61*c4762a1bSJed Brown PorosityDistribution porosityDist; 62*c4762a1bSJed Brown PetscReal inflowState; 63*c4762a1bSJed Brown PetscReal source[3]; 64*c4762a1bSJed Brown /* Monitoring */ 65*c4762a1bSJed Brown PetscInt numMonitorFuncs, maxMonitorFunc; 66*c4762a1bSJed Brown Functional *monitorFuncs; 67*c4762a1bSJed Brown PetscInt errorFunctional; 68*c4762a1bSJed Brown Functional functionalRegistry; 69*c4762a1bSJed Brown } AppCtx; 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown static AppCtx *globalUser; 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 74*c4762a1bSJed Brown { 75*c4762a1bSJed Brown const char *velocityDist[4] = {"zero", "constant", "harmonic", "shear"}; 76*c4762a1bSJed Brown const char *porosityDist[5] = {"zero", "constant", "gaussian", "tilted", "delta"}; 77*c4762a1bSJed Brown PetscInt bd, vd, pd, d; 78*c4762a1bSJed Brown PetscBool flg; 79*c4762a1bSJed Brown PetscErrorCode ierr; 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown PetscFunctionBeginUser; 82*c4762a1bSJed Brown options->dim = 2; 83*c4762a1bSJed Brown options->bd[0] = DM_BOUNDARY_PERIODIC; 84*c4762a1bSJed Brown options->bd[1] = DM_BOUNDARY_PERIODIC; 85*c4762a1bSJed Brown options->filename[0] = '\0'; 86*c4762a1bSJed Brown options->useFV = PETSC_FALSE; 87*c4762a1bSJed Brown options->velocityDist = VEL_HARMONIC; 88*c4762a1bSJed Brown options->porosityDist = ZERO; 89*c4762a1bSJed Brown options->inflowState = -2.0; 90*c4762a1bSJed Brown options->numMonitorFuncs = 0; 91*c4762a1bSJed Brown options->source[0] = 0.5; 92*c4762a1bSJed Brown options->source[1] = 0.5; 93*c4762a1bSJed Brown options->source[2] = 0.5; 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Magma Dynamics Options", "DMPLEX");CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex18.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 97*c4762a1bSJed Brown bd = options->bd[0]; 98*c4762a1bSJed Brown ierr = PetscOptionsEList("-x_bd_type", "The x-boundary type", "ex18.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->bd[0]], &bd, NULL);CHKERRQ(ierr); 99*c4762a1bSJed Brown options->bd[0] = (DMBoundaryType) bd; 100*c4762a1bSJed Brown bd = options->bd[1]; 101*c4762a1bSJed Brown ierr = PetscOptionsEList("-y_bd_type", "The y-boundary type", "ex18.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->bd[1]], &bd, NULL);CHKERRQ(ierr); 102*c4762a1bSJed Brown options->bd[1] = (DMBoundaryType) bd; 103*c4762a1bSJed Brown ierr = PetscOptionsString("-f", "Exodus.II filename to read", "ex18.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = PetscOptionsBool("-use_fv", "Use the finite volume method for advection", "ex18.c", options->useFV, &options->useFV, NULL);CHKERRQ(ierr); 105*c4762a1bSJed Brown vd = options->velocityDist; 106*c4762a1bSJed Brown ierr = PetscOptionsEList("-velocity_dist","Velocity distribution type","ex18.c",velocityDist,4,velocityDist[options->velocityDist],&vd,NULL);CHKERRQ(ierr); 107*c4762a1bSJed Brown options->velocityDist = (VelocityDistribution) vd; 108*c4762a1bSJed Brown pd = options->porosityDist; 109*c4762a1bSJed Brown ierr = PetscOptionsEList("-porosity_dist","Initial porosity distribution type","ex18.c",porosityDist,5,porosityDist[options->porosityDist],&pd,NULL);CHKERRQ(ierr); 110*c4762a1bSJed Brown options->porosityDist = (PorosityDistribution) pd; 111*c4762a1bSJed Brown ierr = PetscOptionsReal("-inflow_state", "The inflow state", "ex18.c", options->inflowState, &options->inflowState, NULL);CHKERRQ(ierr); 112*c4762a1bSJed Brown d = options->dim; 113*c4762a1bSJed Brown ierr = PetscOptionsRealArray("-source_loc", "The source location", "ex18.c", options->source, &d, &flg);CHKERRQ(ierr); 114*c4762a1bSJed Brown if (flg && d != options->dim) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must give dim coordinates for the source location, not %d", d); 115*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown PetscFunctionReturn(0); 118*c4762a1bSJed Brown } 119*c4762a1bSJed Brown 120*c4762a1bSJed Brown static PetscErrorCode ProcessMonitorOptions(MPI_Comm comm, AppCtx *options) 121*c4762a1bSJed Brown { 122*c4762a1bSJed Brown Functional func; 123*c4762a1bSJed Brown char *names[256]; 124*c4762a1bSJed Brown PetscInt f; 125*c4762a1bSJed Brown PetscErrorCode ierr; 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown PetscFunctionBeginUser; 128*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Simulation Monitor Options", "DMPLEX");CHKERRQ(ierr); 129*c4762a1bSJed Brown options->numMonitorFuncs = ALEN(names); 130*c4762a1bSJed Brown ierr = PetscOptionsStringArray("-monitor", "List of functionals to monitor", "", names, &options->numMonitorFuncs, NULL);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = PetscMalloc1(options->numMonitorFuncs, &options->monitorFuncs);CHKERRQ(ierr); 132*c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 133*c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 134*c4762a1bSJed Brown PetscBool match; 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown ierr = PetscStrcasecmp(names[f], func->name, &match);CHKERRQ(ierr); 137*c4762a1bSJed Brown if (match) break; 138*c4762a1bSJed Brown } 139*c4762a1bSJed Brown if (!func) SETERRQ1(comm, PETSC_ERR_USER, "No known functional '%s'", names[f]); 140*c4762a1bSJed Brown options->monitorFuncs[f] = func; 141*c4762a1bSJed Brown /* Jed inserts a de-duplication of functionals here */ 142*c4762a1bSJed Brown ierr = PetscFree(names[f]);CHKERRQ(ierr); 143*c4762a1bSJed Brown } 144*c4762a1bSJed Brown /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ 145*c4762a1bSJed Brown options->maxMonitorFunc = -1; 146*c4762a1bSJed Brown for (func = options->functionalRegistry; func; func = func->next) { 147*c4762a1bSJed Brown for (f = 0; f < options->numMonitorFuncs; ++f) { 148*c4762a1bSJed Brown Functional call = options->monitorFuncs[f]; 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown if (func->func == call->func && func->ctx == call->ctx) options->maxMonitorFunc = PetscMax(options->maxMonitorFunc, func->offset); 151*c4762a1bSJed Brown } 152*c4762a1bSJed Brown } 153*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 154*c4762a1bSJed Brown PetscFunctionReturn(0); 155*c4762a1bSJed Brown } 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown static PetscErrorCode FunctionalRegister(Functional *functionalRegistry, const char name[], PetscInt *offset, FunctionalFunc func, void *ctx) 158*c4762a1bSJed Brown { 159*c4762a1bSJed Brown Functional *ptr, f; 160*c4762a1bSJed Brown PetscInt lastoffset = -1; 161*c4762a1bSJed Brown PetscErrorCode ierr; 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown PetscFunctionBeginUser; 164*c4762a1bSJed Brown for (ptr = functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; 165*c4762a1bSJed Brown ierr = PetscNew(&f);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = PetscStrallocpy(name, &f->name);CHKERRQ(ierr); 167*c4762a1bSJed Brown f->offset = lastoffset + 1; 168*c4762a1bSJed Brown f->func = func; 169*c4762a1bSJed Brown f->ctx = ctx; 170*c4762a1bSJed Brown f->next = NULL; 171*c4762a1bSJed Brown *ptr = f; 172*c4762a1bSJed Brown *offset = f->offset; 173*c4762a1bSJed Brown PetscFunctionReturn(0); 174*c4762a1bSJed Brown } 175*c4762a1bSJed Brown 176*c4762a1bSJed Brown static PetscErrorCode FunctionalDestroy(Functional *link) 177*c4762a1bSJed Brown { 178*c4762a1bSJed Brown Functional next, l; 179*c4762a1bSJed Brown PetscErrorCode ierr; 180*c4762a1bSJed Brown 181*c4762a1bSJed Brown PetscFunctionBeginUser; 182*c4762a1bSJed Brown if (!link) PetscFunctionReturn(0); 183*c4762a1bSJed Brown l = *link; 184*c4762a1bSJed Brown *link = NULL; 185*c4762a1bSJed Brown for (; l; l=next) { 186*c4762a1bSJed Brown next = l->next; 187*c4762a1bSJed Brown ierr = PetscFree(l->name);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = PetscFree(l);CHKERRQ(ierr); 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown PetscFunctionReturn(0); 191*c4762a1bSJed Brown } 192*c4762a1bSJed Brown 193*c4762a1bSJed Brown static void f0_zero_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 194*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 195*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 196*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 197*c4762a1bSJed Brown { 198*c4762a1bSJed Brown PetscInt comp; 199*c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = u[comp]; 200*c4762a1bSJed Brown } 201*c4762a1bSJed Brown 202*c4762a1bSJed Brown static void f0_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 203*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 204*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 205*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 206*c4762a1bSJed Brown { 207*c4762a1bSJed Brown PetscScalar wind[3] = {0.0, 0.0, 0.0}; 208*c4762a1bSJed Brown PetscInt comp; 209*c4762a1bSJed Brown 210*c4762a1bSJed Brown constant_u_2d(dim, t, x, Nf, wind, NULL); 211*c4762a1bSJed Brown for (comp = 0; comp < dim && comp < 3; ++comp) f0[comp] = u[comp] - wind[comp]; 212*c4762a1bSJed Brown } 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown static void f1_constant_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 215*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 216*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 217*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 218*c4762a1bSJed Brown { 219*c4762a1bSJed Brown PetscInt comp; 220*c4762a1bSJed Brown for (comp = 0; comp < dim*dim; ++comp) f1[comp] = 0.0; 221*c4762a1bSJed Brown } 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown static void g0_constant_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 224*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 225*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 226*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 227*c4762a1bSJed Brown { 228*c4762a1bSJed Brown PetscInt d; 229*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d*dim+d] = 1.0; 230*c4762a1bSJed Brown } 231*c4762a1bSJed Brown 232*c4762a1bSJed Brown static void g0_constant_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 233*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 234*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 235*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 236*c4762a1bSJed Brown { 237*c4762a1bSJed Brown g0[0] = 1.0; 238*c4762a1bSJed Brown } 239*c4762a1bSJed Brown 240*c4762a1bSJed Brown static void f0_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 241*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 242*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 243*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 244*c4762a1bSJed Brown { 245*c4762a1bSJed Brown PetscInt comp; 246*c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) f0[comp] = 4.0; 247*c4762a1bSJed Brown } 248*c4762a1bSJed Brown 249*c4762a1bSJed Brown static void f1_lap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 250*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 251*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 252*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 253*c4762a1bSJed Brown { 254*c4762a1bSJed Brown PetscInt comp, d; 255*c4762a1bSJed Brown for (comp = 0; comp < dim; ++comp) { 256*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 257*c4762a1bSJed Brown f1[comp*dim+d] = u_x[comp*dim+d]; 258*c4762a1bSJed Brown } 259*c4762a1bSJed Brown } 260*c4762a1bSJed Brown } 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown static void f0_lap_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 263*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 264*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 265*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 266*c4762a1bSJed Brown { 267*c4762a1bSJed Brown f0[0] = -PetscSinReal(2.0*PETSC_PI*x[0]); 268*c4762a1bSJed Brown f0[1] = 2.0*PETSC_PI*x[1]*PetscCosReal(2.0*PETSC_PI*x[0]); 269*c4762a1bSJed Brown } 270*c4762a1bSJed Brown 271*c4762a1bSJed Brown static void f0_lap_doubly_periodic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 272*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 273*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 274*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 275*c4762a1bSJed Brown { 276*c4762a1bSJed Brown f0[0] = -2.0*PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1]); 277*c4762a1bSJed Brown f0[1] = 2.0*PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0]); 278*c4762a1bSJed Brown } 279*c4762a1bSJed Brown 280*c4762a1bSJed Brown void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 281*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 282*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 283*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 284*c4762a1bSJed Brown { 285*c4762a1bSJed Brown const PetscInt Ncomp = dim; 286*c4762a1bSJed Brown PetscInt compI, d; 287*c4762a1bSJed Brown 288*c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 289*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 290*c4762a1bSJed Brown g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0; 291*c4762a1bSJed Brown } 292*c4762a1bSJed Brown } 293*c4762a1bSJed Brown } 294*c4762a1bSJed Brown 295*c4762a1bSJed Brown /* \frac{\partial\phi}{\partial t} + \nabla\phi \cdot \mathbf{u} + \phi \nabla \cdot \mathbf{u} = 0 */ 296*c4762a1bSJed Brown static void f0_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 297*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 298*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 299*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 300*c4762a1bSJed Brown { 301*c4762a1bSJed Brown PetscInt d; 302*c4762a1bSJed Brown f0[0] = u_t[dim]; 303*c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += u[dim]*u_x[d*dim+d] + u_x[dim*dim+d]*u[d]; 304*c4762a1bSJed Brown } 305*c4762a1bSJed Brown 306*c4762a1bSJed Brown static void f1_advection(PetscInt dim, PetscInt Nf, PetscInt NfAux, 307*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 308*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 309*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 310*c4762a1bSJed Brown { 311*c4762a1bSJed Brown PetscInt d; 312*c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[0] = 0.0; 313*c4762a1bSJed Brown } 314*c4762a1bSJed Brown 315*c4762a1bSJed Brown void g0_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 316*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 317*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 318*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 319*c4762a1bSJed Brown { 320*c4762a1bSJed Brown PetscInt d; 321*c4762a1bSJed Brown g0[0] = u_tShift; 322*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[d*dim+d]; 323*c4762a1bSJed Brown } 324*c4762a1bSJed Brown 325*c4762a1bSJed Brown void g1_adv_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 326*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 327*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 328*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 329*c4762a1bSJed Brown { 330*c4762a1bSJed Brown PetscInt d; 331*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = u[d]; 332*c4762a1bSJed Brown } 333*c4762a1bSJed Brown 334*c4762a1bSJed Brown void g0_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 335*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 336*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 337*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 338*c4762a1bSJed Brown { 339*c4762a1bSJed Brown PetscInt d; 340*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[0] += u_x[dim*dim+d]; 341*c4762a1bSJed Brown } 342*c4762a1bSJed Brown 343*c4762a1bSJed Brown void g1_adv_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 344*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 345*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 346*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 347*c4762a1bSJed Brown { 348*c4762a1bSJed Brown PetscInt d; 349*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = u[dim]; 350*c4762a1bSJed Brown } 351*c4762a1bSJed Brown 352*c4762a1bSJed 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) 353*c4762a1bSJed Brown { 354*c4762a1bSJed Brown PetscReal wind[3] = {0.0, 1.0, 0.0}; 355*c4762a1bSJed Brown PetscReal wn = DMPlex_DotRealD_Internal(PetscMin(dim,3), wind, n); 356*c4762a1bSJed Brown 357*c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 358*c4762a1bSJed Brown } 359*c4762a1bSJed Brown 360*c4762a1bSJed 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) 361*c4762a1bSJed Brown { 362*c4762a1bSJed Brown PetscReal wn = DMPlex_DotD_Internal(dim, uL, n); 363*c4762a1bSJed Brown 364*c4762a1bSJed Brown #if 1 365*c4762a1bSJed Brown flux[0] = (wn > 0 ? uL[dim] : uR[dim]) * wn; 366*c4762a1bSJed Brown #else 367*c4762a1bSJed 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]); */ 368*c4762a1bSJed Brown /* Smear it out */ 369*c4762a1bSJed Brown flux[0] = 0.5*((uL[dim] + uR[dim]) + (uL[dim] - uR[dim])*tanh(1.0e5*wn)) * wn; 370*c4762a1bSJed Brown #endif 371*c4762a1bSJed Brown } 372*c4762a1bSJed Brown 373*c4762a1bSJed Brown static PetscErrorCode zero_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 374*c4762a1bSJed Brown { 375*c4762a1bSJed Brown u[0] = 0.0; 376*c4762a1bSJed Brown u[1] = 0.0; 377*c4762a1bSJed Brown return 0; 378*c4762a1bSJed Brown } 379*c4762a1bSJed Brown 380*c4762a1bSJed Brown static PetscErrorCode constant_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 381*c4762a1bSJed Brown { 382*c4762a1bSJed Brown u[0] = 0.0; 383*c4762a1bSJed Brown u[1] = 1.0; 384*c4762a1bSJed Brown return 0; 385*c4762a1bSJed Brown } 386*c4762a1bSJed Brown 387*c4762a1bSJed Brown /* Coordinates of the point which was at x at t = 0 */ 388*c4762a1bSJed Brown static PetscErrorCode constant_x_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 389*c4762a1bSJed Brown { 390*c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 391*c4762a1bSJed Brown PetscErrorCode ierr; 392*c4762a1bSJed Brown u[0] = x[0]; 393*c4762a1bSJed Brown u[1] = x[1] + t; 394*c4762a1bSJed Brown ierr = DMLocalizeCoordinate(globalUser->dm, u, PETSC_FALSE, u);CHKERRQ(ierr); 395*c4762a1bSJed Brown return 0; 396*c4762a1bSJed Brown } 397*c4762a1bSJed Brown 398*c4762a1bSJed Brown /* 399*c4762a1bSJed Brown In 2D we use the exact solution: 400*c4762a1bSJed Brown 401*c4762a1bSJed Brown u = x^2 + y^2 402*c4762a1bSJed Brown v = 2 x^2 - 2xy 403*c4762a1bSJed Brown phi = h(x + y + (u + v) t) 404*c4762a1bSJed Brown f_x = f_y = 4 405*c4762a1bSJed Brown 406*c4762a1bSJed Brown so that 407*c4762a1bSJed Brown 408*c4762a1bSJed Brown -\Delta u + f = <-4, -4> + <4, 4> = 0 409*c4762a1bSJed Brown {\partial\phi}{\partial t} - \nabla\cdot \phi u = 0 410*c4762a1bSJed Brown h_t(x + y + (u + v) t) - u . grad phi - phi div u 411*c4762a1bSJed Brown = u h' + v h' - u h_x - v h_y 412*c4762a1bSJed Brown = 0 413*c4762a1bSJed Brown 414*c4762a1bSJed Brown We will conserve phi since 415*c4762a1bSJed Brown 416*c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 417*c4762a1bSJed Brown 418*c4762a1bSJed Brown Also try h((x + ut)^2 + (y + vt)^2), so that 419*c4762a1bSJed Brown 420*c4762a1bSJed Brown h_t((x + ut)^2 + (y + vt)^2) - u . grad phi - phi div u 421*c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h_x - v h_y 422*c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt)) - u h' 2 (x + u t) - v h' 2 (y + vt) 423*c4762a1bSJed Brown = 2 h' (u (x + ut) + v (y + vt) - u (x + u t) - v (y + vt)) 424*c4762a1bSJed Brown = 0 425*c4762a1bSJed Brown 426*c4762a1bSJed Brown */ 427*c4762a1bSJed Brown static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 428*c4762a1bSJed Brown { 429*c4762a1bSJed Brown u[0] = x[0]*x[0] + x[1]*x[1]; 430*c4762a1bSJed Brown u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 431*c4762a1bSJed Brown return 0; 432*c4762a1bSJed Brown } 433*c4762a1bSJed Brown 434*c4762a1bSJed Brown /* 435*c4762a1bSJed Brown In 2D we use the exact, periodic solution: 436*c4762a1bSJed Brown 437*c4762a1bSJed Brown u = sin(2 pi x)/4 pi^2 438*c4762a1bSJed Brown v = -y cos(2 pi x)/2 pi 439*c4762a1bSJed Brown phi = h(x + y + (u + v) t) 440*c4762a1bSJed Brown f_x = -sin(2 pi x) 441*c4762a1bSJed Brown f_y = 2 pi y cos(2 pi x) 442*c4762a1bSJed Brown 443*c4762a1bSJed Brown so that 444*c4762a1bSJed Brown 445*c4762a1bSJed Brown -\Delta u + f = <sin(2pi x), -2pi y cos(2pi x)> + <-sin(2pi x), 2pi y cos(2pi x)> = 0 446*c4762a1bSJed Brown 447*c4762a1bSJed Brown We will conserve phi since 448*c4762a1bSJed Brown 449*c4762a1bSJed Brown \nabla \cdot u = cos(2pi x)/2pi - cos(2pi x)/2pi = 0 450*c4762a1bSJed Brown */ 451*c4762a1bSJed Brown static PetscErrorCode periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 452*c4762a1bSJed Brown { 453*c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 454*c4762a1bSJed Brown u[1] = -x[1]*PetscCosReal(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI); 455*c4762a1bSJed Brown return 0; 456*c4762a1bSJed Brown } 457*c4762a1bSJed Brown 458*c4762a1bSJed Brown /* 459*c4762a1bSJed Brown In 2D we use the exact, doubly periodic solution: 460*c4762a1bSJed Brown 461*c4762a1bSJed Brown u = sin(2 pi x) cos(2 pi y)/4 pi^2 462*c4762a1bSJed Brown v = -sin(2 pi y) cos(2 pi x)/4 pi^2 463*c4762a1bSJed Brown phi = h(x + y + (u + v) t) 464*c4762a1bSJed Brown f_x = -2sin(2 pi x) cos(2 pi y) 465*c4762a1bSJed Brown f_y = 2sin(2 pi y) cos(2 pi x) 466*c4762a1bSJed Brown 467*c4762a1bSJed Brown so that 468*c4762a1bSJed Brown 469*c4762a1bSJed 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 470*c4762a1bSJed Brown 471*c4762a1bSJed Brown We will conserve phi since 472*c4762a1bSJed Brown 473*c4762a1bSJed Brown \nabla \cdot u = cos(2pi x) cos(2pi y)/2pi - cos(2pi y) cos(2pi x)/2pi = 0 474*c4762a1bSJed Brown */ 475*c4762a1bSJed Brown static PetscErrorCode doubly_periodic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 476*c4762a1bSJed Brown { 477*c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0])*PetscCosReal(2.0*PETSC_PI*x[1])/PetscSqr(2.0*PETSC_PI); 478*c4762a1bSJed Brown u[1] = -PetscSinReal(2.0*PETSC_PI*x[1])*PetscCosReal(2.0*PETSC_PI*x[0])/PetscSqr(2.0*PETSC_PI); 479*c4762a1bSJed Brown return 0; 480*c4762a1bSJed Brown } 481*c4762a1bSJed Brown 482*c4762a1bSJed Brown static PetscErrorCode shear_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 483*c4762a1bSJed Brown { 484*c4762a1bSJed Brown u[0] = x[1] - 0.5; 485*c4762a1bSJed Brown u[1] = 0.0; 486*c4762a1bSJed Brown return 0; 487*c4762a1bSJed Brown } 488*c4762a1bSJed Brown 489*c4762a1bSJed Brown static PetscErrorCode initialVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 490*c4762a1bSJed Brown { 491*c4762a1bSJed Brown PetscInt d; 492*c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 493*c4762a1bSJed Brown return 0; 494*c4762a1bSJed Brown } 495*c4762a1bSJed Brown 496*c4762a1bSJed Brown static PetscErrorCode zero_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 497*c4762a1bSJed Brown { 498*c4762a1bSJed Brown u[0] = 0.0; 499*c4762a1bSJed Brown return 0; 500*c4762a1bSJed Brown } 501*c4762a1bSJed Brown 502*c4762a1bSJed Brown static PetscErrorCode constant_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 503*c4762a1bSJed Brown { 504*c4762a1bSJed Brown u[0] = 1.0; 505*c4762a1bSJed Brown return 0; 506*c4762a1bSJed Brown } 507*c4762a1bSJed Brown 508*c4762a1bSJed Brown static PetscErrorCode delta_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 509*c4762a1bSJed Brown { 510*c4762a1bSJed Brown PetscReal x0[2]; 511*c4762a1bSJed Brown PetscScalar xn[2]; 512*c4762a1bSJed Brown 513*c4762a1bSJed Brown x0[0] = globalUser->source[0]; 514*c4762a1bSJed Brown x0[1] = globalUser->source[1]; 515*c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 516*c4762a1bSJed Brown { 517*c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 518*c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 519*c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 520*c4762a1bSJed Brown 521*c4762a1bSJed Brown u[0] = r2 < 1.0e-7 ? 1.0 : 0.0; 522*c4762a1bSJed Brown } 523*c4762a1bSJed Brown return 0; 524*c4762a1bSJed Brown } 525*c4762a1bSJed Brown 526*c4762a1bSJed Brown /* 527*c4762a1bSJed Brown Gaussian blob, initially centered on (0.5, 0.5) 528*c4762a1bSJed Brown 529*c4762a1bSJed Brown xi = x(t) - x0, eta = y(t) - y0 530*c4762a1bSJed Brown 531*c4762a1bSJed Brown where x(t), y(t) are the integral curves of v(t), 532*c4762a1bSJed Brown 533*c4762a1bSJed Brown dx/dt . grad f = v . f 534*c4762a1bSJed Brown 535*c4762a1bSJed Brown Check: constant v(t) = {v0, w0}, x(t) = {x0 + v0 t, y0 + w0 t} 536*c4762a1bSJed Brown 537*c4762a1bSJed Brown v0 f_x + w0 f_y = v . f 538*c4762a1bSJed Brown */ 539*c4762a1bSJed Brown static PetscErrorCode gaussian_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 540*c4762a1bSJed Brown { 541*c4762a1bSJed Brown const PetscReal x0[2] = {0.5, 0.5}; 542*c4762a1bSJed Brown const PetscReal sigma = 1.0/6.0; 543*c4762a1bSJed Brown PetscScalar xn[2]; 544*c4762a1bSJed Brown 545*c4762a1bSJed Brown constant_x_2d(dim, time, x0, Nf, xn, ctx); 546*c4762a1bSJed Brown { 547*c4762a1bSJed Brown /* const PetscReal xi = x[0] + (sin(2.0*PETSC_PI*x[0])/(4.0*PETSC_PI*PETSC_PI))*t - x0[0]; */ 548*c4762a1bSJed Brown /* const PetscReal eta = x[1] + (-x[1]*cos(2.0*PETSC_PI*x[0])/(2.0*PETSC_PI))*t - x0[1]; */ 549*c4762a1bSJed Brown const PetscReal xi = x[0] - PetscRealPart(xn[0]); 550*c4762a1bSJed Brown const PetscReal eta = x[1] - PetscRealPart(xn[1]); 551*c4762a1bSJed Brown const PetscReal r2 = xi*xi + eta*eta; 552*c4762a1bSJed Brown 553*c4762a1bSJed Brown u[0] = PetscExpReal(-r2/(2.0*sigma*sigma))/(sigma*PetscSqrtReal(2.0*PETSC_PI)); 554*c4762a1bSJed Brown } 555*c4762a1bSJed Brown return 0; 556*c4762a1bSJed Brown } 557*c4762a1bSJed Brown 558*c4762a1bSJed Brown static PetscErrorCode tilted_phi_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 559*c4762a1bSJed Brown { 560*c4762a1bSJed Brown PetscReal x0[3]; 561*c4762a1bSJed Brown const PetscReal wind[3] = {0.0, 1.0, 0.0}; 562*c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 563*c4762a1bSJed Brown 564*c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, wind, x, x0); 565*c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 566*c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 567*c4762a1bSJed Brown return 0; 568*c4762a1bSJed Brown } 569*c4762a1bSJed Brown 570*c4762a1bSJed Brown static PetscErrorCode tilted_phi_coupled_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 571*c4762a1bSJed Brown { 572*c4762a1bSJed Brown PetscReal ur[3]; 573*c4762a1bSJed Brown PetscReal x0[3]; 574*c4762a1bSJed Brown const PetscReal t = *((PetscReal *) ctx); 575*c4762a1bSJed Brown 576*c4762a1bSJed Brown ur[0] = PetscRealPart(u[0]); ur[1] = PetscRealPart(u[1]); ur[2] = PetscRealPart(u[2]); 577*c4762a1bSJed Brown DMPlex_WaxpyD_Internal(2, -t, ur, x, x0); 578*c4762a1bSJed Brown if (x0[1] > 0) u[0] = 1.0*x[0] + 3.0*x[1]; 579*c4762a1bSJed Brown else u[0] = -2.0; /* Inflow state */ 580*c4762a1bSJed Brown return 0; 581*c4762a1bSJed Brown } 582*c4762a1bSJed Brown 583*c4762a1bSJed Brown static PetscErrorCode advect_inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 584*c4762a1bSJed Brown { 585*c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 586*c4762a1bSJed Brown 587*c4762a1bSJed Brown PetscFunctionBeginUser; 588*c4762a1bSJed Brown xG[0] = user->inflowState; 589*c4762a1bSJed Brown PetscFunctionReturn(0); 590*c4762a1bSJed Brown } 591*c4762a1bSJed Brown 592*c4762a1bSJed Brown static PetscErrorCode advect_outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx) 593*c4762a1bSJed Brown { 594*c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 595*c4762a1bSJed Brown 596*c4762a1bSJed Brown PetscFunctionBeginUser; 597*c4762a1bSJed Brown xG[0] = xI[user->dim]; 598*c4762a1bSJed Brown PetscFunctionReturn(0); 599*c4762a1bSJed Brown } 600*c4762a1bSJed Brown 601*c4762a1bSJed Brown static PetscErrorCode ExactSolution(DM dm, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx) 602*c4762a1bSJed Brown { 603*c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 604*c4762a1bSJed Brown PetscInt dim; 605*c4762a1bSJed Brown PetscErrorCode ierr; 606*c4762a1bSJed Brown 607*c4762a1bSJed Brown PetscFunctionBeginUser; 608*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 609*c4762a1bSJed Brown switch (user->porosityDist) { 610*c4762a1bSJed Brown case TILTED: 611*c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) tilted_phi_2d(dim, time, x, 2, u, (void *) &time); 612*c4762a1bSJed Brown else tilted_phi_coupled_2d(dim, time, x, 2, u, (void *) &time); 613*c4762a1bSJed Brown break; 614*c4762a1bSJed Brown case GAUSSIAN: 615*c4762a1bSJed Brown gaussian_phi_2d(dim, time, x, 2, u, (void *) &time); 616*c4762a1bSJed Brown break; 617*c4762a1bSJed Brown case DELTA: 618*c4762a1bSJed Brown delta_phi_2d(dim, time, x, 2, u, (void *) &time); 619*c4762a1bSJed Brown break; 620*c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); 621*c4762a1bSJed Brown } 622*c4762a1bSJed Brown PetscFunctionReturn(0); 623*c4762a1bSJed Brown } 624*c4762a1bSJed Brown 625*c4762a1bSJed Brown static PetscErrorCode Functional_Error(DM dm, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, void *ctx) 626*c4762a1bSJed Brown { 627*c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 628*c4762a1bSJed Brown PetscScalar yexact[3]={0,0,0}; 629*c4762a1bSJed Brown PetscErrorCode ierr; 630*c4762a1bSJed Brown 631*c4762a1bSJed Brown PetscFunctionBeginUser; 632*c4762a1bSJed Brown ierr = ExactSolution(dm, time, x, yexact, ctx);CHKERRQ(ierr); 633*c4762a1bSJed Brown f[user->errorFunctional] = PetscAbsScalar(y[0] - yexact[0]); 634*c4762a1bSJed Brown PetscFunctionReturn(0); 635*c4762a1bSJed Brown } 636*c4762a1bSJed Brown 637*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 638*c4762a1bSJed Brown { 639*c4762a1bSJed Brown DM distributedMesh = NULL; 640*c4762a1bSJed Brown PetscBool periodic = user->bd[0] == DM_BOUNDARY_PERIODIC || user->bd[0] == DM_BOUNDARY_TWIST || user->bd[1] == DM_BOUNDARY_PERIODIC || user->bd[1] == DM_BOUNDARY_TWIST ? PETSC_TRUE : PETSC_FALSE; 641*c4762a1bSJed Brown const char *filename = user->filename; 642*c4762a1bSJed Brown const PetscInt cells[3] = {3, 3, 3}; 643*c4762a1bSJed Brown const PetscReal L[3] = {1.0, 1.0, 1.0}; 644*c4762a1bSJed Brown PetscReal maxCell[3]; 645*c4762a1bSJed Brown PetscInt d; 646*c4762a1bSJed Brown size_t len; 647*c4762a1bSJed Brown PetscErrorCode ierr; 648*c4762a1bSJed Brown 649*c4762a1bSJed Brown PetscFunctionBeginUser; 650*c4762a1bSJed Brown ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 651*c4762a1bSJed Brown if (!len) { 652*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, user->dim, PETSC_FALSE, cells, NULL, NULL, user->bd, PETSC_TRUE, dm);CHKERRQ(ierr); 653*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 654*c4762a1bSJed Brown } else { 655*c4762a1bSJed Brown ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr); 656*c4762a1bSJed Brown } 657*c4762a1bSJed Brown if (periodic) {for (d = 0; d < 3; ++d) maxCell[d] = 1.1*(L[d]/cells[d]); ierr = DMSetPeriodicity(*dm, PETSC_TRUE, maxCell, L, user->bd);CHKERRQ(ierr);} 658*c4762a1bSJed Brown /* Distribute mesh */ 659*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 660*c4762a1bSJed Brown if (distributedMesh) { 661*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 662*c4762a1bSJed Brown *dm = distributedMesh; 663*c4762a1bSJed Brown } 664*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 665*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 666*c4762a1bSJed Brown PetscFunctionReturn(0); 667*c4762a1bSJed Brown } 668*c4762a1bSJed Brown 669*c4762a1bSJed Brown static PetscErrorCode SetupBC(DM dm, AppCtx *user) 670*c4762a1bSJed Brown { 671*c4762a1bSJed Brown PetscDS prob; 672*c4762a1bSJed Brown DMLabel label; 673*c4762a1bSJed Brown PetscBool check; 674*c4762a1bSJed Brown PetscErrorCode ierr; 675*c4762a1bSJed Brown 676*c4762a1bSJed Brown PetscFunctionBeginUser; 677*c4762a1bSJed Brown /* Set initial guesses and exact solutions */ 678*c4762a1bSJed Brown switch (user->dim) { 679*c4762a1bSJed Brown case 2: 680*c4762a1bSJed Brown user->initialGuess[0] = initialVelocity; 681*c4762a1bSJed Brown switch(user->porosityDist) { 682*c4762a1bSJed Brown case ZERO: user->initialGuess[1] = zero_phi;break; 683*c4762a1bSJed Brown case CONSTANT: user->initialGuess[1] = constant_phi;break; 684*c4762a1bSJed Brown case GAUSSIAN: user->initialGuess[1] = gaussian_phi_2d;break; 685*c4762a1bSJed Brown case DELTA: user->initialGuess[1] = delta_phi_2d;break; 686*c4762a1bSJed Brown case TILTED: 687*c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) user->initialGuess[1] = tilted_phi_2d; 688*c4762a1bSJed Brown else user->initialGuess[1] = tilted_phi_coupled_2d; 689*c4762a1bSJed Brown break; 690*c4762a1bSJed Brown } 691*c4762a1bSJed Brown user->exactFuncs[1] = user->initialGuess[1]; 692*c4762a1bSJed Brown switch (user->velocityDist) { 693*c4762a1bSJed Brown case VEL_ZERO: 694*c4762a1bSJed Brown user->exactFuncs[0] = zero_u_2d; break; 695*c4762a1bSJed Brown case VEL_CONSTANT: 696*c4762a1bSJed Brown user->exactFuncs[0] = constant_u_2d; break; 697*c4762a1bSJed Brown case VEL_HARMONIC: 698*c4762a1bSJed Brown switch (user->bd[0]) { 699*c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 700*c4762a1bSJed Brown switch (user->bd[1]) { 701*c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 702*c4762a1bSJed Brown user->exactFuncs[0] = doubly_periodic_u_2d; break; 703*c4762a1bSJed Brown default: 704*c4762a1bSJed Brown user->exactFuncs[0] = periodic_u_2d; break; 705*c4762a1bSJed Brown } 706*c4762a1bSJed Brown break; 707*c4762a1bSJed Brown default: 708*c4762a1bSJed Brown user->exactFuncs[0] = quadratic_u_2d; break; 709*c4762a1bSJed Brown } 710*c4762a1bSJed Brown break; 711*c4762a1bSJed Brown case VEL_SHEAR: 712*c4762a1bSJed Brown user->exactFuncs[0] = shear_bc; break; 713*c4762a1bSJed Brown break; 714*c4762a1bSJed Brown default: 715*c4762a1bSJed Brown SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", user->dim); 716*c4762a1bSJed Brown } 717*c4762a1bSJed Brown } 718*c4762a1bSJed Brown { 719*c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 720*c4762a1bSJed Brown 721*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit);CHKERRQ(ierr); 722*c4762a1bSJed Brown if (user->velocityDist == VEL_CONSTANT && !isImplicit) user->initialGuess[0] = user->exactFuncs[0]; 723*c4762a1bSJed Brown } 724*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-dmts_check", &check);CHKERRQ(ierr); 725*c4762a1bSJed Brown if (check) { 726*c4762a1bSJed Brown user->initialGuess[0] = user->exactFuncs[0]; 727*c4762a1bSJed Brown user->initialGuess[1] = user->exactFuncs[1]; 728*c4762a1bSJed Brown } 729*c4762a1bSJed Brown /* Set BC */ 730*c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 731*c4762a1bSJed Brown ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 732*c4762a1bSJed Brown if (label) { 733*c4762a1bSJed Brown const PetscInt id = 1; 734*c4762a1bSJed Brown 735*c4762a1bSJed Brown ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) user->exactFuncs[0], 1, &id, user);CHKERRQ(ierr); 736*c4762a1bSJed Brown } 737*c4762a1bSJed Brown ierr = DMGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr); 738*c4762a1bSJed Brown if (label && user->useFV) { 739*c4762a1bSJed Brown const PetscInt inflowids[] = {100,200,300}, outflowids[] = {101}; 740*c4762a1bSJed Brown 741*c4762a1bSJed Brown ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow", "Face Sets", 1, 0, NULL, (void (*)(void)) advect_inflow, ALEN(inflowids), inflowids, user);CHKERRQ(ierr); 742*c4762a1bSJed Brown ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", "Face Sets", 1, 0, NULL, (void (*)(void)) advect_outflow, ALEN(outflowids), outflowids, user);CHKERRQ(ierr); 743*c4762a1bSJed Brown } 744*c4762a1bSJed Brown PetscFunctionReturn(0); 745*c4762a1bSJed Brown } 746*c4762a1bSJed Brown 747*c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 748*c4762a1bSJed Brown { 749*c4762a1bSJed Brown PetscDS prob; 750*c4762a1bSJed Brown PetscErrorCode ierr; 751*c4762a1bSJed Brown 752*c4762a1bSJed Brown PetscFunctionBeginUser; 753*c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 754*c4762a1bSJed Brown switch (user->velocityDist) { 755*c4762a1bSJed Brown case VEL_ZERO: 756*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_zero_u, f1_constant_u);CHKERRQ(ierr); 757*c4762a1bSJed Brown break; 758*c4762a1bSJed Brown case VEL_CONSTANT: 759*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_constant_u, f1_constant_u);CHKERRQ(ierr); 760*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, g0_constant_uu, NULL, NULL, NULL);CHKERRQ(ierr); 761*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 1, 1, g0_constant_pp, NULL, NULL, NULL);CHKERRQ(ierr); 762*c4762a1bSJed Brown break; 763*c4762a1bSJed Brown case VEL_HARMONIC: 764*c4762a1bSJed Brown switch (user->bd[0]) { 765*c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 766*c4762a1bSJed Brown switch (user->bd[1]) { 767*c4762a1bSJed Brown case DM_BOUNDARY_PERIODIC: 768*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_lap_doubly_periodic_u, f1_lap_u);CHKERRQ(ierr); 769*c4762a1bSJed Brown break; 770*c4762a1bSJed Brown default: 771*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_lap_periodic_u, f1_lap_u);CHKERRQ(ierr); 772*c4762a1bSJed Brown break; 773*c4762a1bSJed Brown } 774*c4762a1bSJed Brown break; 775*c4762a1bSJed Brown default: 776*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_lap_u, f1_lap_u);CHKERRQ(ierr); 777*c4762a1bSJed Brown break; 778*c4762a1bSJed Brown } 779*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 780*c4762a1bSJed Brown break; 781*c4762a1bSJed Brown case VEL_SHEAR: 782*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_zero_u, f1_lap_u);CHKERRQ(ierr); 783*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 784*c4762a1bSJed Brown break; 785*c4762a1bSJed Brown } 786*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 1, f0_advection, f1_advection);CHKERRQ(ierr); 787*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 1, 1, g0_adv_pp, g1_adv_pp, NULL, NULL);CHKERRQ(ierr); 788*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 1, 0, g0_adv_pu, g1_adv_pu, NULL, NULL);CHKERRQ(ierr); 789*c4762a1bSJed Brown if (user->velocityDist == VEL_ZERO) {ierr = PetscDSSetRiemannSolver(prob, 1, riemann_advection);CHKERRQ(ierr);} 790*c4762a1bSJed Brown else {ierr = PetscDSSetRiemannSolver(prob, 1, riemann_coupled_advection);CHKERRQ(ierr);} 791*c4762a1bSJed Brown 792*c4762a1bSJed Brown ierr = FunctionalRegister(&user->functionalRegistry, "Error", &user->errorFunctional, Functional_Error, user);CHKERRQ(ierr); 793*c4762a1bSJed Brown PetscFunctionReturn(0); 794*c4762a1bSJed Brown } 795*c4762a1bSJed Brown 796*c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 797*c4762a1bSJed Brown { 798*c4762a1bSJed Brown DM cdm = dm; 799*c4762a1bSJed Brown const PetscInt dim = user->dim; 800*c4762a1bSJed Brown PetscQuadrature q; 801*c4762a1bSJed Brown PetscFE fe[2]; 802*c4762a1bSJed Brown PetscFV fv; 803*c4762a1bSJed Brown MPI_Comm comm; 804*c4762a1bSJed Brown PetscErrorCode ierr; 805*c4762a1bSJed Brown 806*c4762a1bSJed Brown PetscFunctionBeginUser; 807*c4762a1bSJed Brown /* Create finite element */ 808*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 809*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, PETSC_FALSE, "velocity_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 810*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 811*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_FALSE, "porosity_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 812*c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 813*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[1], "porosity");CHKERRQ(ierr); 814*c4762a1bSJed Brown 815*c4762a1bSJed Brown ierr = PetscFVCreate(PetscObjectComm((PetscObject) dm), &fv);CHKERRQ(ierr); 816*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fv, "porosity");CHKERRQ(ierr); 817*c4762a1bSJed Brown ierr = PetscFVSetFromOptions(fv);CHKERRQ(ierr); 818*c4762a1bSJed Brown ierr = PetscFVSetNumComponents(fv, 1);CHKERRQ(ierr); 819*c4762a1bSJed Brown ierr = PetscFVSetSpatialDimension(fv, dim);CHKERRQ(ierr); 820*c4762a1bSJed Brown ierr = PetscFEGetQuadrature(fe[0], &q);CHKERRQ(ierr); 821*c4762a1bSJed Brown ierr = PetscFVSetQuadrature(fv, q);CHKERRQ(ierr); 822*c4762a1bSJed Brown 823*c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 824*c4762a1bSJed Brown if (user->useFV) {ierr = DMSetField(dm, 1, NULL, (PetscObject) fv);CHKERRQ(ierr);} 825*c4762a1bSJed Brown else {ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);} 826*c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 827*c4762a1bSJed Brown ierr = SetupProblem(dm, user);CHKERRQ(ierr); 828*c4762a1bSJed Brown 829*c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 830*c4762a1bSJed Brown while (cdm) { 831*c4762a1bSJed Brown ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 832*c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 833*c4762a1bSJed Brown /* Coordinates were never localized for coarse meshes */ 834*c4762a1bSJed Brown if (cdm) {ierr = DMLocalizeCoordinates(cdm);CHKERRQ(ierr);} 835*c4762a1bSJed Brown } 836*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 837*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 838*c4762a1bSJed Brown ierr = PetscFVDestroy(&fv);CHKERRQ(ierr); 839*c4762a1bSJed Brown PetscFunctionReturn(0); 840*c4762a1bSJed Brown } 841*c4762a1bSJed Brown 842*c4762a1bSJed Brown static PetscErrorCode CreateDM(MPI_Comm comm, AppCtx *user, DM *dm) 843*c4762a1bSJed Brown { 844*c4762a1bSJed Brown PetscErrorCode ierr; 845*c4762a1bSJed Brown 846*c4762a1bSJed Brown PetscFunctionBeginUser; 847*c4762a1bSJed Brown ierr = CreateMesh(comm, user, dm);CHKERRQ(ierr); 848*c4762a1bSJed Brown /* Handle refinement, etc. */ 849*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 850*c4762a1bSJed Brown /* Construct ghost cells */ 851*c4762a1bSJed Brown if (user->useFV) { 852*c4762a1bSJed Brown DM gdm; 853*c4762a1bSJed Brown 854*c4762a1bSJed Brown ierr = DMPlexConstructGhostCells(*dm, NULL, NULL, &gdm);CHKERRQ(ierr); 855*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 856*c4762a1bSJed Brown *dm = gdm; 857*c4762a1bSJed Brown } 858*c4762a1bSJed Brown /* Localize coordinates */ 859*c4762a1bSJed Brown ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 860*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)(*dm),"Mesh");CHKERRQ(ierr); 861*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 862*c4762a1bSJed Brown /* Setup problem */ 863*c4762a1bSJed Brown ierr = SetupDiscretization(*dm, user);CHKERRQ(ierr); 864*c4762a1bSJed Brown /* Setup BC */ 865*c4762a1bSJed Brown ierr = SetupBC(*dm, user);CHKERRQ(ierr); 866*c4762a1bSJed Brown PetscFunctionReturn(0); 867*c4762a1bSJed Brown } 868*c4762a1bSJed Brown 869*c4762a1bSJed Brown static PetscErrorCode SetInitialConditionFVM(DM dm, Vec X, PetscInt field, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx) 870*c4762a1bSJed Brown { 871*c4762a1bSJed Brown PetscDS prob; 872*c4762a1bSJed Brown DM dmCell; 873*c4762a1bSJed Brown Vec cellgeom; 874*c4762a1bSJed Brown const PetscScalar *cgeom; 875*c4762a1bSJed Brown PetscScalar *x; 876*c4762a1bSJed Brown PetscInt dim, Nf, cStart, cEnd, c; 877*c4762a1bSJed Brown PetscErrorCode ierr; 878*c4762a1bSJed Brown 879*c4762a1bSJed Brown PetscFunctionBeginUser; 880*c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 881*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 882*c4762a1bSJed Brown ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 883*c4762a1bSJed Brown ierr = DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);CHKERRQ(ierr); 884*c4762a1bSJed Brown ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr); 885*c4762a1bSJed Brown ierr = DMPlexGetInteriorCellStratum(dm, &cStart, &cEnd);CHKERRQ(ierr); 886*c4762a1bSJed Brown ierr = VecGetArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 887*c4762a1bSJed Brown ierr = VecGetArray(X, &x);CHKERRQ(ierr); 888*c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 889*c4762a1bSJed Brown PetscFVCellGeom *cg; 890*c4762a1bSJed Brown PetscScalar *xc; 891*c4762a1bSJed Brown 892*c4762a1bSJed Brown ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 893*c4762a1bSJed Brown ierr = DMPlexPointGlobalFieldRef(dm, c, field, x, &xc);CHKERRQ(ierr); 894*c4762a1bSJed Brown if (xc) {ierr = (*func)(dim, 0.0, cg->centroid, Nf, xc, ctx);CHKERRQ(ierr);} 895*c4762a1bSJed Brown } 896*c4762a1bSJed Brown ierr = VecRestoreArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 897*c4762a1bSJed Brown ierr = VecRestoreArray(X, &x);CHKERRQ(ierr); 898*c4762a1bSJed Brown PetscFunctionReturn(0); 899*c4762a1bSJed Brown } 900*c4762a1bSJed Brown 901*c4762a1bSJed Brown static PetscErrorCode MonitorFunctionals(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx) 902*c4762a1bSJed Brown { 903*c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 904*c4762a1bSJed Brown char *ftable = NULL; 905*c4762a1bSJed Brown DM dm; 906*c4762a1bSJed Brown PetscSection s; 907*c4762a1bSJed Brown Vec cellgeom; 908*c4762a1bSJed Brown const PetscScalar *x; 909*c4762a1bSJed Brown PetscScalar *a; 910*c4762a1bSJed Brown PetscReal *xnorms; 911*c4762a1bSJed Brown PetscInt pStart, pEnd, p, Nf, f; 912*c4762a1bSJed Brown PetscErrorCode ierr; 913*c4762a1bSJed Brown 914*c4762a1bSJed Brown PetscFunctionBeginUser; 915*c4762a1bSJed Brown ierr = VecViewFromOptions(X, (PetscObject) ts, "-view_solution");CHKERRQ(ierr); 916*c4762a1bSJed Brown ierr = VecGetDM(X, &dm);CHKERRQ(ierr); 917*c4762a1bSJed Brown ierr = DMPlexTSGetGeometryFVM(dm, NULL, &cellgeom, NULL);CHKERRQ(ierr); 918*c4762a1bSJed Brown ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 919*c4762a1bSJed Brown ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr); 920*c4762a1bSJed Brown ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 921*c4762a1bSJed Brown ierr = PetscCalloc1(Nf*2, &xnorms);CHKERRQ(ierr); 922*c4762a1bSJed Brown ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 923*c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 924*c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 925*c4762a1bSJed Brown PetscInt dof, cdof, d; 926*c4762a1bSJed Brown 927*c4762a1bSJed Brown ierr = PetscSectionGetFieldDof(s, p, f, &dof);CHKERRQ(ierr); 928*c4762a1bSJed Brown ierr = PetscSectionGetFieldConstraintDof(s, p, f, &cdof);CHKERRQ(ierr); 929*c4762a1bSJed Brown ierr = DMPlexPointGlobalFieldRead(dm, p, f, x, &a);CHKERRQ(ierr); 930*c4762a1bSJed Brown /* TODO Use constrained indices here */ 931*c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+0] = PetscMax(xnorms[f*2+0], PetscAbsScalar(a[d])); 932*c4762a1bSJed Brown for (d = 0; d < dof-cdof; ++d) xnorms[f*2+1] += PetscAbsScalar(a[d]); 933*c4762a1bSJed Brown } 934*c4762a1bSJed Brown } 935*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 936*c4762a1bSJed Brown if (stepnum >= 0) { /* No summary for final time */ 937*c4762a1bSJed Brown DM dmCell, *fdm; 938*c4762a1bSJed Brown Vec *fv; 939*c4762a1bSJed Brown const PetscScalar *cgeom; 940*c4762a1bSJed Brown PetscScalar **fx; 941*c4762a1bSJed Brown PetscReal *fmin, *fmax, *fint, *ftmp, t; 942*c4762a1bSJed Brown PetscInt cStart, cEnd, c, fcount, f, num; 943*c4762a1bSJed Brown 944*c4762a1bSJed Brown size_t ftableused,ftablealloc; 945*c4762a1bSJed Brown 946*c4762a1bSJed Brown /* Functionals have indices after registering, this is an upper bound */ 947*c4762a1bSJed Brown fcount = user->numMonitorFuncs; 948*c4762a1bSJed Brown ierr = PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fint,fcount,&ftmp);CHKERRQ(ierr); 949*c4762a1bSJed Brown ierr = PetscMalloc3(fcount,&fdm,fcount,&fv,fcount,&fx);CHKERRQ(ierr); 950*c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 951*c4762a1bSJed Brown PetscSection fs; 952*c4762a1bSJed Brown const char *name = user->monitorFuncs[f]->name; 953*c4762a1bSJed Brown 954*c4762a1bSJed Brown fmin[f] = PETSC_MAX_REAL; 955*c4762a1bSJed Brown fmax[f] = PETSC_MIN_REAL; 956*c4762a1bSJed Brown fint[f] = 0; 957*c4762a1bSJed Brown /* Make monitor vecs */ 958*c4762a1bSJed Brown ierr = DMClone(dm, &fdm[f]);CHKERRQ(ierr); 959*c4762a1bSJed Brown ierr = DMGetOutputSequenceNumber(dm, &num, &t);CHKERRQ(ierr); 960*c4762a1bSJed Brown ierr = DMSetOutputSequenceNumber(fdm[f], num, t);CHKERRQ(ierr); 961*c4762a1bSJed Brown ierr = PetscSectionClone(s, &fs);CHKERRQ(ierr); 962*c4762a1bSJed Brown ierr = PetscSectionSetFieldName(fs, 0, NULL);CHKERRQ(ierr); 963*c4762a1bSJed Brown ierr = PetscSectionSetFieldName(fs, 1, name);CHKERRQ(ierr); 964*c4762a1bSJed Brown ierr = DMSetLocalSection(fdm[f], fs);CHKERRQ(ierr); 965*c4762a1bSJed Brown ierr = PetscSectionDestroy(&fs);CHKERRQ(ierr); 966*c4762a1bSJed Brown ierr = DMGetGlobalVector(fdm[f], &fv[f]);CHKERRQ(ierr); 967*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fv[f], name);CHKERRQ(ierr); 968*c4762a1bSJed Brown ierr = VecGetArray(fv[f], &fx[f]);CHKERRQ(ierr); 969*c4762a1bSJed Brown } 970*c4762a1bSJed Brown ierr = DMPlexGetInteriorCellStratum(dm, &cStart, &cEnd);CHKERRQ(ierr); 971*c4762a1bSJed Brown ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr); 972*c4762a1bSJed Brown ierr = VecGetArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 973*c4762a1bSJed Brown ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 974*c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 975*c4762a1bSJed Brown PetscFVCellGeom *cg; 976*c4762a1bSJed Brown PetscScalar *cx; 977*c4762a1bSJed Brown 978*c4762a1bSJed Brown ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 979*c4762a1bSJed Brown ierr = DMPlexPointGlobalFieldRead(dm, c, 1, x, &cx);CHKERRQ(ierr); 980*c4762a1bSJed Brown if (!cx) continue; /* not a global cell */ 981*c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 982*c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 983*c4762a1bSJed Brown PetscScalar *fxc; 984*c4762a1bSJed Brown 985*c4762a1bSJed Brown ierr = DMPlexPointGlobalFieldRef(dm, c, 1, fx[f], &fxc);CHKERRQ(ierr); 986*c4762a1bSJed Brown /* I need to make it easier to get interpolated values here */ 987*c4762a1bSJed Brown ierr = (*func->func)(dm, time, cg->centroid, cx, ftmp, func->ctx);CHKERRQ(ierr); 988*c4762a1bSJed Brown fxc[0] = ftmp[user->monitorFuncs[f]->offset]; 989*c4762a1bSJed Brown } 990*c4762a1bSJed Brown for (f = 0; f < fcount; ++f) { 991*c4762a1bSJed Brown fmin[f] = PetscMin(fmin[f], ftmp[f]); 992*c4762a1bSJed Brown fmax[f] = PetscMax(fmax[f], ftmp[f]); 993*c4762a1bSJed Brown fint[f] += cg->volume * ftmp[f]; 994*c4762a1bSJed Brown } 995*c4762a1bSJed Brown } 996*c4762a1bSJed Brown ierr = VecRestoreArrayRead(cellgeom, &cgeom);CHKERRQ(ierr); 997*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 998*c4762a1bSJed Brown ierr = MPI_Allreduce(MPI_IN_PLACE, fmin, fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 999*c4762a1bSJed Brown ierr = MPI_Allreduce(MPI_IN_PLACE, fmax, fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 1000*c4762a1bSJed Brown ierr = MPI_Allreduce(MPI_IN_PLACE, fint, fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 1001*c4762a1bSJed Brown /* Output functional data */ 1002*c4762a1bSJed Brown ftablealloc = fcount * 100; 1003*c4762a1bSJed Brown ftableused = 0; 1004*c4762a1bSJed Brown ierr = PetscCalloc1(ftablealloc, &ftable);CHKERRQ(ierr); 1005*c4762a1bSJed Brown for (f = 0; f < user->numMonitorFuncs; ++f) { 1006*c4762a1bSJed Brown Functional func = user->monitorFuncs[f]; 1007*c4762a1bSJed Brown PetscInt id = func->offset; 1008*c4762a1bSJed Brown char newline[] = "\n"; 1009*c4762a1bSJed Brown char buffer[256], *p, *prefix; 1010*c4762a1bSJed Brown size_t countused, len; 1011*c4762a1bSJed Brown 1012*c4762a1bSJed Brown /* Create string with functional outputs */ 1013*c4762a1bSJed Brown if (f % 3) { 1014*c4762a1bSJed Brown ierr = PetscArraycpy(buffer, " ", 2);CHKERRQ(ierr); 1015*c4762a1bSJed Brown p = buffer + 2; 1016*c4762a1bSJed Brown } else if (f) { 1017*c4762a1bSJed Brown ierr = PetscArraycpy(buffer, newline, sizeof(newline)-1);CHKERRQ(ierr); 1018*c4762a1bSJed Brown p = buffer + sizeof(newline) - 1; 1019*c4762a1bSJed Brown } else { 1020*c4762a1bSJed Brown p = buffer; 1021*c4762a1bSJed Brown } 1022*c4762a1bSJed Brown ierr = 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]);CHKERRQ(ierr); 1023*c4762a1bSJed Brown countused += p - buffer; 1024*c4762a1bSJed Brown /* reallocate */ 1025*c4762a1bSJed Brown if (countused > ftablealloc-ftableused-1) { 1026*c4762a1bSJed Brown char *ftablenew; 1027*c4762a1bSJed Brown 1028*c4762a1bSJed Brown ftablealloc = 2*ftablealloc + countused; 1029*c4762a1bSJed Brown ierr = PetscMalloc1(ftablealloc, &ftablenew);CHKERRQ(ierr); 1030*c4762a1bSJed Brown ierr = PetscArraycpy(ftablenew, ftable, ftableused);CHKERRQ(ierr); 1031*c4762a1bSJed Brown ierr = PetscFree(ftable);CHKERRQ(ierr); 1032*c4762a1bSJed Brown ftable = ftablenew; 1033*c4762a1bSJed Brown } 1034*c4762a1bSJed Brown ierr = PetscArraycpy(ftable+ftableused, buffer, countused);CHKERRQ(ierr); 1035*c4762a1bSJed Brown ftableused += countused; 1036*c4762a1bSJed Brown ftable[ftableused] = 0; 1037*c4762a1bSJed Brown /* Output vecs */ 1038*c4762a1bSJed Brown ierr = VecRestoreArray(fv[f], &fx[f]);CHKERRQ(ierr); 1039*c4762a1bSJed Brown ierr = PetscStrlen(func->name, &len);CHKERRQ(ierr); 1040*c4762a1bSJed Brown ierr = PetscMalloc1(len+2,&prefix);CHKERRQ(ierr); 1041*c4762a1bSJed Brown ierr = PetscStrcpy(prefix, func->name);CHKERRQ(ierr); 1042*c4762a1bSJed Brown ierr = PetscStrcat(prefix, "_");CHKERRQ(ierr); 1043*c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fv[f], prefix);CHKERRQ(ierr); 1044*c4762a1bSJed Brown ierr = VecViewFromOptions(fv[f], NULL, "-vec_view");CHKERRQ(ierr); 1045*c4762a1bSJed Brown ierr = PetscFree(prefix);CHKERRQ(ierr); 1046*c4762a1bSJed Brown ierr = DMRestoreGlobalVector(fdm[f], &fv[f]);CHKERRQ(ierr); 1047*c4762a1bSJed Brown ierr = DMDestroy(&fdm[f]);CHKERRQ(ierr); 1048*c4762a1bSJed Brown } 1049*c4762a1bSJed Brown ierr = PetscFree4(fmin, fmax, fint, ftmp);CHKERRQ(ierr); 1050*c4762a1bSJed Brown ierr = PetscFree3(fdm, fv, fx);CHKERRQ(ierr); 1051*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), "% 3D time %8.4g |x| (", stepnum, (double) time);CHKERRQ(ierr); 1052*c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 1053*c4762a1bSJed Brown if (f > 0) {ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ", ");CHKERRQ(ierr);} 1054*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+0]);CHKERRQ(ierr); 1055*c4762a1bSJed Brown } 1056*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ") |x|_1 (");CHKERRQ(ierr); 1057*c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 1058*c4762a1bSJed Brown if (f > 0) {ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ", ");CHKERRQ(ierr);} 1059*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), "%8.4g", (double) xnorms[f*2+1]);CHKERRQ(ierr); 1060*c4762a1bSJed Brown } 1061*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject) ts), ") %s\n", ftable ? ftable : "");CHKERRQ(ierr); 1062*c4762a1bSJed Brown ierr = PetscFree(ftable);CHKERRQ(ierr); 1063*c4762a1bSJed Brown } 1064*c4762a1bSJed Brown ierr = PetscFree(xnorms);CHKERRQ(ierr); 1065*c4762a1bSJed Brown PetscFunctionReturn(0); 1066*c4762a1bSJed Brown } 1067*c4762a1bSJed Brown 1068*c4762a1bSJed Brown int main(int argc, char **argv) 1069*c4762a1bSJed Brown { 1070*c4762a1bSJed Brown MPI_Comm comm; 1071*c4762a1bSJed Brown TS ts; 1072*c4762a1bSJed Brown DM dm; 1073*c4762a1bSJed Brown Vec u; 1074*c4762a1bSJed Brown AppCtx user; 1075*c4762a1bSJed Brown PetscReal t0, t = 0.0; 1076*c4762a1bSJed Brown void *ctxs[2]; 1077*c4762a1bSJed Brown PetscErrorCode ierr; 1078*c4762a1bSJed Brown 1079*c4762a1bSJed Brown ctxs[0] = &t; 1080*c4762a1bSJed Brown ctxs[1] = &t; 1081*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, (char*) 0, help);if (ierr) return ierr; 1082*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1083*c4762a1bSJed Brown user.functionalRegistry = NULL; 1084*c4762a1bSJed Brown globalUser = &user; 1085*c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 1086*c4762a1bSJed Brown ierr = TSCreate(comm, &ts);CHKERRQ(ierr); 1087*c4762a1bSJed Brown ierr = TSSetType(ts, TSBEULER);CHKERRQ(ierr); 1088*c4762a1bSJed Brown ierr = CreateDM(comm, &user, &dm);CHKERRQ(ierr); 1089*c4762a1bSJed Brown ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 1090*c4762a1bSJed Brown ierr = ProcessMonitorOptions(comm, &user);CHKERRQ(ierr); 1091*c4762a1bSJed Brown user.dm = dm; 1092*c4762a1bSJed Brown 1093*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 1094*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr); 1095*c4762a1bSJed Brown if (user.useFV) { 1096*c4762a1bSJed Brown PetscBool isImplicit = PETSC_FALSE; 1097*c4762a1bSJed Brown 1098*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,"", "-use_implicit", &isImplicit);CHKERRQ(ierr); 1099*c4762a1bSJed Brown if (isImplicit) { 1100*c4762a1bSJed Brown ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 1101*c4762a1bSJed Brown ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 1102*c4762a1bSJed Brown } 1103*c4762a1bSJed Brown ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 1104*c4762a1bSJed Brown ierr = DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, &user);CHKERRQ(ierr); 1105*c4762a1bSJed Brown } else { 1106*c4762a1bSJed Brown ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 1107*c4762a1bSJed Brown ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 1108*c4762a1bSJed Brown ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 1109*c4762a1bSJed Brown } 1110*c4762a1bSJed Brown if (user.useFV) {ierr = TSMonitorSet(ts, MonitorFunctionals, &user, NULL);CHKERRQ(ierr);} 1111*c4762a1bSJed Brown ierr = TSSetMaxSteps(ts, 1);CHKERRQ(ierr); 1112*c4762a1bSJed Brown ierr = TSSetMaxTime(ts, 2.0);CHKERRQ(ierr); 1113*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr); 1114*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 1115*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1116*c4762a1bSJed Brown 1117*c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 1118*c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, user.initialGuess, ctxs, INSERT_VALUES, u);CHKERRQ(ierr); 1119*c4762a1bSJed Brown if (user.useFV) {ierr = SetInitialConditionFVM(dm, u, 1, user.initialGuess[1], ctxs[1]);CHKERRQ(ierr);} 1120*c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-init_vec_view");CHKERRQ(ierr); 1121*c4762a1bSJed Brown ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 1122*c4762a1bSJed Brown t0 = t; 1123*c4762a1bSJed Brown ierr = DMTSCheckFromOptions(ts, u, user.exactFuncs, ctxs);CHKERRQ(ierr); 1124*c4762a1bSJed Brown ierr = TSSolve(ts, u);CHKERRQ(ierr); 1125*c4762a1bSJed Brown ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 1126*c4762a1bSJed Brown if (t > t0) {ierr = DMTSCheckFromOptions(ts, u, user.exactFuncs, ctxs);CHKERRQ(ierr);} 1127*c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 1128*c4762a1bSJed Brown { 1129*c4762a1bSJed Brown PetscReal ftime; 1130*c4762a1bSJed Brown PetscInt nsteps; 1131*c4762a1bSJed Brown TSConvergedReason reason; 1132*c4762a1bSJed Brown 1133*c4762a1bSJed Brown ierr = TSGetSolveTime(ts, &ftime);CHKERRQ(ierr); 1134*c4762a1bSJed Brown ierr = TSGetStepNumber(ts, &nsteps);CHKERRQ(ierr); 1135*c4762a1bSJed Brown ierr = TSGetConvergedReason(ts, &reason);CHKERRQ(ierr); 1136*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %D steps\n", TSConvergedReasons[reason], (double) ftime, nsteps);CHKERRQ(ierr); 1137*c4762a1bSJed Brown } 1138*c4762a1bSJed Brown 1139*c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 1140*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 1141*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 1142*c4762a1bSJed Brown ierr = PetscFree(user.monitorFuncs);CHKERRQ(ierr); 1143*c4762a1bSJed Brown ierr = FunctionalDestroy(&user.functionalRegistry);CHKERRQ(ierr); 1144*c4762a1bSJed Brown ierr = PetscFinalize(); 1145*c4762a1bSJed Brown return ierr; 1146*c4762a1bSJed Brown } 1147*c4762a1bSJed Brown 1148*c4762a1bSJed Brown /*TEST 1149*c4762a1bSJed Brown 1150*c4762a1bSJed Brown # 2D harmonic velocity, no porosity 1151*c4762a1bSJed Brown test: 1152*c4762a1bSJed Brown suffix: p1p1 1153*c4762a1bSJed Brown requires: !complex !single 1154*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1155*c4762a1bSJed Brown test: 1156*c4762a1bSJed Brown suffix: p1p1_xper 1157*c4762a1bSJed Brown requires: !complex !single 1158*c4762a1bSJed Brown args: -dm_refine 1 -y_bd_type 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 1159*c4762a1bSJed Brown test: 1160*c4762a1bSJed Brown suffix: p1p1_xper_ref 1161*c4762a1bSJed Brown requires: !complex !single 1162*c4762a1bSJed Brown args: -dm_refine 2 -y_bd_type 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 1163*c4762a1bSJed Brown test: 1164*c4762a1bSJed Brown suffix: p1p1_xyper 1165*c4762a1bSJed Brown requires: !complex !single 1166*c4762a1bSJed Brown args: -dm_refine 1 -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 1167*c4762a1bSJed Brown test: 1168*c4762a1bSJed Brown suffix: p1p1_xyper_ref 1169*c4762a1bSJed Brown requires: !complex !single 1170*c4762a1bSJed Brown args: -dm_refine 2 -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 1171*c4762a1bSJed Brown test: 1172*c4762a1bSJed Brown suffix: p2p1 1173*c4762a1bSJed Brown requires: !complex !single 1174*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1175*c4762a1bSJed Brown test: 1176*c4762a1bSJed Brown suffix: p2p1_xyper 1177*c4762a1bSJed Brown requires: !complex !single 1178*c4762a1bSJed Brown args: -dm_refine 1 -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 1179*c4762a1bSJed Brown # Must check that FV BCs propagate to coarse meshes 1180*c4762a1bSJed Brown # Must check that FV BC ids propagate to coarse meshes 1181*c4762a1bSJed Brown # Must check that FE+FV BCs work at the same time 1182*c4762a1bSJed Brown # 2D Advection, matching wind in ex11 8-11 1183*c4762a1bSJed Brown # NOTE implicit solves are limited by accuracy of FD Jacobian 1184*c4762a1bSJed Brown test: 1185*c4762a1bSJed Brown suffix: adv_0 1186*c4762a1bSJed Brown requires: !complex !single exodusii 1187*c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -x_bd_type none -y_bd_type none -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 1188*c4762a1bSJed Brown 1189*c4762a1bSJed Brown test: 1190*c4762a1bSJed Brown suffix: adv_0_im 1191*c4762a1bSJed Brown requires: !complex exodusii 1192*c4762a1bSJed Brown TODO: broken memory corruption 1193*c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -x_bd_type none -y_bd_type none -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 1194*c4762a1bSJed Brown 1195*c4762a1bSJed Brown test: 1196*c4762a1bSJed Brown suffix: adv_0_im_2 1197*c4762a1bSJed Brown requires: !complex exodusii 1198*c4762a1bSJed Brown TODO: broken 1199*c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -x_bd_type none -y_bd_type none -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 1200*c4762a1bSJed Brown 1201*c4762a1bSJed Brown test: 1202*c4762a1bSJed Brown suffix: adv_0_im_3 1203*c4762a1bSJed Brown requires: !complex exodusii 1204*c4762a1bSJed Brown TODO: broken 1205*c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -x_bd_type none -y_bd_type none -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 1206*c4762a1bSJed Brown 1207*c4762a1bSJed Brown test: 1208*c4762a1bSJed Brown suffix: adv_0_im_4 1209*c4762a1bSJed Brown requires: !complex exodusii 1210*c4762a1bSJed Brown TODO: broken 1211*c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -x_bd_type none -y_bd_type none -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 1212*c4762a1bSJed Brown # 2D Advection, misc 1213*c4762a1bSJed Brown 1214*c4762a1bSJed Brown test: 1215*c4762a1bSJed Brown suffix: adv_1 1216*c4762a1bSJed Brown requires: !complex !single 1217*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1218*c4762a1bSJed Brown 1219*c4762a1bSJed Brown test: 1220*c4762a1bSJed Brown suffix: adv_2 1221*c4762a1bSJed Brown requires: !complex 1222*c4762a1bSJed Brown TODO: broken memory corruption 1223*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type 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,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 1224*c4762a1bSJed Brown 1225*c4762a1bSJed Brown test: 1226*c4762a1bSJed Brown suffix: adv_3 1227*c4762a1bSJed Brown requires: !complex 1228*c4762a1bSJed Brown TODO: broken memory corruption 1229*c4762a1bSJed Brown args: -y_bd_type 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 1230*c4762a1bSJed Brown 1231*c4762a1bSJed Brown test: 1232*c4762a1bSJed Brown suffix: adv_3_ex 1233*c4762a1bSJed Brown requires: !complex 1234*c4762a1bSJed Brown args: -y_bd_type 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 1235*c4762a1bSJed Brown 1236*c4762a1bSJed Brown test: 1237*c4762a1bSJed Brown suffix: adv_4 1238*c4762a1bSJed Brown requires: !complex 1239*c4762a1bSJed Brown TODO: broken memory corruption 1240*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type 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 1241*c4762a1bSJed Brown 1242*c4762a1bSJed Brown # 2D Advection, box, delta 1243*c4762a1bSJed Brown test: 1244*c4762a1bSJed Brown suffix: adv_delta_yper_0 1245*c4762a1bSJed Brown requires: !complex 1246*c4762a1bSJed Brown TODO: broken 1247*c4762a1bSJed Brown args: -x_bd_type none -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 1248*c4762a1bSJed Brown 1249*c4762a1bSJed Brown test: 1250*c4762a1bSJed Brown suffix: adv_delta_yper_1 1251*c4762a1bSJed Brown requires: !complex 1252*c4762a1bSJed Brown TODO: broken 1253*c4762a1bSJed Brown args: -x_bd_type none -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 1254*c4762a1bSJed Brown 1255*c4762a1bSJed Brown test: 1256*c4762a1bSJed Brown suffix: adv_delta_yper_2 1257*c4762a1bSJed Brown requires: !complex 1258*c4762a1bSJed Brown TODO: broken 1259*c4762a1bSJed Brown args: -x_bd_type none -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 1260*c4762a1bSJed Brown 1261*c4762a1bSJed Brown test: 1262*c4762a1bSJed Brown suffix: adv_delta_yper_fim_0 1263*c4762a1bSJed Brown requires: !complex 1264*c4762a1bSJed Brown TODO: broken 1265*c4762a1bSJed Brown args: -x_bd_type none -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 1266*c4762a1bSJed Brown 1267*c4762a1bSJed Brown test: 1268*c4762a1bSJed Brown suffix: adv_delta_yper_fim_1 1269*c4762a1bSJed Brown requires: !complex 1270*c4762a1bSJed Brown TODO: broken 1271*c4762a1bSJed Brown args: -x_bd_type none -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 1272*c4762a1bSJed Brown 1273*c4762a1bSJed Brown test: 1274*c4762a1bSJed Brown suffix: adv_delta_yper_fim_2 1275*c4762a1bSJed Brown requires: !complex 1276*c4762a1bSJed Brown TODO: broken 1277*c4762a1bSJed Brown args: -x_bd_type none -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 1278*c4762a1bSJed Brown 1279*c4762a1bSJed Brown test: 1280*c4762a1bSJed Brown suffix: adv_delta_yper_im_0 1281*c4762a1bSJed Brown requires: !complex 1282*c4762a1bSJed Brown TODO: broken 1283*c4762a1bSJed Brown args: -x_bd_type none -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 1284*c4762a1bSJed Brown 1285*c4762a1bSJed Brown test: 1286*c4762a1bSJed Brown suffix: adv_delta_yper_im_1 1287*c4762a1bSJed Brown requires: !complex 1288*c4762a1bSJed Brown TODO: broken 1289*c4762a1bSJed Brown args: -x_bd_type none -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 1290*c4762a1bSJed Brown 1291*c4762a1bSJed Brown test: 1292*c4762a1bSJed Brown suffix: adv_delta_yper_im_2 1293*c4762a1bSJed Brown requires: !complex 1294*c4762a1bSJed Brown TODO: broken 1295*c4762a1bSJed Brown args: -x_bd_type none -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 1296*c4762a1bSJed Brown 1297*c4762a1bSJed Brown test: 1298*c4762a1bSJed Brown suffix: adv_delta_yper_im_3 1299*c4762a1bSJed Brown requires: !complex 1300*c4762a1bSJed Brown TODO: broken 1301*c4762a1bSJed Brown args: -x_bd_type none -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 1302*c4762a1bSJed Brown 1303*c4762a1bSJed Brown # I believe the nullspace is sin(pi y) 1304*c4762a1bSJed Brown test: 1305*c4762a1bSJed Brown suffix: adv_delta_yper_im_4 1306*c4762a1bSJed Brown requires: !complex 1307*c4762a1bSJed Brown TODO: broken 1308*c4762a1bSJed Brown args: -x_bd_type none -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 1309*c4762a1bSJed Brown 1310*c4762a1bSJed Brown test: 1311*c4762a1bSJed Brown suffix: adv_delta_yper_im_5 1312*c4762a1bSJed Brown requires: !complex 1313*c4762a1bSJed Brown TODO: broken 1314*c4762a1bSJed Brown args: -x_bd_type none -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 1315*c4762a1bSJed Brown 1316*c4762a1bSJed Brown test: 1317*c4762a1bSJed Brown suffix: adv_delta_yper_im_6 1318*c4762a1bSJed Brown requires: !complex 1319*c4762a1bSJed Brown TODO: broken 1320*c4762a1bSJed Brown args: -x_bd_type none -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 1321*c4762a1bSJed Brown # 2D Advection, magma benchmark 1 1322*c4762a1bSJed Brown 1323*c4762a1bSJed Brown test: 1324*c4762a1bSJed Brown suffix: adv_delta_shear_im_0 1325*c4762a1bSJed Brown requires: !complex 1326*c4762a1bSJed Brown TODO: broken 1327*c4762a1bSJed Brown args: -y_bd_type 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 1328*c4762a1bSJed Brown # 2D Advection, box, gaussian 1329*c4762a1bSJed Brown 1330*c4762a1bSJed Brown test: 1331*c4762a1bSJed Brown suffix: adv_gauss 1332*c4762a1bSJed Brown requires: !complex 1333*c4762a1bSJed Brown TODO: broken 1334*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1335*c4762a1bSJed Brown 1336*c4762a1bSJed Brown test: 1337*c4762a1bSJed Brown suffix: adv_gauss_im 1338*c4762a1bSJed Brown requires: !complex 1339*c4762a1bSJed Brown TODO: broken 1340*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1341*c4762a1bSJed Brown 1342*c4762a1bSJed Brown test: 1343*c4762a1bSJed Brown suffix: adv_gauss_im_1 1344*c4762a1bSJed Brown requires: !complex 1345*c4762a1bSJed Brown TODO: broken 1346*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1347*c4762a1bSJed Brown 1348*c4762a1bSJed Brown test: 1349*c4762a1bSJed Brown suffix: adv_gauss_im_2 1350*c4762a1bSJed Brown requires: !complex 1351*c4762a1bSJed Brown TODO: broken 1352*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1353*c4762a1bSJed Brown 1354*c4762a1bSJed Brown test: 1355*c4762a1bSJed Brown suffix: adv_gauss_corner 1356*c4762a1bSJed Brown requires: !complex 1357*c4762a1bSJed Brown TODO: broken 1358*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1359*c4762a1bSJed Brown 1360*c4762a1bSJed Brown # 2D Advection+Harmonic 12- 1361*c4762a1bSJed Brown test: 1362*c4762a1bSJed Brown suffix: adv_harm_0 1363*c4762a1bSJed Brown requires: !complex 1364*c4762a1bSJed Brown TODO: broken memory corruption 1365*c4762a1bSJed Brown args: -x_bd_type none -y_bd_type none -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 1366*c4762a1bSJed Brown 1367*c4762a1bSJed Brown TEST*/ 1368