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