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 { 11*30602db0SMatthew G. Knepley PetscInt dummy; 12d94ba093SMatthew G. Knepley } AppCtx; 13d94ba093SMatthew G. Knepley 14d94ba093SMatthew G. Knepley static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 15d94ba093SMatthew G. Knepley { 16d94ba093SMatthew G. Knepley u[0] = x[0]*x[0] + x[1]*x[1]; 17d94ba093SMatthew G. Knepley return 0; 18d94ba093SMatthew G. Knepley } 19d94ba093SMatthew G. Knepley 20d94ba093SMatthew G. Knepley static PetscErrorCode quad_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 21d94ba093SMatthew G. Knepley { 22d94ba093SMatthew G. Knepley u[0] = 2.0; 23d94ba093SMatthew G. Knepley return 0; 24d94ba093SMatthew G. Knepley } 25d94ba093SMatthew G. Knepley 26d94ba093SMatthew G. Knepley static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 27d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 28d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 29d94ba093SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 30d94ba093SMatthew G. Knepley { 31d94ba093SMatthew G. Knepley PetscInt d; 32d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += 2.0; 33d94ba093SMatthew G. Knepley } 34d94ba093SMatthew G. Knepley 35d94ba093SMatthew G. Knepley static void f0_quad_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 36d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 37d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 38d94ba093SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 39d94ba093SMatthew G. Knepley { 40d94ba093SMatthew G. Knepley f0[0] = u[uOff[1]] - 2.0; 41d94ba093SMatthew G. Knepley } 42d94ba093SMatthew G. Knepley 43d94ba093SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 44d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 45d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 46d94ba093SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 47d94ba093SMatthew G. Knepley { 48d94ba093SMatthew G. Knepley PetscInt d; 49d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 50d94ba093SMatthew G. Knepley } 51d94ba093SMatthew G. Knepley 52d94ba093SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 53d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 54d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 55d94ba093SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 56d94ba093SMatthew G. Knepley { 57d94ba093SMatthew G. Knepley PetscInt d; 58d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 59d94ba093SMatthew G. Knepley } 60d94ba093SMatthew G. Knepley 61d94ba093SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 62d94ba093SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 63d94ba093SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 64d94ba093SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 65d94ba093SMatthew G. Knepley { 66d94ba093SMatthew G. Knepley g0[0] = 1.0; 67d94ba093SMatthew G. Knepley } 68d94ba093SMatthew G. Knepley 69d94ba093SMatthew G. Knepley static PetscErrorCode DivideDomain(DM dm, AppCtx *user) 70d94ba093SMatthew G. Knepley { 71d94ba093SMatthew G. Knepley DMLabel top, bottom; 72d94ba093SMatthew G. Knepley PetscReal low[3], high[3], midy; 73d94ba093SMatthew G. Knepley PetscInt cStart, cEnd, c; 74d94ba093SMatthew G. Knepley PetscErrorCode ierr; 75d94ba093SMatthew G. Knepley 76d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 77d94ba093SMatthew G. Knepley ierr = DMCreateLabel(dm, "top");CHKERRQ(ierr); 78d94ba093SMatthew G. Knepley ierr = DMCreateLabel(dm, "bottom");CHKERRQ(ierr); 79d94ba093SMatthew G. Knepley ierr = DMGetLabel(dm, "top", &top);CHKERRQ(ierr); 80d94ba093SMatthew G. Knepley ierr = DMGetLabel(dm, "bottom", &bottom);CHKERRQ(ierr); 81d94ba093SMatthew G. Knepley ierr = DMGetBoundingBox(dm, low, high);CHKERRQ(ierr); 82d94ba093SMatthew G. Knepley midy = 0.5*(high[1] - low[1]); 83d94ba093SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 84d94ba093SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 85d94ba093SMatthew G. Knepley PetscReal centroid[3]; 86d94ba093SMatthew G. Knepley 87d94ba093SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr); 88d94ba093SMatthew G. Knepley if (centroid[1] > midy) {ierr = DMLabelSetValue(top, c, 1);CHKERRQ(ierr);} 89d94ba093SMatthew G. Knepley else {ierr = DMLabelSetValue(bottom, c, 1);CHKERRQ(ierr);} 90d94ba093SMatthew G. Knepley } 91d94ba093SMatthew G. Knepley ierr = DMPlexLabelComplete(dm, top);CHKERRQ(ierr); 92d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 93d94ba093SMatthew G. Knepley } 94d94ba093SMatthew G. Knepley 95d94ba093SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 96d94ba093SMatthew G. Knepley { 97d94ba093SMatthew G. Knepley PetscErrorCode ierr; 98d94ba093SMatthew G. Knepley 99d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 100*30602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 101*30602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 102d94ba093SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 103d94ba093SMatthew G. Knepley ierr = DivideDomain(*dm, user);CHKERRQ(ierr); 104*30602db0SMatthew G. Knepley ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 105d94ba093SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 106d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 107d94ba093SMatthew G. Knepley } 108d94ba093SMatthew G. Knepley 109d94ba093SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 110d94ba093SMatthew G. Knepley { 111d94ba093SMatthew G. Knepley PetscDS ds; 11276fb4698SMatthew G. Knepley PetscWeakForm wf; 11376fb4698SMatthew G. Knepley DMLabel label; 114d94ba093SMatthew G. Knepley const PetscInt id = 1; 115d94ba093SMatthew G. Knepley PetscErrorCode ierr; 116d94ba093SMatthew G. Knepley 117d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 11876fb4698SMatthew G. Knepley ierr = DMGetRegionNumDS(dm, 0, &label, NULL, &ds);CHKERRQ(ierr); 11976fb4698SMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 12076fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, f0_quad_u, 0, f1_u);CHKERRQ(ierr); 12176fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu);CHKERRQ(ierr); 122d94ba093SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr); 12376fb4698SMatthew G. Knepley ierr = DMGetRegionNumDS(dm, 1, &label, NULL, &ds);CHKERRQ(ierr); 12476fb4698SMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 12576fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, f0_quad_u, 0, f1_u);CHKERRQ(ierr); 12676fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu);CHKERRQ(ierr); 12776fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, f0_quad_p, 0, NULL);CHKERRQ(ierr); 12876fb4698SMatthew G. Knepley ierr = PetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL);CHKERRQ(ierr); 129d94ba093SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr); 130d94ba093SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 1, quad_p, user);CHKERRQ(ierr); 13145480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 13245480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL);CHKERRQ(ierr); 133d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 134d94ba093SMatthew G. Knepley } 135d94ba093SMatthew G. Knepley 136d94ba093SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 137d94ba093SMatthew G. Knepley { 138d94ba093SMatthew G. Knepley DM cdm = dm; 139d94ba093SMatthew G. Knepley DMLabel top; 140d94ba093SMatthew G. Knepley PetscFE fe, feTop; 141d94ba093SMatthew G. Knepley PetscQuadrature q; 142*30602db0SMatthew G. Knepley PetscInt dim; 143*30602db0SMatthew G. Knepley PetscBool simplex; 144d94ba093SMatthew G. Knepley const char *nameTop = "pressure"; 145d94ba093SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 146d94ba093SMatthew G. Knepley PetscErrorCode ierr; 147d94ba093SMatthew G. Knepley 148d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 149*30602db0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 150*30602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 151d94ba093SMatthew G. Knepley ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 152*30602db0SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 153d94ba093SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 154d94ba093SMatthew G. Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 155d94ba093SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 156d94ba093SMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 157d94ba093SMatthew G. Knepley ierr = DMGetLabel(dm, "top", &top);CHKERRQ(ierr); 158d94ba093SMatthew G. Knepley ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop);CHKERRQ(ierr); 159*30602db0SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &feTop);CHKERRQ(ierr); 160d94ba093SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) feTop, nameTop);CHKERRQ(ierr); 161d94ba093SMatthew G. Knepley ierr = PetscFESetQuadrature(feTop, q);CHKERRQ(ierr); 162d94ba093SMatthew G. Knepley ierr = DMSetField(dm, 1, top, (PetscObject) feTop);CHKERRQ(ierr); 163d94ba093SMatthew G. Knepley ierr = PetscFEDestroy(&feTop);CHKERRQ(ierr); 164d94ba093SMatthew G. Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 165d94ba093SMatthew G. Knepley ierr = (*setup)(dm, user);CHKERRQ(ierr); 166d94ba093SMatthew G. Knepley while (cdm) { 167d94ba093SMatthew G. Knepley ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 168d94ba093SMatthew G. Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 169d94ba093SMatthew G. Knepley } 170d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 171d94ba093SMatthew G. Knepley } 172d94ba093SMatthew G. Knepley 173d94ba093SMatthew G. Knepley int main(int argc, char **argv) 174d94ba093SMatthew G. Knepley { 175d94ba093SMatthew G. Knepley DM dm; /* Problem specification */ 176d94ba093SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 177d94ba093SMatthew G. Knepley Vec u; /* Solutions */ 178d94ba093SMatthew G. Knepley AppCtx user; /* User-defined work context */ 179d94ba093SMatthew G. Knepley PetscErrorCode ierr; 180d94ba093SMatthew G. Knepley 181d94ba093SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 182d94ba093SMatthew G. Knepley /* Primal system */ 183d94ba093SMatthew G. Knepley ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 184d94ba093SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 185d94ba093SMatthew G. Knepley ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 186d94ba093SMatthew G. Knepley ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr); 187d94ba093SMatthew G. Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 188d94ba093SMatthew G. Knepley ierr = VecSet(u, 0.0);CHKERRQ(ierr); 189d94ba093SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr); 190d94ba093SMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 191d94ba093SMatthew G. Knepley ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 192348a1646SMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 193d94ba093SMatthew G. Knepley ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 194d94ba093SMatthew G. Knepley ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 195d94ba093SMatthew G. Knepley ierr = VecViewFromOptions(u, NULL, "-sol_view");CHKERRQ(ierr); 196d94ba093SMatthew G. Knepley /* Cleanup */ 197d94ba093SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 198d94ba093SMatthew G. Knepley ierr = SNESDestroy(&snes);CHKERRQ(ierr); 199d94ba093SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 200d94ba093SMatthew G. Knepley ierr = PetscFinalize(); 201d94ba093SMatthew G. Knepley return ierr; 202d94ba093SMatthew G. Knepley } 203d94ba093SMatthew G. Knepley 204d94ba093SMatthew G. Knepley /*TEST 205d94ba093SMatthew G. Knepley 206d94ba093SMatthew G. Knepley test: 207d94ba093SMatthew G. Knepley suffix: 2d_p1_0 208d94ba093SMatthew G. Knepley requires: triangle 209d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 2 -pressure_top_petscspace_degree 0 -dm_refine 0 -dmsnes_check 210d94ba093SMatthew G. Knepley 211d94ba093SMatthew G. Knepley test: 212d94ba093SMatthew G. Knepley suffix: 2d_p1_1 213d94ba093SMatthew G. Knepley requires: triangle 214d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 1 -pressure_top_petscspace_degree 0 -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate 215d94ba093SMatthew G. Knepley 216d94ba093SMatthew G. Knepley TEST*/ 217