xref: /petsc/src/snes/tutorials/ex26.c (revision 7a101e5e7ba9859de4c800924a501d6ea3cd325c)
1 static char help[] = "'Good Cop' Helmholtz Problem in 2d and 3d with finite elements.\n\
2 We solve the 'Good Cop' Helmholtz problem in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4 This example supports automatic convergence estimation\n\
5 and coarse space adaptivity.\n\n\n";
6 
7 /*
8    The model problem:
9       Solve "Good Cop" Helmholtz equation on the unit square: (0,1) x (0,1)
10           - \Delta u + u = f,
11            where \Delta = Laplace operator
12       Dirichlet b.c.'s on all sides
13 
14 */
15 
16 #include <petscdmplex.h>
17 #include <petscsnes.h>
18 #include <petscds.h>
19 #include <petscconvest.h>
20 
21 typedef struct {
22   PetscBool trig; /* Use trig function as exact solution */
23 } AppCtx;
24 
25 /*For Primal Problem*/
26 static void g0_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
30 {
31   PetscInt d;
32   for (d = 0; d < dim; ++d) g0[0] = 1.0;
33 }
34 
35 static void g3_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
39 {
40   PetscInt d;
41   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
42 }
43 
44 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
45 {
46   PetscInt d;
47   *u = 0.0;
48   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
49   return 0;
50 }
51 
52 static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
53 {
54   PetscInt d;
55   *u = 1.0;
56   for (d = 0; d < dim; ++d) *u += (d+1)*PetscSqr(x[d]);
57   return 0;
58 }
59 
60 static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
61                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
62                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
63                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
64 {
65   PetscInt d;
66   f0[0] += u[0];
67   for (d = 0; d < dim; ++d) f0[0] -= 4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]) + PetscSinReal(2.0*PETSC_PI*x[d]);
68 }
69 
70 static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
71                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
72                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
73                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
74 {
75     PetscInt d;
76     switch (dim) {
77         case 1:
78             f0[0] = 1.0;
79             break;
80         case 2:
81             f0[0] = 5.0;
82             break;
83         case 3:
84             f0[0] = 11.0;
85             break;
86         default:
87             f0[0] = 5.0;
88             break;
89     }
90     f0[0] += u[0];
91     for (d = 0; d < dim; ++d) f0[0] -= (d+1)*PetscSqr(x[d]);
92 }
93 
94 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
98 {
99   PetscInt d;
100   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
101 }
102 
103 static PetscErrorCode ProcessOptions(DM dm, AppCtx *options)
104 {
105   MPI_Comm       comm;
106   PetscInt       dim;
107 
108   PetscFunctionBeginUser;
109   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
110   PetscCall(DMGetDimension(dm, &dim));
111   options->trig = PETSC_FALSE;
112 
113   PetscOptionsBegin(comm, "", "Helmholtz Problem Options", "DMPLEX");
114   PetscCall(PetscOptionsBool("-exact_trig", "Use trigonometric exact solution (better for more complex finite elements)", "ex26.c", options->trig, &options->trig, NULL));
115   PetscOptionsEnd();
116 
117   PetscFunctionReturn(0);
118 }
119 
120 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
121 {
122   PetscFunctionBeginUser;
123   PetscCall(DMCreate(comm, dm));
124   PetscCall(DMSetType(*dm, DMPLEX));
125   PetscCall(DMSetFromOptions(*dm));
126 
127   PetscCall(PetscObjectSetName((PetscObject) *dm, "Mesh"));
128   PetscCall(DMSetApplicationContext(*dm, user));
129   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
130 
131   PetscFunctionReturn(0);
132 }
133 
134 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
135 {
136   PetscDS        ds;
137   DMLabel        label;
138   const PetscInt id = 1;
139 
140   PetscFunctionBeginUser;
141   PetscCall(DMGetDS(dm, &ds));
142   PetscCall(DMGetLabel(dm, "marker", &label));
143   if (user->trig) {
144     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
145     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu));
146     PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
147     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL));
148     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Trig Exact Solution\n"));
149   } else {
150     PetscCall(PetscDSSetResidual(ds, 0, f0_quad_u, f1_u));
151     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu));
152     PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
153     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL));
154   }
155   PetscFunctionReturn(0);
156 }
157 
158 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
159 {
160   DM             cdm = dm;
161   PetscFE        fe;
162   DMPolytopeType ct;
163   PetscBool      simplex;
164   PetscInt       dim, cStart;
165   char           prefix[PETSC_MAX_PATH_LEN];
166 
167   PetscFunctionBeginUser;
168   PetscCall(DMGetDimension(dm, &dim));
169 
170   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
171   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
172   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
173   /* Create finite element */
174   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
175   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
176   PetscCall(PetscObjectSetName((PetscObject) fe, name));
177   /* Set discretization and boundary conditions for each mesh */
178   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
179   PetscCall(DMCreateDS(dm));
180   PetscCall((*setup)(dm, user));
181   while (cdm) {
182     PetscCall(DMCopyDisc(dm,cdm));
183     PetscCall(DMGetCoarseDM(cdm, &cdm));
184   }
185   PetscCall(PetscFEDestroy(&fe));
186   PetscFunctionReturn(0);
187 }
188 
189 int main(int argc, char **argv)
190 {
191     DM             dm;   /* Problem specification */
192     PetscDS        ds;
193     SNES           snes; /* Nonlinear solver */
194     Vec            u;    /* Solutions */
195     AppCtx         user; /* User-defined work context */
196 
197     PetscFunctionBeginUser;
198   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
199     /* Primal system */
200     PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
201     PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
202     PetscCall(ProcessOptions(dm, &user));
203     PetscCall(SNESSetDM(snes, dm));
204     PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
205     PetscCall(DMCreateGlobalVector(dm, &u));
206     PetscCall(VecSet(u, 0.0));
207     PetscCall(PetscObjectSetName((PetscObject) u, "potential"));
208     PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
209     PetscCall(SNESSetFromOptions(snes));
210     PetscCall(DMSNESCheckFromOptions(snes, u));
211 
212     /*Looking for field error*/
213     PetscInt Nfields;
214     PetscCall(DMGetDS(dm, &ds));
215     PetscCall(PetscDSGetNumFields(ds, &Nfields));
216     PetscCall(SNESSolve(snes, NULL, u));
217     PetscCall(SNESGetSolution(snes, &u));
218     PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
219 
220     /* Cleanup */
221     PetscCall(VecDestroy(&u));
222     PetscCall(SNESDestroy(&snes));
223     PetscCall(DMDestroy(&dm));
224     PetscCall(PetscFinalize());
225     return 0;
226 }
227 
228 /*TEST
229 test:
230   # L_2 convergence rate: 1.9
231   suffix: 2d_p1_conv
232   requires: triangle
233   args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu
234 test:
235   # L_2 convergence rate: 1.9
236   suffix: 2d_q1_conv
237   args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu
238 test:
239   # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5
240   suffix: 3d_p1_conv
241   requires: ctetgen
242   args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu
243 test:
244   # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2
245   suffix: 3d_q1_conv
246   args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu
247 test:
248   # L_2 convergence rate: 1.9
249   suffix: 2d_p1_trig_conv
250   requires: triangle
251   args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig
252 test:
253   # L_2 convergence rate: 1.9
254   suffix: 2d_q1_trig_conv
255   args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig
256 test:
257   # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5
258   suffix: 3d_p1_trig_conv
259   requires: ctetgen
260   args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu -exact_trig
261 test:
262   # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2
263   suffix: 3d_q1_trig_conv
264   args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu -exact_trig
265 test:
266   suffix: 2d_p1_gmg_vcycle
267   requires: triangle
268   args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
269     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \
270     -mg_levels_ksp_max_it 1 \
271     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
272 test:
273   suffix: 2d_p1_gmg_fcycle
274   requires: triangle
275   args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
276     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg -pc_mg_type full \
277     -mg_levels_ksp_max_it 2 \
278     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
279 test:
280   suffix: 2d_p1_gmg_vcycle_trig
281   requires: triangle
282   args: -exact_trig -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
283     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \
284     -mg_levels_ksp_max_it 1 \
285     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
286 test:
287   suffix: 2d_q3_cgns
288   requires: cgns
289   args: -dm_plex_simplex 0 -dm_plex_dim 2 -dm_plex_box_faces 3,3 -snes_view_solution cgns:sol.cgns -potential_petscspace_degree 3 -dm_coord_petscspace_degree 3
290 TEST*/
291