13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3819eb1b3Sjeremylt // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5819eb1b3Sjeremylt // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7819eb1b3Sjeremylt 80c59ef15SJeremy L Thompson // libCEED + PETSc Example: CEED BPs 90c59ef15SJeremy L Thompson // 10*ea61e9acSJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps. 110c59ef15SJeremy L Thompson // 12cb32e2e7SValeria Barra // The code uses higher level communication protocols in DMPlex. 130c59ef15SJeremy L Thompson // 140c59ef15SJeremy L Thompson // Build with: 150c59ef15SJeremy L Thompson // 160c59ef15SJeremy L Thompson // make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 170c59ef15SJeremy L Thompson // 180c59ef15SJeremy L Thompson // Sample runs: 190c59ef15SJeremy L Thompson // 2032d2ee49SValeria Barra // ./bps -problem bp1 -degree 3 2128688798Sjeremylt // ./bps -problem bp2 -degree 3 2228688798Sjeremylt // ./bps -problem bp3 -degree 3 2328688798Sjeremylt // ./bps -problem bp4 -degree 3 2428688798Sjeremylt // ./bps -problem bp5 -degree 3 -ceed /cpu/self 2528688798Sjeremylt // ./bps -problem bp6 -degree 3 -ceed /gpu/cuda 260c59ef15SJeremy L Thompson // 272fbc6e41SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 15,15 286c88e6a2Srezgarshakeri //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -ksp_max_it_clip 50,50 -simplex 290c59ef15SJeremy L Thompson 300c59ef15SJeremy L Thompson /// @file 31cb32e2e7SValeria Barra /// CEED BPs example using PETSc with DMPlex 32cb32e2e7SValeria Barra /// See bpsraw.c for a "raw" implementation using a structured grid. 33cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc with DMPlex\n"; 340c59ef15SJeremy L Thompson 352b730f8bSJeremy L Thompson #include "bps.h" 362b730f8bSJeremy L Thompson 37636cccdbSjeremylt #include <ceed.h> 38636cccdbSjeremylt #include <petsc.h> 39636cccdbSjeremylt #include <petscdmplex.h> 40636cccdbSjeremylt #include <petscksp.h> 41636cccdbSjeremylt #include <petscsys.h> 422b730f8bSJeremy L Thompson #include <stdbool.h> 432b730f8bSJeremy L Thompson #include <string.h> 44636cccdbSjeremylt 45636cccdbSjeremylt #include "include/bpsproblemdata.h" 462b730f8bSJeremy L Thompson #include "include/libceedsetup.h" 472b730f8bSJeremy L Thompson #include "include/matops.h" 48636cccdbSjeremylt #include "include/petscutils.h" 49b8962995SJeremy L Thompson #include "include/petscversion.h" 50636cccdbSjeremylt #include "include/structs.h" 51636cccdbSjeremylt 52636cccdbSjeremylt #if PETSC_VERSION_LT(3, 12, 0) 53636cccdbSjeremylt #ifdef PETSC_HAVE_CUDA 54636cccdbSjeremylt #include <petsccuda.h> 55*ea61e9acSJeremy L Thompson // Note: With PETSc prior to version 3.12.0, providing the source path to include 'cublas_v2.h' will be needed to use 'petsccuda.h'. 56636cccdbSjeremylt #endif 57636cccdbSjeremylt #endif 580c59ef15SJeremy L Thompson 59d5b2ba77SJed Brown // ----------------------------------------------------------------------------- 601e284482Svaleriabarra // Main body of program, called in a loop for performance benchmarking purposes 611e284482Svaleriabarra // ----------------------------------------------------------------------------- 622b730f8bSJeremy L Thompson static PetscErrorCode RunWithDM(RunParams rp, DM dm, const char *ceed_resource) { 631e284482Svaleriabarra double my_rt_start, my_rt, rt_min, rt_max; 649b072555Sjeremylt PetscInt xl_size, l_size, g_size; 651e284482Svaleriabarra PetscScalar *r; 669b072555Sjeremylt Vec X, X_loc, rhs, rhs_loc; 679b072555Sjeremylt Mat mat_O; 68819eb1b3Sjeremylt KSP ksp; 696c88e6a2Srezgarshakeri OperatorApplyContext op_apply_ctx, op_error_ctx; 700c59ef15SJeremy L Thompson Ceed ceed; 719b072555Sjeremylt CeedData ceed_data; 729b072555Sjeremylt CeedQFunction qf_error; 739b072555Sjeremylt CeedOperator op_error; 749b072555Sjeremylt CeedVector rhs_ceed, target; 759b072555Sjeremylt VecType vec_type; 769b072555Sjeremylt PetscMemType mem_type; 771e284482Svaleriabarra 781e284482Svaleriabarra PetscFunctionBeginUser; 791e284482Svaleriabarra // Set up libCEED 809b072555Sjeremylt CeedInit(ceed_resource, &ceed); 819b072555Sjeremylt CeedMemType mem_type_backend; 829b072555Sjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend); 831e284482Svaleriabarra 842b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 859b072555Sjeremylt if (!vec_type) { // Not yet set by user -dm_vec_type 869b072555Sjeremylt switch (mem_type_backend) { 872b730f8bSJeremy L Thompson case CEED_MEM_HOST: 882b730f8bSJeremy L Thompson vec_type = VECSTANDARD; 892b730f8bSJeremy L Thompson break; 90b68a8d79SJed Brown case CEED_MEM_DEVICE: { 91b68a8d79SJed Brown const char *resolved; 92b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 939b072555Sjeremylt if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 942b730f8bSJeremy L Thompson else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 959b072555Sjeremylt else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 969b072555Sjeremylt else vec_type = VECSTANDARD; 971e284482Svaleriabarra } 98b68a8d79SJed Brown } 992b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm, vec_type)); 100b68a8d79SJed Brown } 101b68a8d79SJed Brown 102b68a8d79SJed Brown // Create global and local solution vectors 1032b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(dm, &X)); 1042b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(X, &l_size)); 1052b730f8bSJeremy L Thompson PetscCall(VecGetSize(X, &g_size)); 1062b730f8bSJeremy L Thompson PetscCall(DMCreateLocalVector(dm, &X_loc)); 1072b730f8bSJeremy L Thompson PetscCall(VecGetSize(X_loc, &xl_size)); 1082b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X, &rhs)); 1091e284482Svaleriabarra 1101e284482Svaleriabarra // Operator 1112b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_apply_ctx)); 1122b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_error_ctx)); 1132b730f8bSJeremy L Thompson PetscCall(MatCreateShell(rp->comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O)); 1142b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 1152b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag)); 1162b730f8bSJeremy L Thompson PetscCall(MatShellSetVecType(mat_O, vec_type)); 1171e284482Svaleriabarra 1181e284482Svaleriabarra // Print summary 1191e284482Svaleriabarra if (!rp->test_mode) { 1209b072555Sjeremylt PetscInt P = rp->degree + 1, Q = P + rp->q_extra; 1211e284482Svaleriabarra 1229b072555Sjeremylt const char *used_resource; 1239b072555Sjeremylt CeedGetResource(ceed, &used_resource); 1241e284482Svaleriabarra 1259b072555Sjeremylt VecType vec_type; 1262b730f8bSJeremy L Thompson PetscCall(VecGetType(X, &vec_type)); 1271e284482Svaleriabarra 1289b072555Sjeremylt PetscInt c_start, c_end; 1292b730f8bSJeremy L Thompson PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end)); 130de1229c5Srezgarshakeri DMPolytopeType cell_type; 1312b730f8bSJeremy L Thompson PetscCall(DMPlexGetCellType(dm, c_start, &cell_type)); 132de1229c5Srezgarshakeri CeedElemTopology elem_topo = ElemTopologyP2C(cell_type); 13353a0f73bSJed Brown PetscMPIInt comm_size; 1342b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(rp->comm, &comm_size)); 1352b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, 136990fdeb6SJeremy L Thompson "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n" 13753a0f73bSJed Brown " MPI:\n" 13853a0f73bSJed Brown " Hostname : %s\n" 13953a0f73bSJed Brown " Total ranks : %d\n" 14053a0f73bSJed Brown " Ranks per compute node : %d\n" 1411e284482Svaleriabarra " PETSc:\n" 1421e284482Svaleriabarra " PETSc Vec Type : %s\n" 1431e284482Svaleriabarra " libCEED:\n" 1441e284482Svaleriabarra " libCEED Backend : %s\n" 1451e284482Svaleriabarra " libCEED Backend MemType : %s\n" 1461e284482Svaleriabarra " Mesh:\n" 147751eb813Srezgarshakeri " Solution Order (P) : %" CeedInt_FMT "\n" 148751eb813Srezgarshakeri " Quadrature Order (Q) : %" CeedInt_FMT "\n" 14951ad7d5bSrezgarshakeri " Additional quadrature points (q_extra) : %" CeedInt_FMT "\n" 15008140895SJed Brown " Global nodes : %" PetscInt_FMT "\n" 15108140895SJed Brown " Local Elements : %" PetscInt_FMT "\n" 15251ad7d5bSrezgarshakeri " Element topology : %s\n" 15308140895SJed Brown " Owned nodes : %" PetscInt_FMT "\n" 15408140895SJed Brown " DoF per node : %" PetscInt_FMT "\n", 1552b730f8bSJeremy L Thompson rp->bp_choice + 1, rp->hostname, comm_size, rp->ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, 1562b730f8bSJeremy L Thompson Q, rp->q_extra, g_size / rp->num_comp_u, c_end - c_start, CeedElemTopologies[elem_topo], l_size / rp->num_comp_u, 1572b730f8bSJeremy L Thompson rp->num_comp_u)); 1581e284482Svaleriabarra } 1591e284482Svaleriabarra 1601e284482Svaleriabarra // Create RHS vector 1612b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc, &rhs_loc)); 1622b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs_loc)); 1632b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type)); 1649b072555Sjeremylt CeedVectorCreate(ceed, xl_size, &rhs_ceed); 1659b072555Sjeremylt CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 1661e284482Svaleriabarra 1672b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &ceed_data)); 1682b730f8bSJeremy L Thompson PetscCall(SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->q_extra, rp->dim, rp->num_comp_u, g_size, xl_size, bp_options[rp->bp_choice], 1692b730f8bSJeremy L Thompson ceed_data, true, rhs_ceed, &target)); 1701e284482Svaleriabarra 1711e284482Svaleriabarra // Gather RHS 1729b072555Sjeremylt CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 1732b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r)); 1742b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs)); 1752b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs)); 1769b072555Sjeremylt CeedVectorDestroy(&rhs_ceed); 1771e284482Svaleriabarra 1786a6c615bSJeremy L Thompson // Create the error QFunction 1792b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error, bp_options[rp->bp_choice].error_loc, &qf_error); 1809b072555Sjeremylt CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP); 1819b072555Sjeremylt CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE); 1822b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE); 18338f32c05Srezgarshakeri CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_INTERP); 1841e284482Svaleriabarra 1851e284482Svaleriabarra // Create the error operator 1862b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error); 1872b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 1882b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_COLLOCATED, target); 1892b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 1902b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 1911e284482Svaleriabarra 1926c88e6a2Srezgarshakeri // Set up apply operator context 1932b730f8bSJeremy L Thompson PetscCall(SetupApplyOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_apply_ctx)); 1942b730f8bSJeremy L Thompson PetscCall(KSPCreate(rp->comm, &ksp)); 1951e284482Svaleriabarra { 1961e284482Svaleriabarra PC pc; 1972b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc)); 1989b072555Sjeremylt if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) { 1992b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCJACOBI)); 200345f5e33Srezgarshakeri if (rp->simplex) { 2012b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 202345f5e33Srezgarshakeri } else { 2032b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 204345f5e33Srezgarshakeri } 2051e284482Svaleriabarra } else { 2062b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCNONE)); 2071e284482Svaleriabarra } 2082b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp, KSPCG)); 2092b730f8bSJeremy L Thompson PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 2102b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 2111e284482Svaleriabarra } 2122b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp, mat_O, mat_O)); 2131e284482Svaleriabarra 214da9108adSvaleriabarra // First run's performance log is not considered for benchmarking purposes 2152b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 2161e284482Svaleriabarra my_rt_start = MPI_Wtime(); 2172b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X)); 2181e284482Svaleriabarra my_rt = MPI_Wtime() - my_rt_start; 2192b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm)); 2201e284482Svaleriabarra // Set maxits based on first iteration timing 2211e284482Svaleriabarra if (my_rt > 0.02) { 2222b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[0])); 2231e284482Svaleriabarra } else { 2242b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[1])); 2251e284482Svaleriabarra } 2262b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp)); 2271e284482Svaleriabarra 2281e284482Svaleriabarra // Timed solve 2292b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X)); 2302b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)ksp)); 2311e284482Svaleriabarra 2321e284482Svaleriabarra // -- Performance logging 2332b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(rp->solve_stage)); 2341e284482Svaleriabarra 2351e284482Svaleriabarra // -- Solve 2361e284482Svaleriabarra my_rt_start = MPI_Wtime(); 2372b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X)); 2381e284482Svaleriabarra my_rt = MPI_Wtime() - my_rt_start; 2391e284482Svaleriabarra 2401e284482Svaleriabarra // -- Performance logging 2412b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 2421e284482Svaleriabarra 2431e284482Svaleriabarra // Output results 2441e284482Svaleriabarra { 2459b072555Sjeremylt KSPType ksp_type; 2461e284482Svaleriabarra KSPConvergedReason reason; 2471e284482Svaleriabarra PetscReal rnorm; 2481e284482Svaleriabarra PetscInt its; 2492b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp, &ksp_type)); 2502b730f8bSJeremy L Thompson PetscCall(KSPGetConvergedReason(ksp, &reason)); 2512b730f8bSJeremy L Thompson PetscCall(KSPGetIterationNumber(ksp, &its)); 2522b730f8bSJeremy L Thompson PetscCall(KSPGetResidualNorm(ksp, &rnorm)); 2531e284482Svaleriabarra if (!rp->test_mode || reason < 0 || rnorm > 1e-8) { 2542b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, 2551e284482Svaleriabarra " KSP:\n" 2561e284482Svaleriabarra " KSP Type : %s\n" 2571e284482Svaleriabarra " KSP Convergence : %s\n" 258a9b2c5ddSrezgarshakeri " Total KSP Iterations : %" PetscInt_FMT "\n" 2591e284482Svaleriabarra " Final rnorm : %e\n", 2602b730f8bSJeremy L Thompson ksp_type, KSPConvergedReasons[reason], its, (double)rnorm)); 2611e284482Svaleriabarra } 2621e284482Svaleriabarra if (!rp->test_mode) { 2632b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, " Performance:\n")); 2641e284482Svaleriabarra } 2651e284482Svaleriabarra { 2666c88e6a2Srezgarshakeri // Set up error operator context 2672b730f8bSJeremy L Thompson PetscCall(SetupErrorOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx)); 26838f32c05Srezgarshakeri PetscScalar l2_error; 2692b730f8bSJeremy L Thompson PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx)); 2700a8fc04aSrezgarshakeri PetscReal tol = 5e-2; 27138f32c05Srezgarshakeri if (!rp->test_mode || l2_error > tol) { 2722b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm)); 2732b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm)); 2742b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, 27538f32c05Srezgarshakeri " L2 Error : %e\n" 2761e284482Svaleriabarra " CG Solve Time : %g (%g) sec\n", 2772b730f8bSJeremy L Thompson (double)l2_error, rt_max, rt_min)); 2781e284482Svaleriabarra } 2791e284482Svaleriabarra } 2804c583f1fSvaleriabarra if (!rp->test_mode) { 2812b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g_size * its / rt_max, 2822b730f8bSJeremy L Thompson 1e-6 * g_size * its / rt_min)); 2831e284482Svaleriabarra } 2841e284482Svaleriabarra } 2851e284482Svaleriabarra 2861e284482Svaleriabarra if (rp->write_solution) { 2879b072555Sjeremylt PetscViewer vtk_viewer_soln; 2881e284482Svaleriabarra 2892b730f8bSJeremy L Thompson PetscCall(PetscViewerCreate(rp->comm, &vtk_viewer_soln)); 2902b730f8bSJeremy L Thompson PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK)); 2912b730f8bSJeremy L Thompson PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu")); 2922b730f8bSJeremy L Thompson PetscCall(VecView(X, vtk_viewer_soln)); 2932b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&vtk_viewer_soln)); 2941e284482Svaleriabarra } 2951e284482Svaleriabarra 2961e284482Svaleriabarra // Cleanup 2972b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X)); 2982b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X_loc)); 2992b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_apply_ctx->Y_loc)); 3002b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_error_ctx->Y_loc)); 3012b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_O)); 3022b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx)); 3032b730f8bSJeremy L Thompson PetscCall(PetscFree(op_error_ctx)); 3042b730f8bSJeremy L Thompson PetscCall(CeedDataDestroy(0, ceed_data)); 3051e284482Svaleriabarra 3062b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs)); 3072b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs_loc)); 3082b730f8bSJeremy L Thompson PetscCall(KSPDestroy(&ksp)); 3091e284482Svaleriabarra CeedVectorDestroy(&target); 3109b072555Sjeremylt CeedQFunctionDestroy(&qf_error); 3119b072555Sjeremylt CeedOperatorDestroy(&op_error); 3121e284482Svaleriabarra CeedDestroy(&ceed); 3131e284482Svaleriabarra PetscFunctionReturn(0); 3141e284482Svaleriabarra } 3151e284482Svaleriabarra 3162b730f8bSJeremy L Thompson static PetscErrorCode Run(RunParams rp, PetscInt num_resources, char *const *ceed_resources, PetscInt num_bp_choices, const BPType *bp_choices) { 3175f284d84SJed Brown DM dm; 3185f284d84SJed Brown 3195f284d84SJed Brown PetscFunctionBeginUser; 3205f284d84SJed Brown // Setup DM 3212b730f8bSJeremy L Thompson PetscCall(CreateDistributedDM(rp, &dm)); 3225f284d84SJed Brown 3239b072555Sjeremylt for (PetscInt b = 0; b < num_bp_choices; b++) { 3245f284d84SJed Brown DM dm_deg; 3259b072555Sjeremylt VecType vec_type; 3269b072555Sjeremylt PetscInt q_extra = rp->q_extra; 3279b072555Sjeremylt rp->bp_choice = bp_choices[b]; 3289b072555Sjeremylt rp->num_comp_u = bp_options[rp->bp_choice].num_comp_u; 3299b072555Sjeremylt rp->q_extra = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra; 3302b730f8bSJeremy L Thompson PetscCall(DMClone(dm, &dm_deg)); 3312b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 3322b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_deg, vec_type)); 3335f284d84SJed Brown // Create DM 334e83e87a5Sjeremylt PetscInt dim; 3352b730f8bSJeremy L Thompson PetscCall(DMGetDimension(dm_deg, &dim)); 3362b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(dm_deg, rp->degree, rp->q_extra, rp->num_comp_u, dim, bp_options[rp->bp_choice].enforce_bc)); 3375f284d84SJed Brown for (PetscInt r = 0; r < num_resources; r++) { 3382b730f8bSJeremy L Thompson PetscCall(RunWithDM(rp, dm_deg, ceed_resources[r])); 3395f284d84SJed Brown } 3402b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_deg)); 3419b072555Sjeremylt rp->q_extra = q_extra; 3425f284d84SJed Brown } 3435f284d84SJed Brown 3442b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm)); 3455f284d84SJed Brown PetscFunctionReturn(0); 3465f284d84SJed Brown } 3475f284d84SJed Brown 3481e284482Svaleriabarra int main(int argc, char **argv) { 3492b730f8bSJeremy L Thompson PetscInt comm_size; 3501e284482Svaleriabarra RunParams rp; 35153a0f73bSJed Brown MPI_Comm comm; 35253a0f73bSJed Brown char filename[PETSC_MAX_PATH_LEN]; 3539b072555Sjeremylt char *ceed_resources[30]; 3549b072555Sjeremylt PetscInt num_ceed_resources = 30; 35553a0f73bSJed Brown char hostname[PETSC_MAX_PATH_LEN]; 3561e284482Svaleriabarra 3579b072555Sjeremylt PetscInt dim = 3, mesh_elem[3] = {3, 3, 3}; 3582b730f8bSJeremy L Thompson PetscInt num_degrees = 30, degree[30] = {}, num_local_nodes = 2, local_nodes[2] = {}; 3599b072555Sjeremylt PetscMPIInt ranks_per_node; 360c36f77d8SJed Brown PetscBool degree_set; 3619b072555Sjeremylt BPType bp_choices[10]; 3629b072555Sjeremylt PetscInt num_bp_choices = 10; 3630c59ef15SJeremy L Thompson 3641e284482Svaleriabarra // Initialize PETSc 3652b730f8bSJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3660c59ef15SJeremy L Thompson comm = PETSC_COMM_WORLD; 3672b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(comm, &comm_size)); 36853a0f73bSJed Brown #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 36953a0f73bSJed Brown { 37053a0f73bSJed Brown MPI_Comm splitcomm; 3712b730f8bSJeremy L Thompson PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm)); 3722b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node)); 3732b730f8bSJeremy L Thompson PetscCall(MPI_Comm_free(&splitcomm)); 37453a0f73bSJed Brown } 37553a0f73bSJed Brown #else 3769b072555Sjeremylt ranks_per_node = -1; // Unknown 37753a0f73bSJed Brown #endif 378cb32e2e7SValeria Barra 379c36f77d8SJed Brown // Setup all parameters needed in Run() 3802b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &rp)); 381c36f77d8SJed Brown rp->comm = comm; 382c36f77d8SJed Brown 38332d2ee49SValeria Barra // Read command line options 38467490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 385565a3730SJed Brown { 386565a3730SJed Brown PetscBool set; 3872b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set)); 388565a3730SJed Brown if (!set) { 3899b072555Sjeremylt bp_choices[0] = CEED_BP1; 3909b072555Sjeremylt num_bp_choices = 1; 391565a3730SJed Brown } 392565a3730SJed Brown } 393c36f77d8SJed Brown rp->test_mode = PETSC_FALSE; 3942b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, rp->test_mode, &rp->test_mode, NULL)); 395c36f77d8SJed Brown rp->write_solution = PETSC_FALSE; 3962b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, rp->write_solution, &rp->write_solution, NULL)); 397de1229c5Srezgarshakeri rp->simplex = PETSC_FALSE; 3982b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, rp->simplex, &rp->simplex, NULL)); 399eb6a3b92Srezgarshakeri if ((bp_choices[0] == CEED_BP5 || bp_choices[0] == CEED_BP6) && (rp->simplex)) { 4002b730f8bSJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex"); 401eb6a3b92Srezgarshakeri } 402c36f77d8SJed Brown degree[0] = rp->test_mode ? 3 : 2; 4032b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-degree", "Polynomial degree of tensor product basis", NULL, degree, &num_degrees, °ree_set)); 4042b730f8bSJeremy L Thompson if (!degree_set) num_degrees = 1; 4059b072555Sjeremylt rp->q_extra = PETSC_DECIDE; 4062b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points (-1 for auto)", NULL, rp->q_extra, &rp->q_extra, NULL)); 407565a3730SJed Brown { 408565a3730SJed Brown PetscBool set; 4092b730f8bSJeremy L Thompson PetscCall(PetscOptionsStringArray("-ceed", "CEED resource specifier (comma-separated list)", NULL, ceed_resources, &num_ceed_resources, &set)); 410565a3730SJed Brown if (!set) { 4112b730f8bSJeremy L Thompson PetscCall(PetscStrallocpy("/cpu/self", &ceed_resources[0])); 4129b072555Sjeremylt num_ceed_resources = 1; 413565a3730SJed Brown } 414565a3730SJed Brown } 4152b730f8bSJeremy L Thompson PetscCall(PetscGetHostName(hostname, sizeof hostname)); 4162b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL)); 417c36f77d8SJed Brown rp->read_mesh = PETSC_FALSE; 4182b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &rp->read_mesh)); 419c36f77d8SJed Brown rp->filename = filename; 420c36f77d8SJed Brown if (!rp->read_mesh) { 421cb32e2e7SValeria Barra PetscInt tmp = dim; 4222b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL)); 423cb32e2e7SValeria Barra } 4249b072555Sjeremylt local_nodes[0] = 1000; 4252b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-local_nodes", 42653a0f73bSJed Brown "Target number of locally owned nodes per " 42753a0f73bSJed Brown "process (single value or min,max)", 4282b730f8bSJeremy L Thompson NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes)); 4292b730f8bSJeremy L Thompson if (num_local_nodes < 2) local_nodes[1] = 2 * local_nodes[0]; 430c36f77d8SJed Brown { 431c36f77d8SJed Brown PetscInt two = 2; 432c36f77d8SJed Brown rp->ksp_max_it_clip[0] = 5; 433c36f77d8SJed Brown rp->ksp_max_it_clip[1] = 20; 4342b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during benchmarking", NULL, rp->ksp_max_it_clip, &two, 4352b730f8bSJeremy L Thompson NULL)); 436c36f77d8SJed Brown } 437da9108adSvaleriabarra if (!degree_set) { 4389b072555Sjeremylt PetscInt max_degree = 8; 4392b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-max_degree", "Range of degrees [1, max_degree] to run with", NULL, max_degree, &max_degree, NULL)); 4402b730f8bSJeremy L Thompson for (PetscInt i = 0; i < max_degree; i++) degree[i] = i + 1; 4419b072555Sjeremylt num_degrees = max_degree; 44253a0f73bSJed Brown } 44353a0f73bSJed Brown { 44453a0f73bSJed Brown PetscBool flg; 4459b072555Sjeremylt PetscInt p = ranks_per_node; 4462b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg)); 4479b072555Sjeremylt if (flg) ranks_per_node = p; 44853a0f73bSJed Brown } 4499396343dSjeremylt 45067490bc6SJeremy L Thompson PetscOptionsEnd(); 451cb32e2e7SValeria Barra 4521e284482Svaleriabarra // Register PETSc logging stage 4532b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("Solve Stage", &rp->solve_stage)); 45409a940d7Sjeremylt 45553a0f73bSJed Brown rp->hostname = hostname; 4561e284482Svaleriabarra rp->dim = dim; 4579b072555Sjeremylt rp->mesh_elem = mesh_elem; 4589b072555Sjeremylt rp->ranks_per_node = ranks_per_node; 459cb32e2e7SValeria Barra 46053a0f73bSJed Brown for (PetscInt d = 0; d < num_degrees; d++) { 46153a0f73bSJed Brown PetscInt deg = degree[d]; 4629b072555Sjeremylt for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) { 46353a0f73bSJed Brown rp->degree = deg; 4649b072555Sjeremylt rp->local_nodes = n; 4652b730f8bSJeremy L Thompson PetscCall(Run(rp, num_ceed_resources, ceed_resources, num_bp_choices, bp_choices)); 466565a3730SJed Brown } 467565a3730SJed Brown } 4681e284482Svaleriabarra // Clear memory 4692b730f8bSJeremy L Thompson PetscCall(PetscFree(rp)); 4709b072555Sjeremylt for (PetscInt i = 0; i < num_ceed_resources; i++) { 4712b730f8bSJeremy L Thompson PetscCall(PetscFree(ceed_resources[i])); 472565a3730SJed Brown } 4730c59ef15SJeremy L Thompson return PetscFinalize(); 4740c59ef15SJeremy L Thompson } 475