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
quad_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)14*2a8381b2SBarry Smith static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
15d71ae5a4SJacob Faibussowitsch {
16d94ba093SMatthew G. Knepley u[0] = x[0] * x[0] + x[1] * x[1];
173ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
18d94ba093SMatthew G. Knepley }
19d94ba093SMatthew G. Knepley
quad_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)20*2a8381b2SBarry Smith static PetscErrorCode quad_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
21d71ae5a4SJacob Faibussowitsch {
22d94ba093SMatthew G. Knepley u[0] = 2.0;
233ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
24d94ba093SMatthew G. Knepley }
25d94ba093SMatthew G. Knepley
f0_quad_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])26d71ae5a4SJacob Faibussowitsch static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
27d71ae5a4SJacob Faibussowitsch {
28d94ba093SMatthew G. Knepley PetscInt d;
29d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += 2.0;
30d94ba093SMatthew G. Knepley }
31d94ba093SMatthew G. Knepley
f0_quad_p(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])32d71ae5a4SJacob Faibussowitsch static void f0_quad_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
33d71ae5a4SJacob Faibussowitsch {
34d94ba093SMatthew G. Knepley f0[0] = u[uOff[1]] - 2.0;
35d94ba093SMatthew G. Knepley }
36d94ba093SMatthew G. Knepley
f1_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f1[])37d71ae5a4SJacob Faibussowitsch static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
38d71ae5a4SJacob Faibussowitsch {
39d94ba093SMatthew G. Knepley PetscInt d;
40d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d];
41d94ba093SMatthew G. Knepley }
42d94ba093SMatthew G. Knepley
g3_uu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g3[])43d71ae5a4SJacob Faibussowitsch static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
44d71ae5a4SJacob Faibussowitsch {
45d94ba093SMatthew G. Knepley PetscInt d;
46d94ba093SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
47d94ba093SMatthew G. Knepley }
48d94ba093SMatthew G. Knepley
g0_pp(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g0[])49d71ae5a4SJacob Faibussowitsch static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
50d71ae5a4SJacob Faibussowitsch {
51d94ba093SMatthew G. Knepley g0[0] = 1.0;
52d94ba093SMatthew G. Knepley }
53d94ba093SMatthew G. Knepley
DivideDomain(DM dm,AppCtx * user)54d71ae5a4SJacob Faibussowitsch static PetscErrorCode DivideDomain(DM dm, AppCtx *user)
55d71ae5a4SJacob Faibussowitsch {
56d94ba093SMatthew G. Knepley DMLabel top, bottom;
57d94ba093SMatthew G. Knepley PetscReal low[3], high[3], midy;
58d94ba093SMatthew G. Knepley PetscInt cStart, cEnd, c;
59d94ba093SMatthew G. Knepley
60d94ba093SMatthew G. Knepley PetscFunctionBeginUser;
619566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "top"));
629566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "bottom"));
639566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "top", &top));
649566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "bottom", &bottom));
658fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm));
669566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, low, high));
67d94ba093SMatthew G. Knepley midy = 0.5 * (high[1] - low[1]);
689566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
69d94ba093SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {
70d94ba093SMatthew G. Knepley PetscReal centroid[3];
71d94ba093SMatthew G. Knepley
729566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
739566063dSJacob Faibussowitsch if (centroid[1] > midy) PetscCall(DMLabelSetValue(top, c, 1));
749566063dSJacob Faibussowitsch else PetscCall(DMLabelSetValue(bottom, c, 1));
75d94ba093SMatthew G. Knepley }
769566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(dm, top));
773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
78d94ba093SMatthew G. Knepley }
79d94ba093SMatthew G. Knepley
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)80d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
81d71ae5a4SJacob Faibussowitsch {
82d94ba093SMatthew G. Knepley PetscFunctionBeginUser;
839566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
849566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
859566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
869566063dSJacob Faibussowitsch PetscCall(DivideDomain(*dm, user));
879566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user));
889566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
90d94ba093SMatthew G. Knepley }
91d94ba093SMatthew G. Knepley
SetupPrimalProblem(DM dm,AppCtx * user)92d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
93d71ae5a4SJacob Faibussowitsch {
94d94ba093SMatthew G. Knepley PetscDS ds;
9576fb4698SMatthew G. Knepley PetscWeakForm wf;
9676fb4698SMatthew G. Knepley DMLabel label;
97d94ba093SMatthew G. Knepley const PetscInt id = 1;
98d94ba093SMatthew G. Knepley
99d94ba093SMatthew G. Knepley PetscFunctionBeginUser;
10007218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, 0, &label, NULL, &ds, NULL));
1019566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf));
1029566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
1039566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
1049566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
10507218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, 1, &label, NULL, &ds, NULL));
1069566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf));
1079566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
1089566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
1099566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, 0, f0_quad_p, 0, NULL));
1109566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL));
1119566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
1129566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quad_p, user));
1139566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
11457d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)quad_u, NULL, user, NULL));
1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
116d94ba093SMatthew G. Knepley }
117d94ba093SMatthew G. Knepley
SetupDiscretization(DM dm,const char name[],PetscErrorCode (* setup)(DM,AppCtx *),AppCtx * user)118d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
119d71ae5a4SJacob Faibussowitsch {
120d94ba093SMatthew G. Knepley DM cdm = dm;
121d94ba093SMatthew G. Knepley DMLabel top;
122d94ba093SMatthew G. Knepley PetscFE fe, feTop;
123d94ba093SMatthew G. Knepley PetscQuadrature q;
12430602db0SMatthew G. Knepley PetscInt dim;
12530602db0SMatthew G. Knepley PetscBool simplex;
126d94ba093SMatthew G. Knepley const char *nameTop = "pressure";
127d94ba093SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN];
128d94ba093SMatthew G. Knepley
129d94ba093SMatthew G. Knepley PetscFunctionBeginUser;
1309566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
1319566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex));
1329566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1339566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe));
1349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name));
1359566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1369566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &q));
1379566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
1389566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "top", &top));
1399566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop));
1409566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, name ? prefix : NULL, -1, &feTop));
1419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feTop, nameTop));
1429566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(feTop, q));
1439566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, top, (PetscObject)feTop));
1449566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feTop));
1459566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
1469566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user));
147d94ba093SMatthew G. Knepley while (cdm) {
1489566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
1499566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
150d94ba093SMatthew G. Knepley }
1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
152d94ba093SMatthew G. Knepley }
153d94ba093SMatthew G. Knepley
main(int argc,char ** argv)154d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
155d71ae5a4SJacob Faibussowitsch {
156d94ba093SMatthew G. Knepley DM dm; /* Problem specification */
157d94ba093SMatthew G. Knepley SNES snes; /* Nonlinear solver */
158d94ba093SMatthew G. Knepley Vec u; /* Solutions */
159d94ba093SMatthew G. Knepley AppCtx user; /* User-defined work context */
160d94ba093SMatthew G. Knepley
161327415f7SBarry Smith PetscFunctionBeginUser;
1629566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
163d94ba093SMatthew G. Knepley /* Primal system */
1649566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
1659566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1669566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm));
1679566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
1689566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u));
1699566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0));
1709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
1716493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
1729566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
1739566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u));
1749566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u));
1759566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u));
1769566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_view"));
177d94ba093SMatthew G. Knepley /* Cleanup */
1789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
1799566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
1809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
1819566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
182b122ec5aSJacob Faibussowitsch return 0;
183d94ba093SMatthew G. Knepley }
184d94ba093SMatthew G. Knepley
185d94ba093SMatthew G. Knepley /*TEST
186d94ba093SMatthew G. Knepley
187d94ba093SMatthew G. Knepley test:
188d94ba093SMatthew G. Knepley suffix: 2d_p1_0
189d94ba093SMatthew G. Knepley requires: triangle
190d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 2 -pressure_top_petscspace_degree 0 -dm_refine 0 -dmsnes_check
191d94ba093SMatthew G. Knepley
192d94ba093SMatthew G. Knepley test:
193d94ba093SMatthew G. Knepley suffix: 2d_p1_1
194d94ba093SMatthew G. Knepley requires: triangle
195d94ba093SMatthew G. Knepley args: -potential_petscspace_degree 1 -pressure_top_petscspace_degree 0 -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate
196d94ba093SMatthew G. Knepley
197d94ba093SMatthew G. Knepley TEST*/
198