1d94ba093SMatthew G. Knepley static char help[] = "Poisson Problem with a split domain.\n\ 2d94ba093SMatthew G. Knepley We solve the Poisson problem in two halves of a domain.\n\ 3d94ba093SMatthew G. Knepley In one half, we include an additional field.\n\n\n"; 4d94ba093SMatthew G. Knepley 5d94ba093SMatthew G. Knepley #include <petscdmplex.h> 6d94ba093SMatthew G. Knepley #include <petscsnes.h> 7d94ba093SMatthew G. Knepley #include <petscds.h> 8d94ba093SMatthew G. Knepley #include <petscconvest.h> 9d94ba093SMatthew G. Knepley 10d94ba093SMatthew G. Knepley typedef struct { 11d94ba093SMatthew G. Knepley /* Domain and mesh definition */ 12d94ba093SMatthew G. Knepley PetscInt dim; /* The topological mesh dimension */ 13d94ba093SMatthew G. Knepley PetscBool simplex; /* Simplicial mesh */ 14d94ba093SMatthew G. Knepley } AppCtx; 15d94ba093SMatthew G. Knepley 16d94ba093SMatthew G. Knepley static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 17d94ba093SMatthew G. Knepley { 18d94ba093SMatthew G. Knepley u[0] = x[0]*x[0] + x[1]*x[1]; 19d94ba093SMatthew G. Knepley return 0; 20d94ba093SMatthew G. Knepley } 21d94ba093SMatthew G. Knepley 22d94ba093SMatthew G. Knepley static PetscErrorCode quad_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 23d94ba093SMatthew G. Knepley { 24d94ba093SMatthew G. Knepley u[0] = 2.0; 25d94ba093SMatthew G. Knepley return 0; 26d94ba093SMatthew G. Knepley } 27d94ba093SMatthew G. Knepley 28d94ba093SMatthew G. Knepley static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 29d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 30d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 31d94ba093SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 32d94ba093SMatthew G. Knepley { 33d94ba093SMatthew G. Knepley PetscInt d; 34d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += 2.0; 35d94ba093SMatthew G. Knepley } 36d94ba093SMatthew G. Knepley 37d94ba093SMatthew G. Knepley static void f0_quad_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 38d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 39d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 40d94ba093SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 41d94ba093SMatthew G. Knepley { 42d94ba093SMatthew G. Knepley f0[0] = u[uOff[1]] - 2.0; 43d94ba093SMatthew G. Knepley } 44d94ba093SMatthew G. Knepley 45d94ba093SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 46d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 47d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 48d94ba093SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 49d94ba093SMatthew G. Knepley { 50d94ba093SMatthew G. Knepley PetscInt d; 51d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 52d94ba093SMatthew G. Knepley } 53d94ba093SMatthew G. Knepley 54d94ba093SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 55d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 56d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 57d94ba093SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 58d94ba093SMatthew G. Knepley { 59d94ba093SMatthew G. Knepley PetscInt d; 60d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 61d94ba093SMatthew G. Knepley } 62d94ba093SMatthew G. Knepley 63d94ba093SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 64d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 65d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 66d94ba093SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 67d94ba093SMatthew G. Knepley { 68d94ba093SMatthew G. Knepley g0[0] = 1.0; 69d94ba093SMatthew G. Knepley } 70d94ba093SMatthew G. Knepley 71d94ba093SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 72d94ba093SMatthew G. Knepley { 73d94ba093SMatthew G. Knepley PetscErrorCode ierr; 74d94ba093SMatthew G. Knepley 75d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 76d94ba093SMatthew G. Knepley options->dim = 2; 77d94ba093SMatthew G. Knepley options->simplex = PETSC_TRUE; 78d94ba093SMatthew G. Knepley 79d94ba093SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr); 80d94ba093SMatthew G. Knepley ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex23.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 81d94ba093SMatthew G. Knepley ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex23.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 82d94ba093SMatthew G. Knepley ierr = PetscOptionsEnd(); 83d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 84d94ba093SMatthew G. Knepley } 85d94ba093SMatthew G. Knepley 86d94ba093SMatthew G. Knepley static PetscErrorCode DivideDomain(DM dm, AppCtx *user) 87d94ba093SMatthew G. Knepley { 88d94ba093SMatthew G. Knepley DMLabel top, bottom; 89d94ba093SMatthew G. Knepley PetscReal low[3], high[3], midy; 90d94ba093SMatthew G. Knepley PetscInt cStart, cEnd, c; 91d94ba093SMatthew G. Knepley PetscErrorCode ierr; 92d94ba093SMatthew G. Knepley 93d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 94d94ba093SMatthew G. Knepley ierr = DMCreateLabel(dm, "top");CHKERRQ(ierr); 95d94ba093SMatthew G. Knepley ierr = DMCreateLabel(dm, "bottom");CHKERRQ(ierr); 96d94ba093SMatthew G. Knepley ierr = DMGetLabel(dm, "top", &top);CHKERRQ(ierr); 97d94ba093SMatthew G. Knepley ierr = DMGetLabel(dm, "bottom", &bottom);CHKERRQ(ierr); 98d94ba093SMatthew G. Knepley ierr = DMGetBoundingBox(dm, low, high);CHKERRQ(ierr); 99d94ba093SMatthew G. Knepley midy = 0.5*(high[1] - low[1]); 100d94ba093SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 101d94ba093SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 102d94ba093SMatthew G. Knepley PetscReal centroid[3]; 103d94ba093SMatthew G. Knepley 104d94ba093SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr); 105d94ba093SMatthew G. Knepley if (centroid[1] > midy) {ierr = DMLabelSetValue(top, c, 1);CHKERRQ(ierr);} 106d94ba093SMatthew G. Knepley else {ierr = DMLabelSetValue(bottom, c, 1);CHKERRQ(ierr);} 107d94ba093SMatthew G. Knepley } 108d94ba093SMatthew G. Knepley ierr = DMPlexLabelComplete(dm, top);CHKERRQ(ierr); 109d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 110d94ba093SMatthew G. Knepley } 111d94ba093SMatthew G. Knepley 112d94ba093SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 113d94ba093SMatthew G. Knepley { 114d94ba093SMatthew G. Knepley PetscErrorCode ierr; 115d94ba093SMatthew G. Knepley 116d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 117d94ba093SMatthew G. Knepley /* Create box mesh */ 118d94ba093SMatthew G. Knepley ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 119d94ba093SMatthew G. Knepley /* Distribute mesh over processes */ 120d94ba093SMatthew G. Knepley { 121d94ba093SMatthew G. Knepley DM dmDist = NULL; 122d94ba093SMatthew G. Knepley PetscPartitioner part; 123d94ba093SMatthew G. Knepley 124d94ba093SMatthew G. Knepley ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 125d94ba093SMatthew G. Knepley ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 126d94ba093SMatthew G. Knepley ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr); 127d94ba093SMatthew G. Knepley if (dmDist) { 128d94ba093SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 129d94ba093SMatthew G. Knepley *dm = dmDist; 130d94ba093SMatthew G. Knepley } 131d94ba093SMatthew G. Knepley } 132d94ba093SMatthew G. Knepley /* TODO: This should be pulled into the library */ 133d94ba093SMatthew G. Knepley { 134d94ba093SMatthew G. Knepley char convType[256]; 135d94ba093SMatthew G. Knepley PetscBool flg; 136d94ba093SMatthew G. Knepley 137d94ba093SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 138d94ba093SMatthew G. Knepley ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr); 139d94ba093SMatthew G. Knepley ierr = PetscOptionsEnd(); 140d94ba093SMatthew G. Knepley if (flg) { 141d94ba093SMatthew G. Knepley DM dmConv; 142d94ba093SMatthew G. Knepley 143d94ba093SMatthew G. Knepley ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr); 144d94ba093SMatthew G. Knepley if (dmConv) { 145d94ba093SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 146d94ba093SMatthew G. Knepley *dm = dmConv; 147d94ba093SMatthew G. Knepley } 148d94ba093SMatthew G. Knepley } 149d94ba093SMatthew G. Knepley } 150d94ba093SMatthew G. Knepley ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 151d94ba093SMatthew G. Knepley 152d94ba093SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 153d94ba093SMatthew G. Knepley ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 154d94ba093SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 155d94ba093SMatthew G. Knepley ierr = DivideDomain(*dm, user);CHKERRQ(ierr); 156d94ba093SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 157d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 158d94ba093SMatthew G. Knepley } 159d94ba093SMatthew G. Knepley 160d94ba093SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 161d94ba093SMatthew G. Knepley { 162d94ba093SMatthew G. Knepley PetscDS ds; 163*76fb4698SMatthew G. Knepley PetscWeakForm wf; 164*76fb4698SMatthew G. Knepley DMLabel label; 165d94ba093SMatthew G. Knepley const PetscInt id = 1; 166d94ba093SMatthew G. Knepley PetscErrorCode ierr; 167d94ba093SMatthew G. Knepley 168d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 169*76fb4698SMatthew G. Knepley ierr = DMGetRegionNumDS(dm, 0, &label, NULL, &ds);CHKERRQ(ierr); 170*76fb4698SMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 171*76fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, f0_quad_u, 0, f1_u);CHKERRQ(ierr); 172*76fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu);CHKERRQ(ierr); 173d94ba093SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr); 174*76fb4698SMatthew G. Knepley ierr = DMGetRegionNumDS(dm, 1, &label, NULL, &ds);CHKERRQ(ierr); 175*76fb4698SMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 176*76fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, f0_quad_u, 0, f1_u);CHKERRQ(ierr); 177*76fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu);CHKERRQ(ierr); 178*76fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, f0_quad_p, 0, NULL);CHKERRQ(ierr); 179*76fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL);CHKERRQ(ierr); 180d94ba093SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr); 181d94ba093SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 1, quad_p, user);CHKERRQ(ierr); 18256cf3b9cSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) quad_u, NULL, 1, &id, user);CHKERRQ(ierr); 183d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 184d94ba093SMatthew G. Knepley } 185d94ba093SMatthew G. Knepley 186d94ba093SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 187d94ba093SMatthew G. Knepley { 188d94ba093SMatthew G. Knepley DM cdm = dm; 189d94ba093SMatthew G. Knepley DMLabel top; 190d94ba093SMatthew G. Knepley PetscFE fe, feTop; 191d94ba093SMatthew G. Knepley PetscQuadrature q; 192d94ba093SMatthew G. Knepley const char *nameTop = "pressure"; 193d94ba093SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 194d94ba093SMatthew G. Knepley PetscErrorCode ierr; 195d94ba093SMatthew G. Knepley 196d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 197d94ba093SMatthew G. Knepley ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 198d94ba093SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 199d94ba093SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 200d94ba093SMatthew G. Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 201d94ba093SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 202d94ba093SMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 203d94ba093SMatthew G. Knepley ierr = DMGetLabel(dm, "top", &top);CHKERRQ(ierr); 204d94ba093SMatthew G. Knepley ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop);CHKERRQ(ierr); 205d94ba093SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &feTop);CHKERRQ(ierr); 206d94ba093SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) feTop, nameTop);CHKERRQ(ierr); 207d94ba093SMatthew G. Knepley ierr = PetscFESetQuadrature(feTop, q);CHKERRQ(ierr); 208d94ba093SMatthew G. Knepley ierr = DMSetField(dm, 1, top, (PetscObject) feTop);CHKERRQ(ierr); 209d94ba093SMatthew G. Knepley ierr = PetscFEDestroy(&feTop);CHKERRQ(ierr); 210d94ba093SMatthew G. Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 211d94ba093SMatthew G. Knepley ierr = (*setup)(dm, user);CHKERRQ(ierr); 212d94ba093SMatthew G. Knepley while (cdm) { 213d94ba093SMatthew G. Knepley ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 214d94ba093SMatthew G. Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 215d94ba093SMatthew G. Knepley } 216d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 217d94ba093SMatthew G. Knepley } 218d94ba093SMatthew G. Knepley 219d94ba093SMatthew G. Knepley int main(int argc, char **argv) 220d94ba093SMatthew G. Knepley { 221d94ba093SMatthew G. Knepley DM dm; /* Problem specification */ 222d94ba093SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 223d94ba093SMatthew G. Knepley Vec u; /* Solutions */ 224d94ba093SMatthew G. Knepley AppCtx user; /* User-defined work context */ 225d94ba093SMatthew G. Knepley PetscErrorCode ierr; 226d94ba093SMatthew G. Knepley 227d94ba093SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 228d94ba093SMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 229d94ba093SMatthew G. Knepley /* Primal system */ 230d94ba093SMatthew G. Knepley ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 231d94ba093SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 232d94ba093SMatthew G. Knepley ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 233d94ba093SMatthew G. Knepley ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr); 234d94ba093SMatthew G. Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 235d94ba093SMatthew G. Knepley ierr = VecSet(u, 0.0);CHKERRQ(ierr); 236d94ba093SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr); 237d94ba093SMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 238d94ba093SMatthew G. Knepley ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 239348a1646SMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 240d94ba093SMatthew G. Knepley ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 241d94ba093SMatthew G. Knepley ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 242d94ba093SMatthew G. Knepley ierr = VecViewFromOptions(u, NULL, "-sol_view");CHKERRQ(ierr); 243d94ba093SMatthew G. Knepley /* Cleanup */ 244d94ba093SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 245d94ba093SMatthew G. Knepley ierr = SNESDestroy(&snes);CHKERRQ(ierr); 246d94ba093SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 247d94ba093SMatthew G. Knepley ierr = PetscFinalize(); 248d94ba093SMatthew G. Knepley return ierr; 249d94ba093SMatthew G. Knepley } 250d94ba093SMatthew G. Knepley 251d94ba093SMatthew G. Knepley /*TEST 252d94ba093SMatthew G. Knepley 253d94ba093SMatthew G. Knepley test: 254d94ba093SMatthew G. Knepley suffix: 2d_p1_0 255d94ba093SMatthew G. Knepley requires: triangle 256d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 2 -pressure_top_petscspace_degree 0 -dm_refine 0 -dmsnes_check 257d94ba093SMatthew G. Knepley 258d94ba093SMatthew G. Knepley test: 259d94ba093SMatthew G. Knepley suffix: 2d_p1_1 260d94ba093SMatthew G. Knepley requires: triangle 261d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 1 -pressure_top_petscspace_degree 0 -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate 262d94ba093SMatthew G. Knepley 263d94ba093SMatthew G. Knepley TEST*/ 264