xref: /petsc/src/snes/tutorials/ex20.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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
19467446fbSPierre Jolivet     - Maybe use Barry's "more than 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 
trig_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)34*2a8381b2SBarry Smith static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
35d71ae5a4SJacob Faibussowitsch {
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]);
393ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
40efeeb221SMatthew G. Knepley }
41efeeb221SMatthew G. Knepley 
f0_trig_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[])42d71ae5a4SJacob Faibussowitsch static void f0_trig_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[])
43d71ae5a4SJacob Faibussowitsch {
44efeeb221SMatthew G. Knepley   PetscInt d;
45efeeb221SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
46efeeb221SMatthew G. Knepley }
47efeeb221SMatthew G. Knepley 
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[])48d71ae5a4SJacob Faibussowitsch 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[])
49d71ae5a4SJacob Faibussowitsch {
50efeeb221SMatthew G. Knepley   PetscInt d;
51efeeb221SMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
52efeeb221SMatthew G. Knepley }
53efeeb221SMatthew G. Knepley 
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[])54d71ae5a4SJacob Faibussowitsch 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[])
55d71ae5a4SJacob Faibussowitsch {
56efeeb221SMatthew G. Knepley   PetscInt d;
57efeeb221SMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
58efeeb221SMatthew G. Knepley }
59efeeb221SMatthew G. Knepley 
ProcessOptions(MPI_Comm comm,AppCtx * options)60d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
61d71ae5a4SJacob Faibussowitsch {
62efeeb221SMatthew G. Knepley   PetscFunctionBeginUser;
63efeeb221SMatthew G. Knepley   options->cr = PETSC_FALSE;
64d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cr", "Use compatible relaxarion", "ex20.c", options->cr, &options->cr, NULL));
66d0609cedSBarry Smith   PetscOptionsEnd();
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68efeeb221SMatthew G. Knepley }
69efeeb221SMatthew G. Knepley 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)70d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
71d71ae5a4SJacob Faibussowitsch {
72efeeb221SMatthew G. Knepley   PetscFunctionBeginUser;
739566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
749566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
759566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
769566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
779566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79efeeb221SMatthew G. Knepley }
80efeeb221SMatthew G. Knepley 
SetupPrimalProblem(DM dm,AppCtx * user)81d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
82d71ae5a4SJacob Faibussowitsch {
8330602db0SMatthew G. Knepley   PetscDS        ds;
8445480ffeSMatthew G. Knepley   DMLabel        label;
85efeeb221SMatthew G. Knepley   const PetscInt id = 1;
86efeeb221SMatthew G. Knepley 
87efeeb221SMatthew G. Knepley   PetscFunctionBeginUser;
889566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
899566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
909566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
919566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
929566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
9357d50842SBarry Smith   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)trig_u, NULL, user, NULL));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95efeeb221SMatthew G. Knepley }
96efeeb221SMatthew G. Knepley 
SetupDiscretization(DM dm,const char name[],PetscErrorCode (* setup)(DM,AppCtx *),AppCtx * user)97d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
98d71ae5a4SJacob Faibussowitsch {
99efeeb221SMatthew G. Knepley   DM             cdm = dm;
100efeeb221SMatthew G. Knepley   PetscFE        fe;
101efeeb221SMatthew G. Knepley   DMPolytopeType ct;
102efeeb221SMatthew G. Knepley   PetscBool      simplex;
103efeeb221SMatthew G. Knepley   PetscInt       dim, cStart;
104efeeb221SMatthew G. Knepley   char           prefix[PETSC_MAX_PATH_LEN];
105efeeb221SMatthew G. Knepley 
106efeeb221SMatthew G. Knepley   PetscFunctionBeginUser;
1079566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1089566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1099566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
110efeeb221SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
111efeeb221SMatthew G. Knepley 
1129566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1139566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe));
1149566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
1159566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1169566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1179566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
118efeeb221SMatthew G. Knepley   while (cdm) {
1199566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
1209566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
121efeeb221SMatthew G. Knepley   }
1229566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124efeeb221SMatthew G. Knepley }
125efeeb221SMatthew G. Knepley 
126efeeb221SMatthew G. Knepley /*
127efeeb221SMatthew G. Knepley   How to do CR in PETSc:
128efeeb221SMatthew G. Knepley 
129efeeb221SMatthew G. Knepley Loop over PCMG levels, coarse to fine:
130efeeb221SMatthew G. Knepley   Run smoother for 5 iterates
131efeeb221SMatthew G. Knepley     At each iterate, solve Inj u_f = u_c with LSQR to 1e-15
132efeeb221SMatthew G. Knepley     Suppose that e_k = c^k e_0, which means log e_k = log e_0 + k log c
133efeeb221SMatthew G. Knepley       Fit log of error to look at log c, the slope
134efeeb221SMatthew G. Knepley       Check R^2 for linearity (1 - square residual / variance)
135efeeb221SMatthew G. Knepley   Solve exactly
136efeeb221SMatthew G. Knepley   Prolong to next level
137efeeb221SMatthew G. Knepley */
138efeeb221SMatthew G. Knepley 
main(int argc,char ** argv)139d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
140d71ae5a4SJacob Faibussowitsch {
141efeeb221SMatthew G. Knepley   DM     dm;   /* Problem specification */
142efeeb221SMatthew G. Knepley   SNES   snes; /* Nonlinear solver */
143efeeb221SMatthew G. Knepley   Vec    u;    /* Solutions */
144efeeb221SMatthew G. Knepley   AppCtx user; /* User-defined work context */
145efeeb221SMatthew G. Knepley 
146327415f7SBarry Smith   PetscFunctionBeginUser;
1479566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1489566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
149efeeb221SMatthew G. Knepley   /* Primal system */
1509566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
1519566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1529566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
1539566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
1549566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
1559566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
1569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
1576493148fSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
1589566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
1599566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
1609566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
1619566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
162efeeb221SMatthew G. Knepley   /* Cleanup */
1639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1649566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1659566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1669566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
167b122ec5aSJacob Faibussowitsch   return 0;
168efeeb221SMatthew G. Knepley }
169efeeb221SMatthew G. Knepley 
170efeeb221SMatthew G. Knepley /*TEST
171efeeb221SMatthew G. Knepley 
172efeeb221SMatthew G. Knepley   test:
173efeeb221SMatthew G. Knepley     suffix: 2d_p1_gmg_vcycle_rate
174efeeb221SMatthew G. Knepley     requires: triangle
175efeeb221SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
176efeeb221SMatthew G. Knepley           -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \
177efeeb221SMatthew G. Knepley             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
178efeeb221SMatthew G. Knepley             -mg_levels_esteig_ksp_type cg \
179efeeb221SMatthew G. Knepley             -mg_levels_esteig_ksp_max_it 10 \
180efeeb221SMatthew G. Knepley             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
181efeeb221SMatthew G. Knepley             -mg_levels_pc_type jacobi
182efeeb221SMatthew G. Knepley 
183efeeb221SMatthew G. Knepley   test:
184efeeb221SMatthew G. Knepley     suffix: 2d_p1_gmg_vcycle_cr
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 -pc_type mg  -pc_mg_adapt_cr \
188efeeb221SMatthew G. Knepley             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned \
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_cr_ksp_max_it 5 -mg_levels_cr_ksp_converged_rate -mg_levels_cr_ksp_converged_rate_type error
193efeeb221SMatthew G. Knepley 
194efeeb221SMatthew G. Knepley   test:
195efeeb221SMatthew G. Knepley     suffix: 2d_p1_gmg_fcycle_rate
196efeeb221SMatthew G. Knepley     requires: triangle
197efeeb221SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
198efeeb221SMatthew G. Knepley           -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg -pc_mg_type full \
199efeeb221SMatthew G. Knepley             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
200efeeb221SMatthew G. Knepley             -mg_levels_esteig_ksp_type cg \
201efeeb221SMatthew G. Knepley             -mg_levels_esteig_ksp_max_it 10 \
202efeeb221SMatthew G. Knepley             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
203efeeb221SMatthew G. Knepley             -mg_levels_pc_type jacobi
204efeeb221SMatthew G. Knepley   test:
205efeeb221SMatthew G. Knepley     suffix: 2d_p1_gmg_vcycle_adapt_rate
2062b3cbbdaSStefano Zampini     requires: triangle
207908b9b43SStefano Zampini     args: -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
208efeeb221SMatthew G. Knepley           -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \
2092b3cbbdaSStefano Zampini             -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \
210efeeb221SMatthew G. Knepley             -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
211efeeb221SMatthew G. Knepley             -mg_levels_esteig_ksp_type cg \
212efeeb221SMatthew G. Knepley             -mg_levels_esteig_ksp_max_it 10 \
213efeeb221SMatthew G. Knepley             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
214efeeb221SMatthew G. Knepley             -mg_levels_pc_type jacobi
215efeeb221SMatthew G. Knepley   test:
216efeeb221SMatthew G. Knepley     suffix: 2d_p1_scalable_rate
217efeeb221SMatthew G. Knepley     requires: triangle
218efeeb221SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 3 \
219efeeb221SMatthew G. Knepley       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -ksp_converged_rate \
22073f7197eSJed Brown       -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_esteig_ksp_type cg \
221efeeb221SMatthew G. Knepley         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
222efeeb221SMatthew G. Knepley         -pc_gamg_coarse_eq_limit 1000 \
223efeeb221SMatthew G. Knepley         -pc_gamg_threshold 0.05 \
224efeeb221SMatthew G. Knepley         -pc_gamg_threshold_scale .0 \
225efeeb221SMatthew G. Knepley         -mg_levels_ksp_type chebyshev -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \
226efeeb221SMatthew G. Knepley         -mg_levels_ksp_max_it 5                                                \
227efeeb221SMatthew G. Knepley       -matptap_via scalable
228efeeb221SMatthew G. Knepley 
229efeeb221SMatthew G. Knepley TEST*/
230