xref: /petsc/src/snes/tutorials/ex20.c (revision 8ebe3e4e9e00d86ece2e9fcd0cc84910b0ad437c)
1 static char help[] = "Poisson Problem with finite elements.\n\
2 This example supports automatic convergence estimation for multilevel solvers\n\
3 and solver adaptivity.\n\n\n";
4 
5 #include <petscdmplex.h>
6 #include <petscsnes.h>
7 #include <petscds.h>
8 #include <petscconvest.h>
9 
10 /* Next steps:
11 
12 - Show lowest eigenmodes using SLEPc code from my ex6
13 
14 - Run CR example from Brannick's slides that looks like semicoarsening
15   - Show lowest modes
16   - Show CR convergence rate
17   - Show CR solution to show non-convergence
18   - Refine coarse grid around non-converged dofs
19     - Maybe use Barry's "more then Z% above the average" monitor to label bad dofs
20     - Mark coarse cells that contain bad dofs
21     - Run SBR on coarse grid
22 
23 - Run Helmholtz example from Gander's writeup
24 
25 - Run Low Mach example?
26 
27 - Run subduction example?
28 */
29 
30 typedef struct {
31   PetscBool cr;  /* Use compatible relaxation */
32 } AppCtx;
33 
34 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
35 {
36   PetscInt d;
37   u[0] = 0.0;
38   for (d = 0; d < dim; ++d) u[0] += PetscSinReal(2.0*PETSC_PI*x[d]);
39   return 0;
40 }
41 
42 static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
43                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
44                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
45                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
46 {
47   PetscInt d;
48   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
49 }
50 
51 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
52                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
53                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
54                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
55 {
56   PetscInt d;
57   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
58 }
59 
60 static void g3_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
64 {
65   PetscInt d;
66   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
67 }
68 
69 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
70 {
71   PetscErrorCode ierr;
72 
73   PetscFunctionBeginUser;
74   options->cr = PETSC_FALSE;
75 
76   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
77   ierr = PetscOptionsBool("-cr", "Use compatible relaxarion", "ex20.c", options->cr, &options->cr, NULL);CHKERRQ(ierr);
78   ierr = PetscOptionsEnd();CHKERRQ(ierr);
79   PetscFunctionReturn(0);
80 }
81 
82 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
83 {
84   PetscErrorCode ierr;
85 
86   PetscFunctionBeginUser;
87   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
88   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
89   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
90   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
91   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
96 {
97   PetscDS        ds;
98   DMLabel        label;
99   const PetscInt id = 1;
100   PetscErrorCode ierr;
101 
102   PetscFunctionBeginUser;
103   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
104   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
105   ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr);
106   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
107   ierr = PetscDSSetExactSolution(ds, 0, trig_u, user);CHKERRQ(ierr);
108   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL);CHKERRQ(ierr);
109   PetscFunctionReturn(0);
110 }
111 
112 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
113 {
114   DM             cdm = dm;
115   PetscFE        fe;
116   DMPolytopeType ct;
117   PetscBool      simplex;
118   PetscInt       dim, cStart;
119   char           prefix[PETSC_MAX_PATH_LEN];
120   PetscErrorCode ierr;
121 
122   PetscFunctionBeginUser;
123   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
124   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
125   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
126   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
127 
128   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
129   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
130   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
131   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
132   ierr = DMCreateDS(dm);CHKERRQ(ierr);
133   ierr = (*setup)(dm, user);CHKERRQ(ierr);
134   while (cdm) {
135     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
136     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
137   }
138   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 /*
143   How to do CR in PETSc:
144 
145 Loop over PCMG levels, coarse to fine:
146   Run smoother for 5 iterates
147     At each iterate, solve Inj u_f = u_c with LSQR to 1e-15
148     Suppose that e_k = c^k e_0, which means log e_k = log e_0 + k log c
149       Fit log of error to look at log c, the slope
150       Check R^2 for linearity (1 - square residual / variance)
151   Solve exactly
152   Prolong to next level
153 */
154 
155 int main(int argc, char **argv)
156 {
157   DM             dm;   /* Problem specification */
158   SNES           snes; /* Nonlinear solver */
159   Vec            u;    /* Solutions */
160   AppCtx         user; /* User-defined work context */
161   PetscErrorCode ierr;
162 
163   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
164   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
165   /* Primal system */
166   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
167   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
168   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
169   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
170   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
171   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
172   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
173   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
174   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
175   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
176   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
177   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
178   /* Cleanup */
179   ierr = VecDestroy(&u);CHKERRQ(ierr);
180   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
181   ierr = DMDestroy(&dm);CHKERRQ(ierr);
182   ierr = PetscFinalize();
183   return ierr;
184 }
185 
186 /*TEST
187 
188   test:
189     suffix: 2d_p1_gmg_vcycle_rate
190     requires: triangle
191     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
192           -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \
193             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
194             -mg_levels_esteig_ksp_type cg \
195             -mg_levels_esteig_ksp_max_it 10 \
196             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
197             -mg_levels_pc_type jacobi
198 
199   test:
200     suffix: 2d_p1_gmg_vcycle_cr
201     TODO: broken
202     # cannot MatShift() a MATNORMAL until this MatType inherits from MATSHELL, cf. https://gitlab.com/petsc/petsc/-/issues/972
203     requires: triangle
204     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
205           -ksp_rtol 5e-10 -pc_type mg  -pc_mg_adapt_cr \
206             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned \
207             -mg_levels_esteig_ksp_type cg \
208             -mg_levels_esteig_ksp_max_it 10 \
209             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
210             -mg_levels_cr_ksp_max_it 5 -mg_levels_cr_ksp_converged_rate -mg_levels_cr_ksp_converged_rate_type error
211 
212   test:
213     suffix: 2d_p1_gmg_fcycle_rate
214     requires: triangle
215     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
216           -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg -pc_mg_type full \
217             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
218             -mg_levels_esteig_ksp_type cg \
219             -mg_levels_esteig_ksp_max_it 10 \
220             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
221             -mg_levels_pc_type jacobi
222   test:
223     suffix: 2d_p1_gmg_vcycle_adapt_rate
224     requires: triangle bamg
225     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
226           -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \
227             -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \
228             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
229             -mg_levels_esteig_ksp_type cg \
230             -mg_levels_esteig_ksp_max_it 10 \
231             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
232             -mg_levels_pc_type jacobi
233   test:
234     suffix: 2d_p1_scalable_rate
235     requires: triangle
236     args: -potential_petscspace_degree 1 -dm_refine 3 \
237       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -ksp_converged_rate \
238       -pc_type gamg \
239         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
240         -pc_gamg_coarse_eq_limit 1000 \
241         -pc_gamg_square_graph 1 \
242         -pc_gamg_threshold 0.05 \
243         -pc_gamg_threshold_scale .0 \
244       -mg_levels_ksp_type chebyshev -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
245         -mg_levels_ksp_max_it 5 \
246         -mg_levels_esteig_ksp_type cg \
247         -mg_levels_esteig_ksp_max_it 10 \
248         -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
249         -mg_levels_pc_type jacobi \
250       -matptap_via scalable
251 
252 TEST*/
253