xref: /petsc/src/snes/tutorials/ex23.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
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;
769566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "top"));
779566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "bottom"));
789566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "top", &top));
799566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "bottom", &bottom));
80*8fb5bd83SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
819566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dm, low, high));
82d94ba093SMatthew G. Knepley   midy = 0.5*(high[1] - low[1]);
839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
84d94ba093SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
85d94ba093SMatthew G. Knepley     PetscReal centroid[3];
86d94ba093SMatthew G. Knepley 
879566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
889566063dSJacob Faibussowitsch     if (centroid[1] > midy) PetscCall(DMLabelSetValue(top, c, 1));
899566063dSJacob Faibussowitsch     else                    PetscCall(DMLabelSetValue(bottom, c, 1));
90d94ba093SMatthew G. Knepley   }
919566063dSJacob Faibussowitsch   PetscCall(DMPlexLabelComplete(dm, top));
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   PetscFunctionBeginUser;
989566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
999566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1009566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1019566063dSJacob Faibussowitsch   PetscCall(DivideDomain(*dm, user));
1029566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
1039566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
104d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
105d94ba093SMatthew G. Knepley }
106d94ba093SMatthew G. Knepley 
107d94ba093SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
108d94ba093SMatthew G. Knepley {
109d94ba093SMatthew G. Knepley   PetscDS        ds;
11076fb4698SMatthew G. Knepley   PetscWeakForm  wf;
11176fb4698SMatthew G. Knepley   DMLabel        label;
112d94ba093SMatthew G. Knepley   const PetscInt id = 1;
113d94ba093SMatthew G. Knepley 
114d94ba093SMatthew G. Knepley   PetscFunctionBeginUser;
1159566063dSJacob Faibussowitsch   PetscCall(DMGetRegionNumDS(dm, 0, &label, NULL, &ds));
1169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
1179566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
1189566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
1199566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
1209566063dSJacob Faibussowitsch   PetscCall(DMGetRegionNumDS(dm, 1, &label, NULL, &ds));
1219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
1229566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
1239566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
1249566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, 0, f0_quad_p, 0, NULL));
1259566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL));
1269566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
1279566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 1, quad_p, user));
1289566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
1299566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL));
130d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
131d94ba093SMatthew G. Knepley }
132d94ba093SMatthew G. Knepley 
133d94ba093SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
134d94ba093SMatthew G. Knepley {
135d94ba093SMatthew G. Knepley   DM              cdm = dm;
136d94ba093SMatthew G. Knepley   DMLabel         top;
137d94ba093SMatthew G. Knepley   PetscFE         fe, feTop;
138d94ba093SMatthew G. Knepley   PetscQuadrature q;
13930602db0SMatthew G. Knepley   PetscInt        dim;
14030602db0SMatthew G. Knepley   PetscBool       simplex;
141d94ba093SMatthew G. Knepley   const char     *nameTop = "pressure";
142d94ba093SMatthew G. Knepley   char            prefix[PETSC_MAX_PATH_LEN];
143d94ba093SMatthew G. Knepley 
144d94ba093SMatthew G. Knepley   PetscFunctionBeginUser;
1459566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1469566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1479566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1489566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe));
1499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, name));
1509566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &q));
1529566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1539566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "top", &top));
1549566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &feTop));
1569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) feTop, nameTop));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(feTop, q));
1589566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 1, top, (PetscObject) feTop));
1599566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&feTop));
1609566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1619566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
162d94ba093SMatthew G. Knepley   while (cdm) {
1639566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm,cdm));
1649566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
165d94ba093SMatthew G. Knepley   }
166d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
167d94ba093SMatthew G. Knepley }
168d94ba093SMatthew G. Knepley 
169d94ba093SMatthew G. Knepley int main(int argc, char **argv)
170d94ba093SMatthew G. Knepley {
171d94ba093SMatthew G. Knepley   DM             dm;   /* Problem specification */
172d94ba093SMatthew G. Knepley   SNES           snes; /* Nonlinear solver */
173d94ba093SMatthew G. Knepley   Vec            u;    /* Solutions */
174d94ba093SMatthew G. Knepley   AppCtx         user; /* User-defined work context */
175d94ba093SMatthew G. Knepley 
1769566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
177d94ba093SMatthew G. Knepley   /* Primal system */
1789566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
1799566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1809566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
1819566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
1829566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
1839566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) u, "solution"));
1859566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
1869566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
1879566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
1889566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
1899566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
1909566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_view"));
191d94ba093SMatthew G. Knepley   /* Cleanup */
1929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1939566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1949566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1959566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
196b122ec5aSJacob Faibussowitsch   return 0;
197d94ba093SMatthew G. Knepley }
198d94ba093SMatthew G. Knepley 
199d94ba093SMatthew G. Knepley /*TEST
200d94ba093SMatthew G. Knepley 
201d94ba093SMatthew G. Knepley   test:
202d94ba093SMatthew G. Knepley     suffix: 2d_p1_0
203d94ba093SMatthew G. Knepley     requires: triangle
204d94ba093SMatthew G. Knepley     args: -potential_petscspace_degree 2 -pressure_top_petscspace_degree 0 -dm_refine 0 -dmsnes_check
205d94ba093SMatthew G. Knepley 
206d94ba093SMatthew G. Knepley   test:
207d94ba093SMatthew G. Knepley     suffix: 2d_p1_1
208d94ba093SMatthew G. Knepley     requires: triangle
209d94ba093SMatthew G. Knepley     args: -potential_petscspace_degree 1 -pressure_top_petscspace_degree 0 -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate
210d94ba093SMatthew G. Knepley 
211d94ba093SMatthew G. Knepley TEST*/
212