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 // 100c59ef15SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the 110c59ef15SJeremy L Thompson // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps. 120c59ef15SJeremy L Thompson // 13cb32e2e7SValeria Barra // The code uses higher level communication protocols in DMPlex. 140c59ef15SJeremy L Thompson // 150c59ef15SJeremy L Thompson // Build with: 160c59ef15SJeremy L Thompson // 170c59ef15SJeremy L Thompson // make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 180c59ef15SJeremy L Thompson // 190c59ef15SJeremy L Thompson // Sample runs: 200c59ef15SJeremy L Thompson // 2132d2ee49SValeria Barra // ./bps -problem bp1 -degree 3 2228688798Sjeremylt // ./bps -problem bp2 -degree 3 2328688798Sjeremylt // ./bps -problem bp3 -degree 3 2428688798Sjeremylt // ./bps -problem bp4 -degree 3 2528688798Sjeremylt // ./bps -problem bp5 -degree 3 -ceed /cpu/self 2628688798Sjeremylt // ./bps -problem bp6 -degree 3 -ceed /gpu/cuda 270c59ef15SJeremy L Thompson // 282fbc6e41SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 15,15 296c88e6a2Srezgarshakeri //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -ksp_max_it_clip 50,50 -simplex 300c59ef15SJeremy L Thompson 310c59ef15SJeremy L Thompson /// @file 32cb32e2e7SValeria Barra /// CEED BPs example using PETSc with DMPlex 33cb32e2e7SValeria Barra /// See bpsraw.c for a "raw" implementation using a structured grid. 34cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc with DMPlex\n"; 350c59ef15SJeremy L Thompson 36*2b730f8bSJeremy L Thompson #include "bps.h" 37*2b730f8bSJeremy L Thompson 38636cccdbSjeremylt #include <ceed.h> 39636cccdbSjeremylt #include <petsc.h> 40636cccdbSjeremylt #include <petscdmplex.h> 41636cccdbSjeremylt #include <petscksp.h> 42636cccdbSjeremylt #include <petscsys.h> 43*2b730f8bSJeremy L Thompson #include <stdbool.h> 44*2b730f8bSJeremy L Thompson #include <string.h> 45636cccdbSjeremylt 46636cccdbSjeremylt #include "include/bpsproblemdata.h" 47*2b730f8bSJeremy L Thompson #include "include/libceedsetup.h" 48*2b730f8bSJeremy L Thompson #include "include/matops.h" 49636cccdbSjeremylt #include "include/petscutils.h" 50b8962995SJeremy L Thompson #include "include/petscversion.h" 51636cccdbSjeremylt #include "include/structs.h" 52636cccdbSjeremylt 53636cccdbSjeremylt #if PETSC_VERSION_LT(3, 12, 0) 54636cccdbSjeremylt #ifdef PETSC_HAVE_CUDA 55636cccdbSjeremylt #include <petsccuda.h> 56636cccdbSjeremylt // Note: With PETSc prior to version 3.12.0, providing the source path to 57636cccdbSjeremylt // include 'cublas_v2.h' will be needed to use 'petsccuda.h'. 58636cccdbSjeremylt #endif 59636cccdbSjeremylt #endif 600c59ef15SJeremy L Thompson 61d5b2ba77SJed Brown // ----------------------------------------------------------------------------- 621e284482Svaleriabarra // Main body of program, called in a loop for performance benchmarking purposes 631e284482Svaleriabarra // ----------------------------------------------------------------------------- 64*2b730f8bSJeremy L Thompson static PetscErrorCode RunWithDM(RunParams rp, DM dm, const char *ceed_resource) { 651e284482Svaleriabarra double my_rt_start, my_rt, rt_min, rt_max; 669b072555Sjeremylt PetscInt xl_size, l_size, g_size; 671e284482Svaleriabarra PetscScalar *r; 689b072555Sjeremylt Vec X, X_loc, rhs, rhs_loc; 699b072555Sjeremylt Mat mat_O; 70819eb1b3Sjeremylt KSP ksp; 716c88e6a2Srezgarshakeri OperatorApplyContext op_apply_ctx, op_error_ctx; 720c59ef15SJeremy L Thompson Ceed ceed; 739b072555Sjeremylt CeedData ceed_data; 749b072555Sjeremylt CeedQFunction qf_error; 759b072555Sjeremylt CeedOperator op_error; 769b072555Sjeremylt CeedVector rhs_ceed, target; 779b072555Sjeremylt VecType vec_type; 789b072555Sjeremylt PetscMemType mem_type; 791e284482Svaleriabarra 801e284482Svaleriabarra PetscFunctionBeginUser; 811e284482Svaleriabarra // Set up libCEED 829b072555Sjeremylt CeedInit(ceed_resource, &ceed); 839b072555Sjeremylt CeedMemType mem_type_backend; 849b072555Sjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend); 851e284482Svaleriabarra 86*2b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 879b072555Sjeremylt if (!vec_type) { // Not yet set by user -dm_vec_type 889b072555Sjeremylt switch (mem_type_backend) { 89*2b730f8bSJeremy L Thompson case CEED_MEM_HOST: 90*2b730f8bSJeremy L Thompson vec_type = VECSTANDARD; 91*2b730f8bSJeremy L Thompson break; 92b68a8d79SJed Brown case CEED_MEM_DEVICE: { 93b68a8d79SJed Brown const char *resolved; 94b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 959b072555Sjeremylt if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 96*2b730f8bSJeremy L Thompson else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 979b072555Sjeremylt else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 989b072555Sjeremylt else vec_type = VECSTANDARD; 991e284482Svaleriabarra } 100b68a8d79SJed Brown } 101*2b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm, vec_type)); 102b68a8d79SJed Brown } 103b68a8d79SJed Brown 104b68a8d79SJed Brown // Create global and local solution vectors 105*2b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(dm, &X)); 106*2b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(X, &l_size)); 107*2b730f8bSJeremy L Thompson PetscCall(VecGetSize(X, &g_size)); 108*2b730f8bSJeremy L Thompson PetscCall(DMCreateLocalVector(dm, &X_loc)); 109*2b730f8bSJeremy L Thompson PetscCall(VecGetSize(X_loc, &xl_size)); 110*2b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X, &rhs)); 1111e284482Svaleriabarra 1121e284482Svaleriabarra // Operator 113*2b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_apply_ctx)); 114*2b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_error_ctx)); 115*2b730f8bSJeremy L Thompson PetscCall(MatCreateShell(rp->comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O)); 116*2b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 117*2b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag)); 118*2b730f8bSJeremy L Thompson PetscCall(MatShellSetVecType(mat_O, vec_type)); 1191e284482Svaleriabarra 1201e284482Svaleriabarra // Print summary 1211e284482Svaleriabarra if (!rp->test_mode) { 1229b072555Sjeremylt PetscInt P = rp->degree + 1, Q = P + rp->q_extra; 1231e284482Svaleriabarra 1249b072555Sjeremylt const char *used_resource; 1259b072555Sjeremylt CeedGetResource(ceed, &used_resource); 1261e284482Svaleriabarra 1279b072555Sjeremylt VecType vec_type; 128*2b730f8bSJeremy L Thompson PetscCall(VecGetType(X, &vec_type)); 1291e284482Svaleriabarra 1309b072555Sjeremylt PetscInt c_start, c_end; 131*2b730f8bSJeremy L Thompson PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end)); 132de1229c5Srezgarshakeri DMPolytopeType cell_type; 133*2b730f8bSJeremy L Thompson PetscCall(DMPlexGetCellType(dm, c_start, &cell_type)); 134de1229c5Srezgarshakeri CeedElemTopology elem_topo = ElemTopologyP2C(cell_type); 13553a0f73bSJed Brown PetscMPIInt comm_size; 136*2b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(rp->comm, &comm_size)); 137*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, 138990fdeb6SJeremy L Thompson "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n" 13953a0f73bSJed Brown " MPI:\n" 14053a0f73bSJed Brown " Hostname : %s\n" 14153a0f73bSJed Brown " Total ranks : %d\n" 14253a0f73bSJed Brown " Ranks per compute node : %d\n" 1431e284482Svaleriabarra " PETSc:\n" 1441e284482Svaleriabarra " PETSc Vec Type : %s\n" 1451e284482Svaleriabarra " libCEED:\n" 1461e284482Svaleriabarra " libCEED Backend : %s\n" 1471e284482Svaleriabarra " libCEED Backend MemType : %s\n" 1481e284482Svaleriabarra " Mesh:\n" 149751eb813Srezgarshakeri " Solution Order (P) : %" CeedInt_FMT "\n" 150751eb813Srezgarshakeri " Quadrature Order (Q) : %" CeedInt_FMT "\n" 15151ad7d5bSrezgarshakeri " Additional quadrature points (q_extra) : %" CeedInt_FMT "\n" 15208140895SJed Brown " Global nodes : %" PetscInt_FMT "\n" 15308140895SJed Brown " Local Elements : %" PetscInt_FMT "\n" 15451ad7d5bSrezgarshakeri " Element topology : %s\n" 15508140895SJed Brown " Owned nodes : %" PetscInt_FMT "\n" 15608140895SJed Brown " DoF per node : %" PetscInt_FMT "\n", 157*2b730f8bSJeremy L Thompson rp->bp_choice + 1, rp->hostname, comm_size, rp->ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, 158*2b730f8bSJeremy 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, 159*2b730f8bSJeremy L Thompson rp->num_comp_u)); 1601e284482Svaleriabarra } 1611e284482Svaleriabarra 1621e284482Svaleriabarra // Create RHS vector 163*2b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc, &rhs_loc)); 164*2b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs_loc)); 165*2b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type)); 1669b072555Sjeremylt CeedVectorCreate(ceed, xl_size, &rhs_ceed); 1679b072555Sjeremylt CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 1681e284482Svaleriabarra 169*2b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &ceed_data)); 170*2b730f8bSJeremy 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], 171*2b730f8bSJeremy L Thompson ceed_data, true, rhs_ceed, &target)); 1721e284482Svaleriabarra 1731e284482Svaleriabarra // Gather RHS 1749b072555Sjeremylt CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 175*2b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r)); 176*2b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs)); 177*2b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs)); 1789b072555Sjeremylt CeedVectorDestroy(&rhs_ceed); 1791e284482Svaleriabarra 1806a6c615bSJeremy L Thompson // Create the error QFunction 181*2b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error, bp_options[rp->bp_choice].error_loc, &qf_error); 1829b072555Sjeremylt CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP); 1839b072555Sjeremylt CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE); 184*2b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE); 18538f32c05Srezgarshakeri CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_INTERP); 1861e284482Svaleriabarra 1871e284482Svaleriabarra // Create the error operator 188*2b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error); 189*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 190*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_COLLOCATED, target); 191*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 192*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 1931e284482Svaleriabarra 1946c88e6a2Srezgarshakeri // Set up apply operator context 195*2b730f8bSJeremy L Thompson PetscCall(SetupApplyOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_apply_ctx)); 196*2b730f8bSJeremy L Thompson PetscCall(KSPCreate(rp->comm, &ksp)); 1971e284482Svaleriabarra { 1981e284482Svaleriabarra PC pc; 199*2b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc)); 2009b072555Sjeremylt if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) { 201*2b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCJACOBI)); 202345f5e33Srezgarshakeri if (rp->simplex) { 203*2b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 204345f5e33Srezgarshakeri } else { 205*2b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 206345f5e33Srezgarshakeri } 2071e284482Svaleriabarra } else { 208*2b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCNONE)); 2091e284482Svaleriabarra } 210*2b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp, KSPCG)); 211*2b730f8bSJeremy L Thompson PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 212*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 2131e284482Svaleriabarra } 214*2b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp, mat_O, mat_O)); 2151e284482Svaleriabarra 216da9108adSvaleriabarra // First run's performance log is not considered for benchmarking purposes 217*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 2181e284482Svaleriabarra my_rt_start = MPI_Wtime(); 219*2b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X)); 2201e284482Svaleriabarra my_rt = MPI_Wtime() - my_rt_start; 221*2b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm)); 2221e284482Svaleriabarra // Set maxits based on first iteration timing 2231e284482Svaleriabarra if (my_rt > 0.02) { 224*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[0])); 2251e284482Svaleriabarra } else { 226*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[1])); 2271e284482Svaleriabarra } 228*2b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp)); 2291e284482Svaleriabarra 2301e284482Svaleriabarra // Timed solve 231*2b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X)); 232*2b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)ksp)); 2331e284482Svaleriabarra 2341e284482Svaleriabarra // -- Performance logging 235*2b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(rp->solve_stage)); 2361e284482Svaleriabarra 2371e284482Svaleriabarra // -- Solve 2381e284482Svaleriabarra my_rt_start = MPI_Wtime(); 239*2b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X)); 2401e284482Svaleriabarra my_rt = MPI_Wtime() - my_rt_start; 2411e284482Svaleriabarra 2421e284482Svaleriabarra // -- Performance logging 243*2b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 2441e284482Svaleriabarra 2451e284482Svaleriabarra // Output results 2461e284482Svaleriabarra { 2479b072555Sjeremylt KSPType ksp_type; 2481e284482Svaleriabarra KSPConvergedReason reason; 2491e284482Svaleriabarra PetscReal rnorm; 2501e284482Svaleriabarra PetscInt its; 251*2b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp, &ksp_type)); 252*2b730f8bSJeremy L Thompson PetscCall(KSPGetConvergedReason(ksp, &reason)); 253*2b730f8bSJeremy L Thompson PetscCall(KSPGetIterationNumber(ksp, &its)); 254*2b730f8bSJeremy L Thompson PetscCall(KSPGetResidualNorm(ksp, &rnorm)); 2551e284482Svaleriabarra if (!rp->test_mode || reason < 0 || rnorm > 1e-8) { 256*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, 2571e284482Svaleriabarra " KSP:\n" 2581e284482Svaleriabarra " KSP Type : %s\n" 2591e284482Svaleriabarra " KSP Convergence : %s\n" 260a9b2c5ddSrezgarshakeri " Total KSP Iterations : %" PetscInt_FMT "\n" 2611e284482Svaleriabarra " Final rnorm : %e\n", 262*2b730f8bSJeremy L Thompson ksp_type, KSPConvergedReasons[reason], its, (double)rnorm)); 2631e284482Svaleriabarra } 2641e284482Svaleriabarra if (!rp->test_mode) { 265*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, " Performance:\n")); 2661e284482Svaleriabarra } 2671e284482Svaleriabarra { 2686c88e6a2Srezgarshakeri // Set up error operator context 269*2b730f8bSJeremy L Thompson PetscCall(SetupErrorOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx)); 27038f32c05Srezgarshakeri PetscScalar l2_error; 271*2b730f8bSJeremy L Thompson PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx)); 2720a8fc04aSrezgarshakeri PetscReal tol = 5e-2; 27338f32c05Srezgarshakeri if (!rp->test_mode || l2_error > tol) { 274*2b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm)); 275*2b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm)); 276*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, 27738f32c05Srezgarshakeri " L2 Error : %e\n" 2781e284482Svaleriabarra " CG Solve Time : %g (%g) sec\n", 279*2b730f8bSJeremy L Thompson (double)l2_error, rt_max, rt_min)); 2801e284482Svaleriabarra } 2811e284482Svaleriabarra } 2824c583f1fSvaleriabarra if (!rp->test_mode) { 283*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(rp->comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g_size * its / rt_max, 284*2b730f8bSJeremy L Thompson 1e-6 * g_size * its / rt_min)); 2851e284482Svaleriabarra } 2861e284482Svaleriabarra } 2871e284482Svaleriabarra 2881e284482Svaleriabarra if (rp->write_solution) { 2899b072555Sjeremylt PetscViewer vtk_viewer_soln; 2901e284482Svaleriabarra 291*2b730f8bSJeremy L Thompson PetscCall(PetscViewerCreate(rp->comm, &vtk_viewer_soln)); 292*2b730f8bSJeremy L Thompson PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK)); 293*2b730f8bSJeremy L Thompson PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu")); 294*2b730f8bSJeremy L Thompson PetscCall(VecView(X, vtk_viewer_soln)); 295*2b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&vtk_viewer_soln)); 2961e284482Svaleriabarra } 2971e284482Svaleriabarra 2981e284482Svaleriabarra // Cleanup 299*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X)); 300*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X_loc)); 301*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_apply_ctx->Y_loc)); 302*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_error_ctx->Y_loc)); 303*2b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_O)); 304*2b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx)); 305*2b730f8bSJeremy L Thompson PetscCall(PetscFree(op_error_ctx)); 306*2b730f8bSJeremy L Thompson PetscCall(CeedDataDestroy(0, ceed_data)); 3071e284482Svaleriabarra 308*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs)); 309*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs_loc)); 310*2b730f8bSJeremy L Thompson PetscCall(KSPDestroy(&ksp)); 3111e284482Svaleriabarra CeedVectorDestroy(&target); 3129b072555Sjeremylt CeedQFunctionDestroy(&qf_error); 3139b072555Sjeremylt CeedOperatorDestroy(&op_error); 3141e284482Svaleriabarra CeedDestroy(&ceed); 3151e284482Svaleriabarra PetscFunctionReturn(0); 3161e284482Svaleriabarra } 3171e284482Svaleriabarra 318*2b730f8bSJeremy L Thompson static PetscErrorCode Run(RunParams rp, PetscInt num_resources, char *const *ceed_resources, PetscInt num_bp_choices, const BPType *bp_choices) { 3195f284d84SJed Brown DM dm; 3205f284d84SJed Brown 3215f284d84SJed Brown PetscFunctionBeginUser; 3225f284d84SJed Brown // Setup DM 323*2b730f8bSJeremy L Thompson PetscCall(CreateDistributedDM(rp, &dm)); 3245f284d84SJed Brown 3259b072555Sjeremylt for (PetscInt b = 0; b < num_bp_choices; b++) { 3265f284d84SJed Brown DM dm_deg; 3279b072555Sjeremylt VecType vec_type; 3289b072555Sjeremylt PetscInt q_extra = rp->q_extra; 3299b072555Sjeremylt rp->bp_choice = bp_choices[b]; 3309b072555Sjeremylt rp->num_comp_u = bp_options[rp->bp_choice].num_comp_u; 3319b072555Sjeremylt rp->q_extra = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra; 332*2b730f8bSJeremy L Thompson PetscCall(DMClone(dm, &dm_deg)); 333*2b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 334*2b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_deg, vec_type)); 3355f284d84SJed Brown // Create DM 336e83e87a5Sjeremylt PetscInt dim; 337*2b730f8bSJeremy L Thompson PetscCall(DMGetDimension(dm_deg, &dim)); 338*2b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(dm_deg, rp->degree, rp->q_extra, rp->num_comp_u, dim, bp_options[rp->bp_choice].enforce_bc)); 3395f284d84SJed Brown for (PetscInt r = 0; r < num_resources; r++) { 340*2b730f8bSJeremy L Thompson PetscCall(RunWithDM(rp, dm_deg, ceed_resources[r])); 3415f284d84SJed Brown } 342*2b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_deg)); 3439b072555Sjeremylt rp->q_extra = q_extra; 3445f284d84SJed Brown } 3455f284d84SJed Brown 346*2b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm)); 3475f284d84SJed Brown PetscFunctionReturn(0); 3485f284d84SJed Brown } 3495f284d84SJed Brown 3501e284482Svaleriabarra int main(int argc, char **argv) { 351*2b730f8bSJeremy L Thompson PetscInt comm_size; 3521e284482Svaleriabarra RunParams rp; 35353a0f73bSJed Brown MPI_Comm comm; 35453a0f73bSJed Brown char filename[PETSC_MAX_PATH_LEN]; 3559b072555Sjeremylt char *ceed_resources[30]; 3569b072555Sjeremylt PetscInt num_ceed_resources = 30; 35753a0f73bSJed Brown char hostname[PETSC_MAX_PATH_LEN]; 3581e284482Svaleriabarra 3599b072555Sjeremylt PetscInt dim = 3, mesh_elem[3] = {3, 3, 3}; 360*2b730f8bSJeremy L Thompson PetscInt num_degrees = 30, degree[30] = {}, num_local_nodes = 2, local_nodes[2] = {}; 3619b072555Sjeremylt PetscMPIInt ranks_per_node; 362c36f77d8SJed Brown PetscBool degree_set; 3639b072555Sjeremylt BPType bp_choices[10]; 3649b072555Sjeremylt PetscInt num_bp_choices = 10; 3650c59ef15SJeremy L Thompson 3661e284482Svaleriabarra // Initialize PETSc 367*2b730f8bSJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3680c59ef15SJeremy L Thompson comm = PETSC_COMM_WORLD; 369*2b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(comm, &comm_size)); 37053a0f73bSJed Brown #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 37153a0f73bSJed Brown { 37253a0f73bSJed Brown MPI_Comm splitcomm; 373*2b730f8bSJeremy L Thompson PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm)); 374*2b730f8bSJeremy L Thompson PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node)); 375*2b730f8bSJeremy L Thompson PetscCall(MPI_Comm_free(&splitcomm)); 37653a0f73bSJed Brown } 37753a0f73bSJed Brown #else 3789b072555Sjeremylt ranks_per_node = -1; // Unknown 37953a0f73bSJed Brown #endif 380cb32e2e7SValeria Barra 381c36f77d8SJed Brown // Setup all parameters needed in Run() 382*2b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &rp)); 383c36f77d8SJed Brown rp->comm = comm; 384c36f77d8SJed Brown 38532d2ee49SValeria Barra // Read command line options 38667490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 387565a3730SJed Brown { 388565a3730SJed Brown PetscBool set; 389*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set)); 390565a3730SJed Brown if (!set) { 3919b072555Sjeremylt bp_choices[0] = CEED_BP1; 3929b072555Sjeremylt num_bp_choices = 1; 393565a3730SJed Brown } 394565a3730SJed Brown } 395c36f77d8SJed Brown rp->test_mode = PETSC_FALSE; 396*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, rp->test_mode, &rp->test_mode, NULL)); 397c36f77d8SJed Brown rp->write_solution = PETSC_FALSE; 398*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, rp->write_solution, &rp->write_solution, NULL)); 399de1229c5Srezgarshakeri rp->simplex = PETSC_FALSE; 400*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, rp->simplex, &rp->simplex, NULL)); 401eb6a3b92Srezgarshakeri if ((bp_choices[0] == CEED_BP5 || bp_choices[0] == CEED_BP6) && (rp->simplex)) { 402*2b730f8bSJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex"); 403eb6a3b92Srezgarshakeri } 404c36f77d8SJed Brown degree[0] = rp->test_mode ? 3 : 2; 405*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-degree", "Polynomial degree of tensor product basis", NULL, degree, &num_degrees, °ree_set)); 406*2b730f8bSJeremy L Thompson if (!degree_set) num_degrees = 1; 4079b072555Sjeremylt rp->q_extra = PETSC_DECIDE; 408*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points (-1 for auto)", NULL, rp->q_extra, &rp->q_extra, NULL)); 409565a3730SJed Brown { 410565a3730SJed Brown PetscBool set; 411*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsStringArray("-ceed", "CEED resource specifier (comma-separated list)", NULL, ceed_resources, &num_ceed_resources, &set)); 412565a3730SJed Brown if (!set) { 413*2b730f8bSJeremy L Thompson PetscCall(PetscStrallocpy("/cpu/self", &ceed_resources[0])); 4149b072555Sjeremylt num_ceed_resources = 1; 415565a3730SJed Brown } 416565a3730SJed Brown } 417*2b730f8bSJeremy L Thompson PetscCall(PetscGetHostName(hostname, sizeof hostname)); 418*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL)); 419c36f77d8SJed Brown rp->read_mesh = PETSC_FALSE; 420*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &rp->read_mesh)); 421c36f77d8SJed Brown rp->filename = filename; 422c36f77d8SJed Brown if (!rp->read_mesh) { 423cb32e2e7SValeria Barra PetscInt tmp = dim; 424*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL)); 425cb32e2e7SValeria Barra } 4269b072555Sjeremylt local_nodes[0] = 1000; 427*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-local_nodes", 42853a0f73bSJed Brown "Target number of locally owned nodes per " 42953a0f73bSJed Brown "process (single value or min,max)", 430*2b730f8bSJeremy L Thompson NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes)); 431*2b730f8bSJeremy L Thompson if (num_local_nodes < 2) local_nodes[1] = 2 * local_nodes[0]; 432c36f77d8SJed Brown { 433c36f77d8SJed Brown PetscInt two = 2; 434c36f77d8SJed Brown rp->ksp_max_it_clip[0] = 5; 435c36f77d8SJed Brown rp->ksp_max_it_clip[1] = 20; 436*2b730f8bSJeremy 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, 437*2b730f8bSJeremy L Thompson NULL)); 438c36f77d8SJed Brown } 439da9108adSvaleriabarra if (!degree_set) { 4409b072555Sjeremylt PetscInt max_degree = 8; 441*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-max_degree", "Range of degrees [1, max_degree] to run with", NULL, max_degree, &max_degree, NULL)); 442*2b730f8bSJeremy L Thompson for (PetscInt i = 0; i < max_degree; i++) degree[i] = i + 1; 4439b072555Sjeremylt num_degrees = max_degree; 44453a0f73bSJed Brown } 44553a0f73bSJed Brown { 44653a0f73bSJed Brown PetscBool flg; 4479b072555Sjeremylt PetscInt p = ranks_per_node; 448*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg)); 4499b072555Sjeremylt if (flg) ranks_per_node = p; 45053a0f73bSJed Brown } 4519396343dSjeremylt 45267490bc6SJeremy L Thompson PetscOptionsEnd(); 453cb32e2e7SValeria Barra 4541e284482Svaleriabarra // Register PETSc logging stage 455*2b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("Solve Stage", &rp->solve_stage)); 45609a940d7Sjeremylt 45753a0f73bSJed Brown rp->hostname = hostname; 4581e284482Svaleriabarra rp->dim = dim; 4599b072555Sjeremylt rp->mesh_elem = mesh_elem; 4609b072555Sjeremylt rp->ranks_per_node = ranks_per_node; 461cb32e2e7SValeria Barra 46253a0f73bSJed Brown for (PetscInt d = 0; d < num_degrees; d++) { 46353a0f73bSJed Brown PetscInt deg = degree[d]; 4649b072555Sjeremylt for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) { 46553a0f73bSJed Brown rp->degree = deg; 4669b072555Sjeremylt rp->local_nodes = n; 467*2b730f8bSJeremy L Thompson PetscCall(Run(rp, num_ceed_resources, ceed_resources, num_bp_choices, bp_choices)); 468565a3730SJed Brown } 469565a3730SJed Brown } 4701e284482Svaleriabarra // Clear memory 471*2b730f8bSJeremy L Thompson PetscCall(PetscFree(rp)); 4729b072555Sjeremylt for (PetscInt i = 0; i < num_ceed_resources; i++) { 473*2b730f8bSJeremy L Thompson PetscCall(PetscFree(ceed_resources[i])); 474565a3730SJed Brown } 4750c59ef15SJeremy L Thompson return PetscFinalize(); 4760c59ef15SJeremy L Thompson } 477