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); 77efeeb221SMatthew G. Knepley ierr = PetscOptionsBool("-cr", "Use compatible relaxarion", "ex20.c", options->cr, &options->cr, NULL);CHKERRQ(ierr); 78*1e1ea65dSPierre 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 PetscErrorCode ierr; 85efeeb221SMatthew G. Knepley 86efeeb221SMatthew G. Knepley PetscFunctionBeginUser; 8730602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 8830602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 89efeeb221SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 9030602db0SMatthew G. Knepley ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 91efeeb221SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 92efeeb221SMatthew G. Knepley PetscFunctionReturn(0); 93efeeb221SMatthew G. Knepley } 94efeeb221SMatthew G. Knepley 95efeeb221SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 96efeeb221SMatthew G. Knepley { 9730602db0SMatthew G. Knepley PetscDS ds; 9845480ffeSMatthew G. Knepley DMLabel label; 99efeeb221SMatthew G. Knepley const PetscInt id = 1; 100efeeb221SMatthew G. Knepley PetscErrorCode ierr; 101efeeb221SMatthew G. Knepley 102efeeb221SMatthew G. Knepley PetscFunctionBeginUser; 10330602db0SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 10445480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 10530602db0SMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr); 10630602db0SMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 10730602db0SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, trig_u, user);CHKERRQ(ierr); 10845480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL);CHKERRQ(ierr); 109efeeb221SMatthew G. Knepley PetscFunctionReturn(0); 110efeeb221SMatthew G. Knepley } 111efeeb221SMatthew G. Knepley 112efeeb221SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 113efeeb221SMatthew G. Knepley { 114efeeb221SMatthew G. Knepley DM cdm = dm; 115efeeb221SMatthew G. Knepley PetscFE fe; 116efeeb221SMatthew G. Knepley DMPolytopeType ct; 117efeeb221SMatthew G. Knepley PetscBool simplex; 118efeeb221SMatthew G. Knepley PetscInt dim, cStart; 119efeeb221SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 120efeeb221SMatthew G. Knepley PetscErrorCode ierr; 121efeeb221SMatthew G. Knepley 122efeeb221SMatthew G. Knepley PetscFunctionBeginUser; 123efeeb221SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 124efeeb221SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 125efeeb221SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 126efeeb221SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 127efeeb221SMatthew G. Knepley 128efeeb221SMatthew G. Knepley ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 129efeeb221SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 130efeeb221SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 131efeeb221SMatthew G. Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 132efeeb221SMatthew G. Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 133efeeb221SMatthew G. Knepley ierr = (*setup)(dm, user);CHKERRQ(ierr); 134efeeb221SMatthew G. Knepley while (cdm) { 135efeeb221SMatthew G. Knepley ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 136efeeb221SMatthew G. Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 137efeeb221SMatthew G. Knepley } 138efeeb221SMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 139efeeb221SMatthew G. Knepley PetscFunctionReturn(0); 140efeeb221SMatthew G. Knepley } 141efeeb221SMatthew G. Knepley 142efeeb221SMatthew G. Knepley /* 143efeeb221SMatthew G. Knepley How to do CR in PETSc: 144efeeb221SMatthew G. Knepley 145efeeb221SMatthew G. Knepley Loop over PCMG levels, coarse to fine: 146efeeb221SMatthew G. Knepley Run smoother for 5 iterates 147efeeb221SMatthew G. Knepley At each iterate, solve Inj u_f = u_c with LSQR to 1e-15 148efeeb221SMatthew G. Knepley Suppose that e_k = c^k e_0, which means log e_k = log e_0 + k log c 149efeeb221SMatthew G. Knepley Fit log of error to look at log c, the slope 150efeeb221SMatthew G. Knepley Check R^2 for linearity (1 - square residual / variance) 151efeeb221SMatthew G. Knepley Solve exactly 152efeeb221SMatthew G. Knepley Prolong to next level 153efeeb221SMatthew G. Knepley */ 154efeeb221SMatthew G. Knepley 155efeeb221SMatthew G. Knepley int main(int argc, char **argv) 156efeeb221SMatthew G. Knepley { 157efeeb221SMatthew G. Knepley DM dm; /* Problem specification */ 158efeeb221SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 159efeeb221SMatthew G. Knepley Vec u; /* Solutions */ 160efeeb221SMatthew G. Knepley AppCtx user; /* User-defined work context */ 161efeeb221SMatthew G. Knepley PetscErrorCode ierr; 162efeeb221SMatthew G. Knepley 163efeeb221SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 164efeeb221SMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 165efeeb221SMatthew G. Knepley /* Primal system */ 166efeeb221SMatthew G. Knepley ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 167efeeb221SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 168efeeb221SMatthew G. Knepley ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 169efeeb221SMatthew G. Knepley ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr); 170efeeb221SMatthew G. Knepley ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 171efeeb221SMatthew G. Knepley ierr = VecSet(u, 0.0);CHKERRQ(ierr); 172efeeb221SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr); 173efeeb221SMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 174efeeb221SMatthew G. Knepley ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 175efeeb221SMatthew G. Knepley ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 176efeeb221SMatthew G. Knepley ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 177efeeb221SMatthew G. Knepley ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr); 178efeeb221SMatthew G. Knepley /* Cleanup */ 179efeeb221SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 180efeeb221SMatthew G. Knepley ierr = SNESDestroy(&snes);CHKERRQ(ierr); 181efeeb221SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 182efeeb221SMatthew G. Knepley ierr = PetscFinalize(); 183efeeb221SMatthew G. Knepley return ierr; 184efeeb221SMatthew G. Knepley } 185efeeb221SMatthew G. Knepley 186efeeb221SMatthew G. Knepley /*TEST 187efeeb221SMatthew G. Knepley 188efeeb221SMatthew G. Knepley test: 189efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_rate 190efeeb221SMatthew G. Knepley requires: triangle 191efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 192efeeb221SMatthew G. Knepley -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \ 193efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 194efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 195efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 196efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 197efeeb221SMatthew G. Knepley -mg_levels_pc_type jacobi 198efeeb221SMatthew G. Knepley 199efeeb221SMatthew G. Knepley test: 200efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_cr 201efeeb221SMatthew G. Knepley requires: triangle 202efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 203efeeb221SMatthew G. Knepley -ksp_rtol 5e-10 -pc_type mg -pc_mg_adapt_cr \ 204efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned \ 205efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 206efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 207efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 208efeeb221SMatthew G. Knepley -mg_levels_cr_ksp_max_it 5 -mg_levels_cr_ksp_converged_rate -mg_levels_cr_ksp_converged_rate_type error 209efeeb221SMatthew G. Knepley 210efeeb221SMatthew G. Knepley test: 211efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_fcycle_rate 212efeeb221SMatthew G. Knepley requires: triangle 213efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 214efeeb221SMatthew G. Knepley -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg -pc_mg_type full \ 215efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 216efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 217efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 218efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 219efeeb221SMatthew G. Knepley -mg_levels_pc_type jacobi 220efeeb221SMatthew G. Knepley test: 221efeeb221SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_adapt_rate 222efeeb221SMatthew G. Knepley requires: triangle bamg 223efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 224efeeb221SMatthew G. Knepley -ksp_rtol 5e-10 -ksp_converged_rate -pc_type mg \ 225efeeb221SMatthew G. Knepley -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \ 226efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 227efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 228efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 229efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 230efeeb221SMatthew G. Knepley -mg_levels_pc_type jacobi 231efeeb221SMatthew G. Knepley test: 232efeeb221SMatthew G. Knepley suffix: 2d_p1_scalable_rate 233efeeb221SMatthew G. Knepley requires: triangle 234efeeb221SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 3 \ 235efeeb221SMatthew G. Knepley -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -ksp_converged_rate \ 236efeeb221SMatthew G. Knepley -pc_type gamg \ 237efeeb221SMatthew G. Knepley -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 238efeeb221SMatthew G. Knepley -pc_gamg_coarse_eq_limit 1000 \ 239efeeb221SMatthew G. Knepley -pc_gamg_square_graph 1 \ 240efeeb221SMatthew G. Knepley -pc_gamg_threshold 0.05 \ 241efeeb221SMatthew G. Knepley -pc_gamg_threshold_scale .0 \ 242efeeb221SMatthew G. Knepley -mg_levels_ksp_type chebyshev -mg_levels_ksp_norm_type preconditioned -mg_levels_ksp_converged_rate \ 243efeeb221SMatthew G. Knepley -mg_levels_ksp_max_it 5 \ 244efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 245efeeb221SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 246efeeb221SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 247efeeb221SMatthew G. Knepley -mg_levels_pc_type jacobi \ 248efeeb221SMatthew G. Knepley -matptap_via scalable 249efeeb221SMatthew G. Knepley 250efeeb221SMatthew G. Knepley TEST*/ 251