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