1 static char help[] = "'Good Cop' Helmholtz Problem in 2d and 3d with finite elements.\n\ 2 We solve the 'Good Cop' Helmholtz problem in a rectangular\n\ 3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4 This example supports automatic convergence estimation\n\ 5 and coarse space adaptivity.\n\n\n"; 6 7 /* 8 The model problem: 9 Solve "Good Cop" Helmholtz equation on the unit square: (0,1) x (0,1) 10 - \Delta u + u = f, 11 where \Delta = Laplace operator 12 Dirichlet b.c.'s on all sides 13 14 */ 15 16 #include <petscdmplex.h> 17 #include <petscsnes.h> 18 #include <petscds.h> 19 #include <petscconvest.h> 20 21 typedef struct { 22 PetscBool trig; /* Use trig function as exact solution */ 23 } AppCtx; 24 25 /*For Primal Problem*/ 26 static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 27 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 28 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 29 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 30 { 31 PetscInt d; 32 for (d = 0; d < dim; ++d) g0[0] = 1.0; 33 } 34 35 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 36 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 37 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 38 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 39 { 40 PetscInt d; 41 for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 42 } 43 44 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 45 { 46 PetscInt d; 47 *u = 0.0; 48 for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 49 return 0; 50 } 51 52 static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 53 { 54 PetscInt d; 55 *u = 1.0; 56 for (d = 0; d < dim; ++d) *u += (d+1)*PetscSqr(x[d]); 57 return 0; 58 } 59 60 static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 61 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 62 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 63 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 64 { 65 PetscInt d; 66 f0[0] += u[0]; 67 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]); 68 } 69 70 static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 71 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 72 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 73 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 74 { 75 PetscInt d; 76 switch (dim) { 77 case 1: 78 f0[0] = 1.0; 79 break; 80 case 2: 81 f0[0] = 5.0; 82 break; 83 case 3: 84 f0[0] = 11.0; 85 break; 86 default: 87 f0[0] = 5.0; 88 break; 89 } 90 f0[0] += u[0]; 91 for (d = 0; d < dim; ++d) f0[0] -= (d+1)*PetscSqr(x[d]); 92 } 93 94 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 98 { 99 PetscInt d; 100 for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 101 } 102 103 static PetscErrorCode ProcessOptions(DM dm, AppCtx *options) 104 { 105 MPI_Comm comm; 106 PetscInt dim; 107 PetscErrorCode ierr; 108 109 PetscFunctionBeginUser; 110 CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 111 CHKERRQ(DMGetDimension(dm, &dim)); 112 options->trig = PETSC_FALSE; 113 114 ierr = PetscOptionsBegin(comm, "", "Helmholtz Problem Options", "DMPLEX");CHKERRQ(ierr); 115 CHKERRQ(PetscOptionsBool("-exact_trig", "Use trigonometric exact solution (better for more complex finite elements)", "ex26.c", options->trig, &options->trig, NULL)); 116 ierr = PetscOptionsEnd();CHKERRQ(ierr); 117 118 PetscFunctionReturn(0); 119 } 120 121 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 122 { 123 PetscFunctionBeginUser; 124 CHKERRQ(DMCreate(comm, dm)); 125 CHKERRQ(DMSetType(*dm, DMPLEX)); 126 CHKERRQ(DMSetFromOptions(*dm)); 127 128 CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Mesh")); 129 CHKERRQ(DMSetApplicationContext(*dm, user)); 130 CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 131 132 PetscFunctionReturn(0); 133 } 134 135 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 136 { 137 PetscDS ds; 138 DMLabel label; 139 const PetscInt id = 1; 140 141 PetscFunctionBeginUser; 142 CHKERRQ(DMGetDS(dm, &ds)); 143 CHKERRQ(DMGetLabel(dm, "marker", &label)); 144 if (user->trig) { 145 CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 146 CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu)); 147 CHKERRQ(PetscDSSetExactSolution(ds, 0, trig_u, user)); 148 CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL)); 149 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Trig Exact Solution\n")); 150 } else { 151 CHKERRQ(PetscDSSetResidual(ds, 0, f0_quad_u, f1_u)); 152 CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu)); 153 CHKERRQ(PetscDSSetExactSolution(ds, 0, quad_u, user)); 154 CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL)); 155 } 156 PetscFunctionReturn(0); 157 } 158 159 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 160 { 161 DM cdm = dm; 162 PetscFE fe; 163 DMPolytopeType ct; 164 PetscBool simplex; 165 PetscInt dim, cStart; 166 char prefix[PETSC_MAX_PATH_LEN]; 167 168 PetscFunctionBeginUser; 169 CHKERRQ(DMGetDimension(dm, &dim)); 170 171 CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 172 CHKERRQ(DMPlexGetCellType(dm, cStart, &ct)); 173 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 174 /* Create finite element */ 175 CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 176 CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 177 CHKERRQ(PetscObjectSetName((PetscObject) fe, name)); 178 /* Set discretization and boundary conditions for each mesh */ 179 CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 180 CHKERRQ(DMCreateDS(dm)); 181 CHKERRQ((*setup)(dm, user)); 182 while (cdm) { 183 CHKERRQ(DMCopyDisc(dm,cdm)); 184 CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 185 } 186 CHKERRQ(PetscFEDestroy(&fe)); 187 PetscFunctionReturn(0); 188 } 189 190 int main(int argc, char **argv) 191 { 192 DM dm; /* Problem specification */ 193 PetscDS ds; 194 SNES snes; /* Nonlinear solver */ 195 Vec u; /* Solutions */ 196 AppCtx user; /* User-defined work context */ 197 198 CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 199 /* Primal system */ 200 CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 201 CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 202 CHKERRQ(ProcessOptions(dm, &user)); 203 CHKERRQ(SNESSetDM(snes, dm)); 204 CHKERRQ(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 205 CHKERRQ(DMCreateGlobalVector(dm, &u)); 206 CHKERRQ(VecSet(u, 0.0)); 207 CHKERRQ(PetscObjectSetName((PetscObject) u, "potential")); 208 CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 209 CHKERRQ(SNESSetFromOptions(snes)); 210 CHKERRQ(DMSNESCheckFromOptions(snes, u)); 211 212 /*Looking for field error*/ 213 PetscInt Nfields; 214 CHKERRQ(DMGetDS(dm, &ds)); 215 CHKERRQ(PetscDSGetNumFields(ds, &Nfields)); 216 CHKERRQ(SNESSolve(snes, NULL, u)); 217 CHKERRQ(SNESGetSolution(snes, &u)); 218 CHKERRQ(VecViewFromOptions(u, NULL, "-potential_view")); 219 220 /* Cleanup */ 221 CHKERRQ(VecDestroy(&u)); 222 CHKERRQ(SNESDestroy(&snes)); 223 CHKERRQ(DMDestroy(&dm)); 224 CHKERRQ(PetscFinalize()); 225 return 0; 226 } 227 228 /*TEST 229 test: 230 # L_2 convergence rate: 1.9 231 suffix: 2d_p1_conv 232 requires: triangle 233 args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu 234 test: 235 # L_2 convergence rate: 1.9 236 suffix: 2d_q1_conv 237 args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu 238 test: 239 # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5 240 suffix: 3d_p1_conv 241 requires: ctetgen 242 args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 243 test: 244 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2 245 suffix: 3d_q1_conv 246 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 247 test: 248 # L_2 convergence rate: 1.9 249 suffix: 2d_p1_trig_conv 250 requires: triangle 251 args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig 252 test: 253 # L_2 convergence rate: 1.9 254 suffix: 2d_q1_trig_conv 255 args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig 256 test: 257 # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5 258 suffix: 3d_p1_trig_conv 259 requires: ctetgen 260 args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu -exact_trig 261 test: 262 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2 263 suffix: 3d_q1_trig_conv 264 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 265 test: 266 suffix: 2d_p1_gmg_vcycle 267 requires: triangle 268 args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 269 -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \ 270 -mg_levels_ksp_max_it 1 \ 271 -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 272 test: 273 suffix: 2d_p1_gmg_fcycle 274 requires: triangle 275 args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 276 -ksp_type cg -ksp_rtol 1e-10 -pc_type mg -pc_mg_type full \ 277 -mg_levels_ksp_max_it 2 \ 278 -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 279 test: 280 suffix: 2d_p1_gmg_vcycle_trig 281 requires: triangle 282 args: -exact_trig -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 283 -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \ 284 -mg_levels_ksp_max_it 1 \ 285 -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 286 TEST*/ 287