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