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 { 1130602db0SMatthew 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 75d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 76*9566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "top")); 77*9566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "bottom")); 78*9566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "top", &top)); 79*9566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "bottom", &bottom)); 80*9566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, low, high)); 81d94ba093SMatthew G. Knepley midy = 0.5*(high[1] - low[1]); 82*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 83d94ba093SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 84d94ba093SMatthew G. Knepley PetscReal centroid[3]; 85d94ba093SMatthew G. Knepley 86*9566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 87*9566063dSJacob Faibussowitsch if (centroid[1] > midy) PetscCall(DMLabelSetValue(top, c, 1)); 88*9566063dSJacob Faibussowitsch else PetscCall(DMLabelSetValue(bottom, c, 1)); 89d94ba093SMatthew G. Knepley } 90*9566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(dm, top)); 91d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 92d94ba093SMatthew G. Knepley } 93d94ba093SMatthew G. Knepley 94d94ba093SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 95d94ba093SMatthew G. Knepley { 96d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 97*9566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 98*9566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 99*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 100*9566063dSJacob Faibussowitsch PetscCall(DivideDomain(*dm, user)); 101*9566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 102*9566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 103d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 104d94ba093SMatthew G. Knepley } 105d94ba093SMatthew G. Knepley 106d94ba093SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 107d94ba093SMatthew G. Knepley { 108d94ba093SMatthew G. Knepley PetscDS ds; 10976fb4698SMatthew G. Knepley PetscWeakForm wf; 11076fb4698SMatthew G. Knepley DMLabel label; 111d94ba093SMatthew G. Knepley const PetscInt id = 1; 112d94ba093SMatthew G. Knepley 113d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 114*9566063dSJacob Faibussowitsch PetscCall(DMGetRegionNumDS(dm, 0, &label, NULL, &ds)); 115*9566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 116*9566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u)); 117*9566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu)); 118*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user)); 119*9566063dSJacob Faibussowitsch PetscCall(DMGetRegionNumDS(dm, 1, &label, NULL, &ds)); 120*9566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 121*9566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u)); 122*9566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu)); 123*9566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, 0, f0_quad_p, 0, NULL)); 124*9566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL)); 125*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user)); 126*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quad_p, user)); 127*9566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 128*9566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL)); 129d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 130d94ba093SMatthew G. Knepley } 131d94ba093SMatthew G. Knepley 132d94ba093SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 133d94ba093SMatthew G. Knepley { 134d94ba093SMatthew G. Knepley DM cdm = dm; 135d94ba093SMatthew G. Knepley DMLabel top; 136d94ba093SMatthew G. Knepley PetscFE fe, feTop; 137d94ba093SMatthew G. Knepley PetscQuadrature q; 13830602db0SMatthew G. Knepley PetscInt dim; 13930602db0SMatthew G. Knepley PetscBool simplex; 140d94ba093SMatthew G. Knepley const char *nameTop = "pressure"; 141d94ba093SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 142d94ba093SMatthew G. Knepley 143d94ba093SMatthew G. Knepley PetscFunctionBeginUser; 144*9566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 145*9566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 146*9566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 147*9566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 148*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 149*9566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 150*9566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &q)); 151*9566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 152*9566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "top", &top)); 153*9566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop)); 154*9566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &feTop)); 155*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) feTop, nameTop)); 156*9566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(feTop, q)); 157*9566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, top, (PetscObject) feTop)); 158*9566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feTop)); 159*9566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 160*9566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 161d94ba093SMatthew G. Knepley while (cdm) { 162*9566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 163*9566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 164d94ba093SMatthew G. Knepley } 165d94ba093SMatthew G. Knepley PetscFunctionReturn(0); 166d94ba093SMatthew G. Knepley } 167d94ba093SMatthew G. Knepley 168d94ba093SMatthew G. Knepley int main(int argc, char **argv) 169d94ba093SMatthew G. Knepley { 170d94ba093SMatthew G. Knepley DM dm; /* Problem specification */ 171d94ba093SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 172d94ba093SMatthew G. Knepley Vec u; /* Solutions */ 173d94ba093SMatthew G. Knepley AppCtx user; /* User-defined work context */ 174d94ba093SMatthew G. Knepley 175*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 176d94ba093SMatthew G. Knepley /* Primal system */ 177*9566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 178*9566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 179*9566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 180*9566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 181*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 182*9566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 183*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "solution")); 184*9566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 185*9566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 186*9566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 187*9566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 188*9566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 189*9566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_view")); 190d94ba093SMatthew G. Knepley /* Cleanup */ 191*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 192*9566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 193*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 194*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 195b122ec5aSJacob Faibussowitsch return 0; 196d94ba093SMatthew G. Knepley } 197d94ba093SMatthew G. Knepley 198d94ba093SMatthew G. Knepley /*TEST 199d94ba093SMatthew G. Knepley 200d94ba093SMatthew G. Knepley test: 201d94ba093SMatthew G. Knepley suffix: 2d_p1_0 202d94ba093SMatthew G. Knepley requires: triangle 203d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 2 -pressure_top_petscspace_degree 0 -dm_refine 0 -dmsnes_check 204d94ba093SMatthew G. Knepley 205d94ba093SMatthew G. Knepley test: 206d94ba093SMatthew G. Knepley suffix: 2d_p1_1 207d94ba093SMatthew G. Knepley requires: triangle 208d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 1 -pressure_top_petscspace_degree 0 -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate 209d94ba093SMatthew G. Knepley 210d94ba093SMatthew G. Knepley TEST*/ 211