Lines Matching +full:style +full:- +full:c

1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
4 // SPDX-License-Identifier: BSD-2-Clause
12 // The code is intentionally "raw", using only low-level communication primitives.
20 // ./bpsraw -problem bp1
21 // ./bpsraw -problem bp2
22 // ./bpsraw -problem bp3
23 // ./bpsraw -problem bp4
24 // ./bpsraw -problem bp5 -ceed /cpu/self
25 // ./bpsraw -problem bp6 -ceed /gpu/cuda
27 //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -q_extra 1 -ksp_max_it_clip 15,15
31 /// See bps.c for an implementation using DMPlex unstructured grids.
50 PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1. / (3 - d))); in Split3()
52 m[reverse ? 2 - d : d] = try; in Split3()
60 for (int d = 0; d < 3; d++) m_nodes[d] = degree * mesh_elem[d] + (i_rank[d] == p[d] - 1); in GlobalNodes()
64 // Dumb brute-force is easier to read in GlobalStart()
75 return -1; in GlobalStart()
83 for (int d = 0; d < 3; d++) m_nodes[d] = mesh_elem[d] * (P - 1) + 1; in CreateRestriction()
91 if (0) { // This is the C-style (i,j,k) ordering that I prefer in CreateRestriction()
92 …P + jj) * P + kk] = num_comp * (((i * (P - 1) + ii) * m_nodes[1] + (j * (P - 1) + jj)) * m_nodes[2… in CreateRestriction()
94 … * (jj + P * kk)] = num_comp * (((i * (P - 1) + ii) * m_nodes[1] + (j * (P - 1) + jj)) * m_nodes[2… in CreateRestriction()
115 VecScatter g_to_g_D; // global-to-global; only Dirichlet values
232 // Global-to-local in MatMult_Mass()
233 …PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REV… in MatMult_Mass()
234 …PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVER… in MatMult_Mass()
237 PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type)); in MatMult_Mass()
238 PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type)); in MatMult_Mass()
239 CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x); in MatMult_Mass()
240 CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y); in MatMult_Mass()
243 …CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMED… in MatMult_Mass()
246 CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL); in MatMult_Mass()
247 CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL); in MatMult_Mass()
248 PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x)); in MatMult_Mass()
249 PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y)); in MatMult_Mass()
251 // Local-to-global in MatMult_Mass()
254 …PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWAR… in MatMult_Mass()
255 …PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWARD)… in MatMult_Mass()
269 // Global-to-local in MatMult_Diff()
270 …PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_R… in MatMult_Diff()
271 …PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REV… in MatMult_Diff()
274 PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &x_mem_type)); in MatMult_Diff()
275 PetscCall(VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type)); in MatMult_Diff()
276 CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x); in MatMult_Diff()
277 CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y); in MatMult_Diff()
280 …CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, CEED_REQUEST_IMMED… in MatMult_Diff()
283 CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL); in MatMult_Diff()
284 CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL); in MatMult_Diff()
285 PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x)); in MatMult_Diff()
286 PetscCall(VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y)); in MatMult_Diff()
288 // Local-to-global in MatMult_Diff()
290 PetscCall(VecScatterBegin(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD)); in MatMult_Diff()
291 PetscCall(VecScatterEnd(op_apply_ctx->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD)); in MatMult_Diff()
292 …PetscCall(VecScatterBegin(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORW… in MatMult_Diff()
293 …PetscCall(VecScatterEnd(op_apply_ctx->l_to_g_0, op_apply_ctx->Y_loc, Y, ADD_VALUES, SCATTER_FORWAR… in MatMult_Diff()
306 CeedVectorCreate(op_apply_ctx->ceed, length, &collocated_error); in ComputeErrorMax()
308 // Global-to-local in ComputeErrorMax()
309 …PetscCall(VecScatterBegin(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REV… in ComputeErrorMax()
310 …PetscCall(VecScatterEnd(op_apply_ctx->l_to_g, X, op_apply_ctx->X_loc, INSERT_VALUES, SCATTER_REVER… in ComputeErrorMax()
313 PetscCall(VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, &mem_type)); in ComputeErrorMax()
314 CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x); in ComputeErrorMax()
317 CeedOperatorApply(op_error, op_apply_ctx->x_ceed, collocated_error, CEED_REQUEST_IMMEDIATE); in ComputeErrorMax()
320 CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), NULL); in ComputeErrorMax()
321 PetscCall(VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x)); in ComputeErrorMax()
331 PetscCall(MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX, op_apply_ctx->comm)); in ComputeErrorMax()
370 …PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEn… in main()
373 …PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, tes… in main()
375 …PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, … in main()
377 …PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solu… in main()
379 …PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &d… in main()
381 …PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra… in main()
382 …PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resourc… in main()
384 …PetscCall(PetscOptionsInt("-local", "Target number of locally owned nodes per process", NULL, loca… in main()
388 …PetscCall(PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during… in main()
429 rank_left -= i_rank[d] * pstride[d]; in main()
452 "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n" in main()
483 // Create local-to-global scatter in main()
502 …for (PetscInt i = 0, ir, ii; ir = i >= m_nodes[0], ii = i - ir * m_nodes[0], i < l_nodes[0]; i++) { in main()
503 …for (PetscInt j = 0, jr, jj; jr = j >= m_nodes[1], jj = j - jr * m_nodes[1], j < l_nodes[1]; j++) { in main()
504 …for (PetscInt k = 0, kr, kk; kr = k >= m_nodes[2], kk = k - kr * m_nodes[2], k < l_nodes[2]; k++) { in main()
523 …// Create global-to-global scatter for Dirichlet values (everything not in l_to_g_is_0, which is t… in main()
532 PetscCall(PetscMalloc1(x_end - x_start, &ind_D)); in main()
534 for (PetscInt i = 0; i < x_end - x_start; i++) { in main()
644 op_apply_ctx->comm = comm; in main()
645 op_apply_ctx->l_to_g = l_to_g; in main()
647 op_apply_ctx->l_to_g_0 = l_to_g_0; in main()
648 op_apply_ctx->g_to_g_D = g_to_g_D; in main()
650 op_apply_ctx->X_loc = X_loc; in main()
651 PetscCall(VecDuplicate(X_loc, &op_apply_ctx->Y_loc)); in main()
652 CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->x_ceed); in main()
653 CeedVectorCreate(ceed, l_size * num_comp_u, &op_apply_ctx->y_ceed); in main()
654 op_apply_ctx->op = op_apply; in main()
655 op_apply_ctx->q_data = q_data; in main()
656 op_apply_ctx->ceed = ceed; in main()
700 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); in main()
704 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1)); in main()
707 my_rt = MPI_Wtime() - my_rt_start; in main()
711 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[0])); in main()
713 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, ksp_max_it_clip[1])); in main()
721 // -- Performance logging in main()
725 // -- Solve in main()
728 my_rt = MPI_Wtime() - my_rt_start; in main()
730 // -- Performance logging in main()
743 if (!test_mode || reason < 0 || rnorm > 1e-8) { in main()
758 PetscReal tol = 5e-2; in main()
769 …PetscCall(PetscPrintf(comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g… in main()
770 1e-6 * gsize * its / rt_min)); in main()
787 PetscCall(VecDestroy(&op_apply_ctx->X_loc)); in main()
788 PetscCall(VecDestroy(&op_apply_ctx->Y_loc)); in main()
795 CeedVectorDestroy(&op_apply_ctx->x_ceed); in main()
796 CeedVectorDestroy(&op_apply_ctx->y_ceed); in main()
797 CeedVectorDestroy(&op_apply_ctx->q_data); in main()