xref: /petsc/src/snes/tutorials/ex23.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLabel(dm, "top"));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLabel(dm, "bottom"));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "top", &top));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "bottom", &bottom));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetBoundingBox(dm, low, high));
81d94ba093SMatthew G. Knepley   midy = 0.5*(high[1] - low[1]);
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
83d94ba093SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
84d94ba093SMatthew G. Knepley     PetscReal centroid[3];
85d94ba093SMatthew G. Knepley 
865f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
875f80ce2aSJacob Faibussowitsch     if (centroid[1] > midy) CHKERRQ(DMLabelSetValue(top, c, 1));
885f80ce2aSJacob Faibussowitsch     else                    CHKERRQ(DMLabelSetValue(bottom, c, 1));
89d94ba093SMatthew G. Knepley   }
905f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(DivideDomain(*dm, user));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(*dm, user));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetRegionNumDS(dm, 0, &label, NULL, &ds));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetWeakForm(ds, &wf));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 0, quad_u, user));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetRegionNumDS(dm, 1, &label, NULL, &ds));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetWeakForm(ds, &wf));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, 0, f0_quad_p, 0, NULL));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 0, quad_u, user));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 1, quad_p, user));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsSimplex(dm, &simplex));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, name));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetQuadrature(fe, &q));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "top", &top));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &feTop));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) feTop, nameTop));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetQuadrature(feTop, q));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 1, top, (PetscObject) feTop));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&feTop));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ((*setup)(dm, user));
161d94ba093SMatthew G. Knepley   while (cdm) {
1625f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm,cdm));
1635f80ce2aSJacob Faibussowitsch     CHKERRQ(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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL,help));
176d94ba093SMatthew G. Knepley   /* Primal system */
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes, dm));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &u));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u, 0.0));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) u, "solution"));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESCheckFromOptions(snes, u));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes, NULL, u));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes, &u));
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-sol_view"));
190d94ba093SMatthew G. Knepley   /* Cleanup */
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
194*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
195*b122ec5aSJacob 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