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 PetscFunctionBeginUser; 72efeeb221SMatthew G. Knepley options->cr = PETSC_FALSE; 73d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cr", "Use compatible relaxarion", "ex20.c", options->cr, &options->cr, NULL)); 75d0609cedSBarry Smith PetscOptionsEnd(); 76efeeb221SMatthew G. Knepley PetscFunctionReturn(0); 77efeeb221SMatthew G. Knepley } 78efeeb221SMatthew G. Knepley 79efeeb221SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 80efeeb221SMatthew G. Knepley { 81efeeb221SMatthew G. Knepley PetscFunctionBeginUser; 829566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 839566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 849566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 859566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 869566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 87efeeb221SMatthew G. Knepley PetscFunctionReturn(0); 88efeeb221SMatthew G. Knepley } 89efeeb221SMatthew G. Knepley 90efeeb221SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 91efeeb221SMatthew G. Knepley { 9230602db0SMatthew G. Knepley PetscDS ds; 9345480ffeSMatthew G. Knepley DMLabel label; 94efeeb221SMatthew G. Knepley const PetscInt id = 1; 95efeeb221SMatthew G. Knepley 96efeeb221SMatthew G. Knepley PetscFunctionBeginUser; 979566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 989566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 999566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 1009566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1019566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user)); 1029566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL)); 103efeeb221SMatthew G. Knepley PetscFunctionReturn(0); 104efeeb221SMatthew G. Knepley } 105efeeb221SMatthew G. Knepley 106efeeb221SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 107efeeb221SMatthew G. Knepley { 108efeeb221SMatthew G. Knepley DM cdm = dm; 109efeeb221SMatthew G. Knepley PetscFE fe; 110efeeb221SMatthew G. Knepley DMPolytopeType ct; 111efeeb221SMatthew G. Knepley PetscBool simplex; 112efeeb221SMatthew G. Knepley PetscInt dim, cStart; 113efeeb221SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 114efeeb221SMatthew G. Knepley 115efeeb221SMatthew G. Knepley PetscFunctionBeginUser; 1169566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1179566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1189566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 119efeeb221SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 120efeeb221SMatthew G. Knepley 1219566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 1229566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 1239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 1249566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 1259566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1269566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 127efeeb221SMatthew G. Knepley while (cdm) { 1289566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 1299566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 130efeeb221SMatthew G. Knepley } 1319566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 132efeeb221SMatthew G. Knepley PetscFunctionReturn(0); 133efeeb221SMatthew G. Knepley } 134efeeb221SMatthew G. Knepley 135efeeb221SMatthew G. Knepley /* 136efeeb221SMatthew G. Knepley How to do CR in PETSc: 137efeeb221SMatthew G. Knepley 138efeeb221SMatthew G. Knepley Loop over PCMG levels, coarse to fine: 139efeeb221SMatthew G. Knepley Run smoother for 5 iterates 140efeeb221SMatthew G. Knepley At each iterate, solve Inj u_f = u_c with LSQR to 1e-15 141efeeb221SMatthew G. Knepley Suppose that e_k = c^k e_0, which means log e_k = log e_0 + k log c 142efeeb221SMatthew G. Knepley Fit log of error to look at log c, the slope 143efeeb221SMatthew G. Knepley Check R^2 for linearity (1 - square residual / variance) 144efeeb221SMatthew G. Knepley Solve exactly 145efeeb221SMatthew G. Knepley Prolong to next level 146efeeb221SMatthew G. Knepley */ 147efeeb221SMatthew G. Knepley 148efeeb221SMatthew G. Knepley int main(int argc, char **argv) 149efeeb221SMatthew G. Knepley { 150efeeb221SMatthew G. Knepley DM dm; /* Problem specification */ 151efeeb221SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 152efeeb221SMatthew G. Knepley Vec u; /* Solutions */ 153efeeb221SMatthew G. Knepley AppCtx user; /* User-defined work context */ 154efeeb221SMatthew G. Knepley 1559566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 1569566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 157efeeb221SMatthew G. Knepley /* Primal system */ 1589566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 1599566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1609566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 1619566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 1629566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 1639566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 1649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 1659566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 1669566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 1679566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 1689566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 1699566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 170efeeb221SMatthew G. Knepley /* Cleanup */ 1719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1729566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1739566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 175b122ec5aSJacob Faibussowitsch return 0; 176efeeb221SMatthew G. Knepley } 177efeeb221SMatthew G. Knepley 178efeeb221SMatthew G. Knepley /*TEST 179efeeb221SMatthew G. Knepley 180efeeb221SMatthew G. Knepley test: 181efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_rate 182efeeb221SMatthew G. Knepley requires: triangle 183efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 184efeeb221SMatthew G. Knepley -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \ 185efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 186efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 187efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 188efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 189efeeb221SMatthew G. Knepley -mg_levels_pc_type jacobi 190efeeb221SMatthew G. Knepley 191efeeb221SMatthew G. Knepley test: 192efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_cr 193dcc276ddSPierre Jolivet TODO: broken 194dcc276ddSPierre Jolivet # cannot MatShift() a MATNORMAL until this MatType inherits from MATSHELL, cf. https://gitlab.com/petsc/petsc/-/issues/972 195efeeb221SMatthew G. Knepley requires: triangle 196efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 197efeeb221SMatthew G. Knepley -ksp_rtol 5e-10 -pc_type mg -pc_mg_adapt_cr \ 198efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned \ 199efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 200efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 201efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 202efeeb221SMatthew G. Knepley -mg_levels_cr_ksp_max_it 5 -mg_levels_cr_ksp_converged_rate -mg_levels_cr_ksp_converged_rate_type error 203efeeb221SMatthew G. Knepley 204efeeb221SMatthew G. Knepley test: 205efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_fcycle_rate 206efeeb221SMatthew G. Knepley requires: triangle 207efeeb221SMatthew G. Knepley args: -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 -pc_mg_type full \ 209efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 210efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 211efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 212efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 213efeeb221SMatthew G. Knepley -mg_levels_pc_type jacobi 214efeeb221SMatthew G. Knepley test: 215efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_adapt_rate 216*2b3cbbdaSStefano Zampini requires: triangle 217efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 218efeeb221SMatthew G. Knepley -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \ 219*2b3cbbdaSStefano Zampini -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \ 220efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 221efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 222efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 223efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 224efeeb221SMatthew G. Knepley -mg_levels_pc_type jacobi 225efeeb221SMatthew G. Knepley test: 226efeeb221SMatthew G. Knepley suffix: 2d_p1_scalable_rate 227efeeb221SMatthew G. Knepley requires: triangle 228efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 3 \ 229efeeb221SMatthew G. Knepley -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -ksp_converged_rate \ 23073f7197eSJed Brown -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_esteig_ksp_type cg \ 231efeeb221SMatthew G. Knepley -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 232efeeb221SMatthew G. Knepley -pc_gamg_coarse_eq_limit 1000 \ 233efeeb221SMatthew G. Knepley -pc_gamg_square_graph 1 \ 234efeeb221SMatthew G. Knepley -pc_gamg_threshold 0.05 \ 235efeeb221SMatthew G. Knepley -pc_gamg_threshold_scale .0 \ 236efeeb221SMatthew G. Knepley -mg_levels_ksp_type chebyshev -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 237efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 \ 238efeeb221SMatthew G. Knepley -matptap_via scalable 239efeeb221SMatthew G. Knepley 240efeeb221SMatthew G. Knepley TEST*/ 241