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