138be53adSDaniel Finn static char help[] = "'Good Cop' Helmholtz Problem in 2d and 3d with finite elements.\n\ 238be53adSDaniel Finn We solve the 'Good Cop' Helmholtz problem in a rectangular\n\ 338be53adSDaniel Finn domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 438be53adSDaniel Finn This example supports automatic convergence estimation\n\ 538be53adSDaniel Finn and coarse space adaptivity.\n\n\n"; 638be53adSDaniel Finn 738be53adSDaniel Finn /* 838be53adSDaniel Finn The model problem: 938be53adSDaniel Finn Solve "Good Cop" Helmholtz equation on the unit square: (0,1) x (0,1) 1038be53adSDaniel Finn - \Delta u + u = f, 1138be53adSDaniel Finn where \Delta = Laplace operator 1238be53adSDaniel Finn Dirichlet b.c.'s on all sides 1338be53adSDaniel Finn 1438be53adSDaniel Finn */ 1538be53adSDaniel Finn 1638be53adSDaniel Finn #include <petscdmplex.h> 1738be53adSDaniel Finn #include <petscsnes.h> 1838be53adSDaniel Finn #include <petscds.h> 1938be53adSDaniel Finn #include <petscconvest.h> 2038be53adSDaniel Finn 2138be53adSDaniel Finn typedef struct { 2238be53adSDaniel Finn PetscBool trig; /* Use trig function as exact solution */ 2338be53adSDaniel Finn } AppCtx; 2438be53adSDaniel Finn 2538be53adSDaniel Finn /*For Primal Problem*/ 2638be53adSDaniel Finn static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2738be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2838be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2938be53adSDaniel Finn PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 3038be53adSDaniel Finn { 3138be53adSDaniel Finn PetscInt d; 3238be53adSDaniel Finn for (d = 0; d < dim; ++d) g0[0] = 1.0; 3338be53adSDaniel Finn } 3438be53adSDaniel Finn 3538be53adSDaniel Finn static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 3638be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 3738be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 3838be53adSDaniel Finn PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 3938be53adSDaniel Finn { 4038be53adSDaniel Finn PetscInt d; 4138be53adSDaniel Finn for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 4238be53adSDaniel Finn } 4338be53adSDaniel Finn 4438be53adSDaniel Finn static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 4538be53adSDaniel Finn { 4638be53adSDaniel Finn PetscInt d; 4738be53adSDaniel Finn *u = 0.0; 4838be53adSDaniel Finn for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 4938be53adSDaniel Finn return 0; 5038be53adSDaniel Finn } 5138be53adSDaniel Finn 5238be53adSDaniel Finn static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 5338be53adSDaniel Finn { 5438be53adSDaniel Finn PetscInt d; 5538be53adSDaniel Finn *u = 1.0; 5638be53adSDaniel Finn for (d = 0; d < dim; ++d) *u += (d+1)*PetscSqr(x[d]); 5738be53adSDaniel Finn return 0; 5838be53adSDaniel Finn } 5938be53adSDaniel Finn 6038be53adSDaniel Finn static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 6138be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 6238be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 6338be53adSDaniel Finn PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 6438be53adSDaniel Finn { 6538be53adSDaniel Finn PetscInt d; 6638be53adSDaniel Finn f0[0] += u[0]; 6738be53adSDaniel Finn for (d = 0; d < dim; ++d) f0[0] -= 4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]) + PetscSinReal(2.0*PETSC_PI*x[d]); 6838be53adSDaniel Finn } 6938be53adSDaniel Finn 7038be53adSDaniel Finn static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 7138be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 7238be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 7338be53adSDaniel Finn PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 7438be53adSDaniel Finn { 7538be53adSDaniel Finn PetscInt d; 7638be53adSDaniel Finn switch (dim) { 7738be53adSDaniel Finn case 1: 7838be53adSDaniel Finn f0[0] = 1.0; 7938be53adSDaniel Finn break; 8038be53adSDaniel Finn case 2: 8138be53adSDaniel Finn f0[0] = 5.0; 8238be53adSDaniel Finn break; 8338be53adSDaniel Finn case 3: 8438be53adSDaniel Finn f0[0] = 11.0; 8538be53adSDaniel Finn break; 8638be53adSDaniel Finn default: 8738be53adSDaniel Finn f0[0] = 5.0; 8838be53adSDaniel Finn break; 8938be53adSDaniel Finn } 9038be53adSDaniel Finn f0[0] += u[0]; 9138be53adSDaniel Finn for (d = 0; d < dim; ++d) f0[0] -= (d+1)*PetscSqr(x[d]); 9238be53adSDaniel Finn } 9338be53adSDaniel Finn 9438be53adSDaniel Finn static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 9538be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 9638be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 9738be53adSDaniel Finn PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 9838be53adSDaniel Finn { 9938be53adSDaniel Finn PetscInt d; 10038be53adSDaniel Finn for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 10138be53adSDaniel Finn } 10238be53adSDaniel Finn 10338be53adSDaniel Finn static PetscErrorCode ProcessOptions(DM dm, AppCtx *options) 10438be53adSDaniel Finn { 10538be53adSDaniel Finn MPI_Comm comm; 10638be53adSDaniel Finn PetscInt dim; 10738be53adSDaniel Finn PetscErrorCode ierr; 10838be53adSDaniel Finn 10938be53adSDaniel Finn PetscFunctionBeginUser; 110*9566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 111*9566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 11238be53adSDaniel Finn options->trig = PETSC_FALSE; 11338be53adSDaniel Finn 114*9566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "Helmholtz Problem Options", "DMPLEX");PetscCall(ierr); 115*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact_trig", "Use trigonometric exact solution (better for more complex finite elements)", "ex26.c", options->trig, &options->trig, NULL)); 116*9566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 11738be53adSDaniel Finn 11838be53adSDaniel Finn PetscFunctionReturn(0); 11938be53adSDaniel Finn } 12038be53adSDaniel Finn 12138be53adSDaniel Finn static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 12238be53adSDaniel Finn { 12338be53adSDaniel Finn PetscFunctionBeginUser; 124*9566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 125*9566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 126*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 12738be53adSDaniel Finn 128*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Mesh")); 129*9566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 130*9566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 13138be53adSDaniel Finn 13238be53adSDaniel Finn PetscFunctionReturn(0); 13338be53adSDaniel Finn } 13438be53adSDaniel Finn 13538be53adSDaniel Finn static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 13638be53adSDaniel Finn { 13738be53adSDaniel Finn PetscDS ds; 13838be53adSDaniel Finn DMLabel label; 13938be53adSDaniel Finn const PetscInt id = 1; 14038be53adSDaniel Finn 14138be53adSDaniel Finn PetscFunctionBeginUser; 142*9566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 143*9566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 14438be53adSDaniel Finn if (user->trig) { 145*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 146*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu)); 147*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user)); 148*9566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL)); 149*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Trig Exact Solution\n")); 15038be53adSDaniel Finn } else { 151*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quad_u, f1_u)); 152*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu)); 153*9566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user)); 154*9566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL)); 15538be53adSDaniel Finn } 15638be53adSDaniel Finn PetscFunctionReturn(0); 15738be53adSDaniel Finn } 15838be53adSDaniel Finn 15938be53adSDaniel Finn static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 16038be53adSDaniel Finn { 16138be53adSDaniel Finn DM cdm = dm; 16238be53adSDaniel Finn PetscFE fe; 16338be53adSDaniel Finn DMPolytopeType ct; 16438be53adSDaniel Finn PetscBool simplex; 16538be53adSDaniel Finn PetscInt dim, cStart; 16638be53adSDaniel Finn char prefix[PETSC_MAX_PATH_LEN]; 16738be53adSDaniel Finn 16838be53adSDaniel Finn PetscFunctionBeginUser; 169*9566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 17038be53adSDaniel Finn 171*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 172*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 17338be53adSDaniel Finn simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 17438be53adSDaniel Finn /* Create finite element */ 175*9566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 176*9566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 177*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 17838be53adSDaniel Finn /* Set discretization and boundary conditions for each mesh */ 179*9566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 180*9566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 181*9566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 18238be53adSDaniel Finn while (cdm) { 183*9566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 184*9566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 18538be53adSDaniel Finn } 186*9566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 18738be53adSDaniel Finn PetscFunctionReturn(0); 18838be53adSDaniel Finn } 18938be53adSDaniel Finn 19038be53adSDaniel Finn int main(int argc, char **argv) 19138be53adSDaniel Finn { 19238be53adSDaniel Finn DM dm; /* Problem specification */ 19338be53adSDaniel Finn PetscDS ds; 19438be53adSDaniel Finn SNES snes; /* Nonlinear solver */ 19538be53adSDaniel Finn Vec u; /* Solutions */ 19638be53adSDaniel Finn AppCtx user; /* User-defined work context */ 19738be53adSDaniel Finn 198*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 19938be53adSDaniel Finn /* Primal system */ 200*9566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 201*9566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 202*9566063dSJacob Faibussowitsch PetscCall(ProcessOptions(dm, &user)); 203*9566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 204*9566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 205*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 206*9566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 207*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 208*9566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 209*9566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 210*9566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 21138be53adSDaniel Finn 21238be53adSDaniel Finn /*Looking for field error*/ 21338be53adSDaniel Finn PetscInt Nfields; 214*9566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 215*9566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nfields)); 216*9566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 217*9566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 218*9566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 21938be53adSDaniel Finn 22038be53adSDaniel Finn /* Cleanup */ 221*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 222*9566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 223*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 224*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 225b122ec5aSJacob Faibussowitsch return 0; 22638be53adSDaniel Finn } 22738be53adSDaniel Finn 22838be53adSDaniel Finn /*TEST 22938be53adSDaniel Finn test: 23038be53adSDaniel Finn # L_2 convergence rate: 1.9 23138be53adSDaniel Finn suffix: 2d_p1_conv 23238be53adSDaniel Finn requires: triangle 23338be53adSDaniel Finn args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu 23438be53adSDaniel Finn test: 23538be53adSDaniel Finn # L_2 convergence rate: 1.9 23638be53adSDaniel Finn suffix: 2d_q1_conv 23738be53adSDaniel Finn args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu 23838be53adSDaniel Finn test: 23938be53adSDaniel Finn # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5 24038be53adSDaniel Finn suffix: 3d_p1_conv 24138be53adSDaniel Finn requires: ctetgen 24238be53adSDaniel Finn args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 24338be53adSDaniel Finn test: 24438be53adSDaniel Finn # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2 24538be53adSDaniel Finn suffix: 3d_q1_conv 24638be53adSDaniel Finn args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 24738be53adSDaniel Finn test: 24838be53adSDaniel Finn # L_2 convergence rate: 1.9 24938be53adSDaniel Finn suffix: 2d_p1_trig_conv 25038be53adSDaniel Finn requires: triangle 25138be53adSDaniel Finn args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig 25238be53adSDaniel Finn test: 25338be53adSDaniel Finn # L_2 convergence rate: 1.9 25438be53adSDaniel Finn suffix: 2d_q1_trig_conv 25538be53adSDaniel Finn args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig 25638be53adSDaniel Finn test: 25738be53adSDaniel Finn # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5 25838be53adSDaniel Finn suffix: 3d_p1_trig_conv 25938be53adSDaniel Finn requires: ctetgen 26038be53adSDaniel Finn args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu -exact_trig 26138be53adSDaniel Finn test: 26238be53adSDaniel Finn # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2 26338be53adSDaniel Finn suffix: 3d_q1_trig_conv 26438be53adSDaniel Finn args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu -exact_trig 26538be53adSDaniel Finn test: 26638be53adSDaniel Finn suffix: 2d_p1_gmg_vcycle 26738be53adSDaniel Finn requires: triangle 26838be53adSDaniel Finn args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 26938be53adSDaniel Finn -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \ 27038be53adSDaniel Finn -mg_levels_ksp_max_it 1 \ 27138be53adSDaniel Finn -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 27238be53adSDaniel Finn test: 27338be53adSDaniel Finn suffix: 2d_p1_gmg_fcycle 27438be53adSDaniel Finn requires: triangle 27538be53adSDaniel Finn args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 27638be53adSDaniel Finn -ksp_type cg -ksp_rtol 1e-10 -pc_type mg -pc_mg_type full \ 27738be53adSDaniel Finn -mg_levels_ksp_max_it 2 \ 27838be53adSDaniel Finn -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 27938be53adSDaniel Finn test: 28038be53adSDaniel Finn suffix: 2d_p1_gmg_vcycle_trig 28138be53adSDaniel Finn requires: triangle 28238be53adSDaniel Finn args: -exact_trig -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 28338be53adSDaniel Finn -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \ 28438be53adSDaniel Finn -mg_levels_ksp_max_it 1 \ 28538be53adSDaniel Finn -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 28638be53adSDaniel Finn TEST*/ 287