1 static char help[] = "Benchmark Poisson Problem in 2d and 3d with finite elements.\n\ 2 We solve the Poisson problem in a rectangular domain\n\ 3 using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4 5 #include <petscdmplex.h> 6 #include <petscsnes.h> 7 #include <petscds.h> 8 #include <petscconvest.h> 9 10 typedef struct { 11 PetscInt nit; /* Number of benchmark iterations */ 12 PetscBool strong; /* Do not integrate the Laplacian by parts */ 13 } AppCtx; 14 15 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 16 { 17 PetscInt d; 18 *u = 0.0; 19 for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 20 return 0; 21 } 22 23 static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 24 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 25 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 26 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 27 { 28 PetscInt d; 29 for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 30 } 31 32 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 33 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 34 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 35 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 36 { 37 PetscInt d; 38 for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 39 } 40 41 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 42 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 43 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 44 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 45 { 46 PetscInt d; 47 for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 48 } 49 50 static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 51 { 52 *u = PetscSqr(x[0]) + PetscSqr(x[1]); 53 return 0; 54 } 55 56 static void f0_strong_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 57 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 58 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 59 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 60 { 61 PetscInt d; 62 for (d = 0; d < dim; ++d) f0[0] -= u_x[dim + d*dim+d]; 63 f0[0] += 4.0; 64 } 65 66 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 67 { 68 PetscErrorCode ierr; 69 70 PetscFunctionBeginUser; 71 options->nit = 10; 72 options->strong = PETSC_FALSE; 73 ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");PetscCall(ierr); 74 PetscCall(PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL)); 75 PetscCall(PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL)); 76 ierr = PetscOptionsEnd();PetscCall(ierr); 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 81 { 82 PetscFunctionBeginUser; 83 PetscCall(DMCreate(comm, dm)); 84 PetscCall(DMSetType(*dm, DMPLEX)); 85 PetscCall(DMSetFromOptions(*dm)); 86 PetscCall(DMSetApplicationContext(*dm, user)); 87 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 88 PetscFunctionReturn(0); 89 } 90 91 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 92 { 93 PetscDS ds; 94 DMLabel label; 95 const PetscInt id = 1; 96 97 PetscFunctionBeginUser; 98 PetscCall(DMGetDS(dm, &ds)); 99 PetscCall(DMGetLabel(dm, "marker", &label)); 100 if (user->strong) { 101 PetscCall(PetscDSSetResidual(ds, 0, f0_strong_u, NULL)); 102 PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, user)); 103 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quadratic_u, NULL, user, NULL)); 104 } else { 105 PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 106 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 107 PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user)); 108 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL)); 109 } 110 PetscFunctionReturn(0); 111 } 112 113 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 114 { 115 DM cdm = dm; 116 PetscFE fe; 117 DMPolytopeType ct; 118 PetscBool simplex; 119 PetscInt dim, cStart; 120 char prefix[PETSC_MAX_PATH_LEN]; 121 122 PetscFunctionBeginUser; 123 PetscCall(DMGetDimension(dm, &dim)); 124 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 125 PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 126 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; // false 127 /* Create finite element */ 128 PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 129 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 130 PetscCall(PetscObjectSetName((PetscObject) fe, name)); 131 /* Set discretization and boundary conditions for each mesh */ 132 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 133 PetscCall(DMCreateDS(dm)); 134 PetscCall((*setup)(dm, user)); 135 while (cdm) { 136 PetscCall(DMCopyDisc(dm,cdm)); 137 /* TODO: Check whether the boundary of coarse meshes is marked */ 138 PetscCall(DMGetCoarseDM(cdm, &cdm)); 139 } 140 PetscCall(PetscFEDestroy(&fe)); 141 PetscFunctionReturn(0); 142 } 143 144 int main(int argc, char **argv) 145 { 146 DM dm; /* Problem specification */ 147 SNES snes; /* Nonlinear solver */ 148 Vec u; /* Solutions */ 149 AppCtx user; /* User-defined work context */ 150 PetscLogDouble time; 151 Mat Amat; 152 153 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 154 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 155 /* system */ 156 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 157 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 158 PetscCall(SNESSetDM(snes, dm)); 159 PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 160 PetscCall(DMCreateGlobalVector(dm, &u)); 161 PetscCall(SNESSetFromOptions(snes)); 162 PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 163 PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 164 PetscCall(DMSNESCheckFromOptions(snes, u)); 165 PetscCall(PetscTime(&time)); 166 PetscCall(SNESSetUp(snes)); 167 PetscCall(SNESGetJacobian(snes, &Amat, NULL, NULL, NULL)); 168 PetscCall(MatSetOption(Amat,MAT_SPD,PETSC_TRUE)); 169 PetscCall(SNESSolve(snes, NULL, u)); 170 PetscCall(PetscTimeSubtract(&time)); 171 // PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First Solve time: %g\n",-time)); 172 /* Benchmark system */ 173 if (user.nit) { 174 Vec b; 175 PetscInt i; 176 #if defined(PETSC_USE_LOG) 177 PetscLogStage kspstage; 178 #endif 179 PetscCall(PetscLogStageRegister("Solve only", &kspstage)); 180 PetscCall(PetscLogStagePush(kspstage)); 181 PetscCall(SNESGetSolution(snes, &u)); 182 PetscCall(SNESGetFunction(snes, &b, NULL, NULL)); 183 for (i=0;i<user.nit;i++) { 184 PetscCall(VecZeroEntries(u)); 185 PetscCall(SNESSolve(snes, NULL, u)); 186 } 187 PetscCall(PetscLogStagePop()); 188 } 189 PetscCall(SNESGetSolution(snes, &u)); 190 PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 191 /* Cleanup */ 192 PetscCall(VecDestroy(&u)); 193 PetscCall(SNESDestroy(&snes)); 194 PetscCall(DMDestroy(&dm)); 195 PetscCall(PetscFinalize()); 196 return 0; 197 } 198 199 /*TEST 200 201 test: 202 suffix: strong 203 requires: triangle 204 args: -dm_plex_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check \ 205 -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong 206 207 test: 208 suffix: bench 209 nsize: 4 210 args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,1 -dm_refine 2 -dm_view -ksp_monitor \ 211 -benchmark_it 1 \ 212 -dm_plex_box_upper 2,2,1 \ 213 -dm_plex_box_lower 0,0,0 \ 214 -dm_plex_dim 3 \ 215 -ksp_converged_reason \ 216 -ksp_norm_type unpreconditioned \ 217 -ksp_rtol 1.e-6 \ 218 -ksp_type cg \ 219 -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.1 \ 220 -mg_levels_ksp_max_it 1 \ 221 -mg_levels_ksp_type chebyshev \ 222 -mg_levels_pc_type jacobi \ 223 -pc_gamg_coarse_eq_limit 200 \ 224 -pc_gamg_coarse_grid_layout_type compact \ 225 -pc_gamg_esteig_ksp_max_it 5 \ 226 -pc_gamg_process_eq_limit 200 \ 227 -pc_gamg_repartition false \ 228 -pc_gamg_reuse_interpolation true \ 229 -pc_gamg_square_graph 0 \ 230 -pc_gamg_threshold 0.001 -pc_gamg_threshold_scale .5\ 231 -pc_gamg_type agg \ 232 -pc_type gamg \ 233 -petscpartitioner_simple_node_grid 1,2,1 \ 234 -petscpartitioner_simple_process_grid 2,1,1 \ 235 -petscpartitioner_type simple \ 236 -potential_petscspace_degree 2 \ 237 -snes_lag_jacobian -2 \ 238 -snes_max_it 1 \ 239 -snes_rtol 1.e-8 \ 240 -snes_type ksponly \ 241 -use_gpu_aware_mpi true 242 243 test: 244 suffix: comparison 245 nsize: 4 246 args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 247 -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \ 248 -dm_plex_simplex 0 -snes_type ksponly -dm_view -ksp_type cg -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 249 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 250 251 test: 252 suffix: cuda 253 nsize: 4 254 requires: cuda 255 output_file: output/ex13_comparison.out 256 args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 257 -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \ 258 -dm_plex_simplex 0 -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 259 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijcusparse -dm_vec_type cuda 260 261 test: 262 suffix: kokkos_comp 263 nsize: 4 264 requires: !sycl kokkos_kernels 265 output_file: output/ex13_comparison.out 266 args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 267 -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \ 268 -dm_plex_simplex 0 -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 269 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijkokkos -dm_vec_type kokkos 270 271 test: 272 nsize: 4 273 requires: !sycl kokkos_kernels 274 suffix: kokkos 275 args: -dm_plex_dim 2 -dm_plex_box_faces 2,8 -petscpartitioner_type simple -petscpartitioner_simple_process_grid 2,1 \ 276 -petscpartitioner_simple_node_grid 2,1 -dm_plex_simplex 0 -potential_petscspace_degree 1 -dm_refine 1 -ksp_type cg -pc_type gamg -ksp_norm_type unpreconditioned \ 277 -pc_gamg_esteig_ksp_type cg -ksp_converged_reason -snes_rtol 1.e-4 -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 278 279 test: 280 suffix: aijmkl_comp 281 nsize: 4 282 requires: mkl_sparse 283 output_file: output/ex13_comparison.out 284 args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 285 -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \ 286 -dm_plex_simplex 0 -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 287 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl 288 289 test: 290 suffix: aijmkl_seq 291 nsize: 1 292 requires: mkl_sparse 293 TODO: broken (INDEFINITE PC) 294 args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_refine 1 -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_simplex 0 \ 295 -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_sym_graph 0 -pc_gamg_threshold -1 -pc_gamg_square_graph 10 -pc_gamg_process_eq_limit 400 \ 296 -pc_gamg_reuse_interpolation -pc_gamg_coarse_eq_limit 10 -pc_gamg_esteig_ksp_type cg -ksp_type cg -ksp_norm_type unpreconditioned \ 297 -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl -dm_vec_type standard 298 299 TEST*/ 300