xref: /petsc/src/snes/tutorials/ex23.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Poisson Problem with a split domain.\n\
2 We solve the Poisson problem in two halves of a domain.\n\
3 In one half, we include an additional field.\n\n\n";
4 
5 #include <petscdmplex.h>
6 #include <petscsnes.h>
7 #include <petscds.h>
8 #include <petscconvest.h>
9 
10 typedef struct {
11   PetscInt dummy;
12 } AppCtx;
13 
14 static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
15   u[0] = x[0] * x[0] + x[1] * x[1];
16   return 0;
17 }
18 
19 static PetscErrorCode quad_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
20   u[0] = 2.0;
21   return 0;
22 }
23 
24 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[]) {
25   PetscInt d;
26   for (d = 0; d < dim; ++d) f0[0] += 2.0;
27 }
28 
29 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[]) {
30   f0[0] = u[uOff[1]] - 2.0;
31 }
32 
33 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[]) {
34   PetscInt d;
35   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
36 }
37 
38 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[]) {
39   PetscInt d;
40   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
41 }
42 
43 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[]) {
44   g0[0] = 1.0;
45 }
46 
47 static PetscErrorCode DivideDomain(DM dm, AppCtx *user) {
48   DMLabel   top, bottom;
49   PetscReal low[3], high[3], midy;
50   PetscInt  cStart, cEnd, c;
51 
52   PetscFunctionBeginUser;
53   PetscCall(DMCreateLabel(dm, "top"));
54   PetscCall(DMCreateLabel(dm, "bottom"));
55   PetscCall(DMGetLabel(dm, "top", &top));
56   PetscCall(DMGetLabel(dm, "bottom", &bottom));
57   PetscCall(DMGetCoordinatesLocalSetUp(dm));
58   PetscCall(DMGetBoundingBox(dm, low, high));
59   midy = 0.5 * (high[1] - low[1]);
60   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
61   for (c = cStart; c < cEnd; ++c) {
62     PetscReal centroid[3];
63 
64     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
65     if (centroid[1] > midy) PetscCall(DMLabelSetValue(top, c, 1));
66     else PetscCall(DMLabelSetValue(bottom, c, 1));
67   }
68   PetscCall(DMPlexLabelComplete(dm, top));
69   PetscFunctionReturn(0);
70 }
71 
72 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
73   PetscFunctionBeginUser;
74   PetscCall(DMCreate(comm, dm));
75   PetscCall(DMSetType(*dm, DMPLEX));
76   PetscCall(DMSetFromOptions(*dm));
77   PetscCall(DivideDomain(*dm, user));
78   PetscCall(DMSetApplicationContext(*dm, user));
79   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
80   PetscFunctionReturn(0);
81 }
82 
83 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) {
84   PetscDS        ds;
85   PetscWeakForm  wf;
86   DMLabel        label;
87   const PetscInt id = 1;
88 
89   PetscFunctionBeginUser;
90   PetscCall(DMGetRegionNumDS(dm, 0, &label, NULL, &ds));
91   PetscCall(PetscDSGetWeakForm(ds, &wf));
92   PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
93   PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
94   PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
95   PetscCall(DMGetRegionNumDS(dm, 1, &label, NULL, &ds));
96   PetscCall(PetscDSGetWeakForm(ds, &wf));
97   PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u));
98   PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
99   PetscCall(PetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, 0, f0_quad_p, 0, NULL));
100   PetscCall(PetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL));
101   PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
102   PetscCall(PetscDSSetExactSolution(ds, 1, quad_p, user));
103   PetscCall(DMGetLabel(dm, "marker", &label));
104   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))quad_u, NULL, user, NULL));
105   PetscFunctionReturn(0);
106 }
107 
108 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) {
109   DM              cdm = dm;
110   DMLabel         top;
111   PetscFE         fe, feTop;
112   PetscQuadrature q;
113   PetscInt        dim;
114   PetscBool       simplex;
115   const char     *nameTop = "pressure";
116   char            prefix[PETSC_MAX_PATH_LEN];
117 
118   PetscFunctionBeginUser;
119   PetscCall(DMGetDimension(dm, &dim));
120   PetscCall(DMPlexIsSimplex(dm, &simplex));
121   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
122   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe));
123   PetscCall(PetscObjectSetName((PetscObject)fe, name));
124   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
125   PetscCall(PetscFEGetQuadrature(fe, &q));
126   PetscCall(PetscFEDestroy(&fe));
127   PetscCall(DMGetLabel(dm, "top", &top));
128   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop));
129   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, name ? prefix : NULL, -1, &feTop));
130   PetscCall(PetscObjectSetName((PetscObject)feTop, nameTop));
131   PetscCall(PetscFESetQuadrature(feTop, q));
132   PetscCall(DMSetField(dm, 1, top, (PetscObject)feTop));
133   PetscCall(PetscFEDestroy(&feTop));
134   PetscCall(DMCreateDS(dm));
135   PetscCall((*setup)(dm, user));
136   while (cdm) {
137     PetscCall(DMCopyDisc(dm, cdm));
138     PetscCall(DMGetCoarseDM(cdm, &cdm));
139   }
140   PetscFunctionReturn(0);
141 }
142 
143 int main(int argc, char **argv) {
144   DM     dm;   /* Problem specification */
145   SNES   snes; /* Nonlinear solver */
146   Vec    u;    /* Solutions */
147   AppCtx user; /* User-defined work context */
148 
149   PetscFunctionBeginUser;
150   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
151   /* Primal system */
152   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
153   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
154   PetscCall(SNESSetDM(snes, dm));
155   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
156   PetscCall(DMCreateGlobalVector(dm, &u));
157   PetscCall(VecSet(u, 0.0));
158   PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
159   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
160   PetscCall(SNESSetFromOptions(snes));
161   PetscCall(DMSNESCheckFromOptions(snes, u));
162   PetscCall(SNESSolve(snes, NULL, u));
163   PetscCall(SNESGetSolution(snes, &u));
164   PetscCall(VecViewFromOptions(u, NULL, "-sol_view"));
165   /* Cleanup */
166   PetscCall(VecDestroy(&u));
167   PetscCall(SNESDestroy(&snes));
168   PetscCall(DMDestroy(&dm));
169   PetscCall(PetscFinalize());
170   return 0;
171 }
172 
173 /*TEST
174 
175   test:
176     suffix: 2d_p1_0
177     requires: triangle
178     args: -potential_petscspace_degree 2 -pressure_top_petscspace_degree 0 -dm_refine 0 -dmsnes_check
179 
180   test:
181     suffix: 2d_p1_1
182     requires: triangle
183     args: -potential_petscspace_degree 1 -pressure_top_petscspace_degree 0 -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate
184 
185 TEST*/
186