xref: /libCEED/examples/petsc/bps.c (revision 3c11f1fc4579930f14e261f40b0038a1553ab9e9)
1d275d636SJeremy L Thompson // Copyright (c) 2017-2025, 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 //
10ea61e9acSJeremy 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 //
278ad9a6ffSJames Wright //TESTARGS(name="BP3, tet elements") -ceed {ceed_resource} -test -problem bp3 -degree 3 -ksp_max_it_clip 50,50 -simplex
2896093a69SJeremy L Thompson //TESTARGS(name="BP5, hex elements") -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 18,18
29*3c11f1fcSJeremy L Thompson //TESTARGS(name="BP1+3, hex elements") -ceed {ceed_resource} -test -problem bp1_3 -degree 3 -ksp_max_it_clip 18,18
30*3c11f1fcSJeremy L Thompson //TESTARGS(name="BP2+4, hex elements") -ceed {ceed_resource} -test -problem bp2_4 -degree 3 -ksp_max_it_clip 18,18
310c59ef15SJeremy L Thompson 
320c59ef15SJeremy L Thompson /// @file
33cb32e2e7SValeria Barra /// CEED BPs example using PETSc with DMPlex
34cb32e2e7SValeria Barra /// See bpsraw.c for a "raw" implementation using a structured grid.
35cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc with DMPlex\n";
360c59ef15SJeremy L Thompson 
372b730f8bSJeremy L Thompson #include "bps.h"
382b730f8bSJeremy L Thompson 
39636cccdbSjeremylt #include <ceed.h>
40636cccdbSjeremylt #include <petscdmplex.h>
41636cccdbSjeremylt #include <petscksp.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 
52d5b2ba77SJed Brown // -----------------------------------------------------------------------------
531e284482Svaleriabarra // Main body of program, called in a loop for performance benchmarking purposes
541e284482Svaleriabarra // -----------------------------------------------------------------------------
552b730f8bSJeremy L Thompson static PetscErrorCode RunWithDM(RunParams rp, DM dm, const char *ceed_resource) {
561e284482Svaleriabarra   double               my_rt_start, my_rt, rt_min, rt_max;
579b072555Sjeremylt   PetscInt             xl_size, l_size, g_size;
589b072555Sjeremylt   Vec                  X, X_loc, rhs, rhs_loc;
599b072555Sjeremylt   Mat                  mat_O;
60819eb1b3Sjeremylt   KSP                  ksp;
616c88e6a2Srezgarshakeri   OperatorApplyContext op_apply_ctx, op_error_ctx;
620c59ef15SJeremy L Thompson   Ceed                 ceed;
639b072555Sjeremylt   CeedData             ceed_data;
649b072555Sjeremylt   CeedQFunction        qf_error;
659b072555Sjeremylt   CeedOperator         op_error;
669b072555Sjeremylt   CeedVector           rhs_ceed, target;
678c03e814SJeremy L Thompson   VecType              vec_type = VECSTANDARD;
689b072555Sjeremylt   PetscMemType         mem_type;
691e284482Svaleriabarra 
701e284482Svaleriabarra   PetscFunctionBeginUser;
711e284482Svaleriabarra   // Set up libCEED
729b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
739b072555Sjeremylt   CeedMemType mem_type_backend;
749b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
751e284482Svaleriabarra 
768c03e814SJeremy L Thompson   // Set mesh vec_type
779b072555Sjeremylt   switch (mem_type_backend) {
782b730f8bSJeremy L Thompson     case CEED_MEM_HOST:
792b730f8bSJeremy L Thompson       vec_type = VECSTANDARD;
802b730f8bSJeremy L Thompson       break;
81b68a8d79SJed Brown     case CEED_MEM_DEVICE: {
82b68a8d79SJed Brown       const char *resolved;
838c03e814SJeremy L Thompson 
84b68a8d79SJed Brown       CeedGetResource(ceed, &resolved);
859b072555Sjeremylt       if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
862b730f8bSJeremy L Thompson       else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
879b072555Sjeremylt       else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
889b072555Sjeremylt       else vec_type = VECSTANDARD;
891e284482Svaleriabarra     }
90b68a8d79SJed Brown   }
912b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(dm, vec_type));
92072672edSJeremy L Thompson   PetscCall(DMSetFromOptions(dm));
93b68a8d79SJed Brown 
94b68a8d79SJed Brown   // Create global and local solution vectors
952b730f8bSJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &X));
962b730f8bSJeremy L Thompson   PetscCall(VecGetLocalSize(X, &l_size));
972b730f8bSJeremy L Thompson   PetscCall(VecGetSize(X, &g_size));
982b730f8bSJeremy L Thompson   PetscCall(DMCreateLocalVector(dm, &X_loc));
992b730f8bSJeremy L Thompson   PetscCall(VecGetSize(X_loc, &xl_size));
1002b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X, &rhs));
1011e284482Svaleriabarra 
1021e284482Svaleriabarra   // Operator
1032b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &op_apply_ctx));
1042b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &op_error_ctx));
1052b730f8bSJeremy L Thompson   PetscCall(MatCreateShell(rp->comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
1062b730f8bSJeremy L Thompson   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
1072b730f8bSJeremy L Thompson   PetscCall(MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag));
1082b730f8bSJeremy L Thompson   PetscCall(MatShellSetVecType(mat_O, vec_type));
1091e284482Svaleriabarra 
1101e284482Svaleriabarra   // Print summary
1111e284482Svaleriabarra   if (!rp->test_mode) {
1129b072555Sjeremylt     PetscInt P = rp->degree + 1, Q = P + rp->q_extra;
1131e284482Svaleriabarra 
1149b072555Sjeremylt     const char *used_resource;
1159b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
1161e284482Svaleriabarra 
1179b072555Sjeremylt     VecType vec_type;
1182b730f8bSJeremy L Thompson     PetscCall(VecGetType(X, &vec_type));
1191e284482Svaleriabarra 
1209b072555Sjeremylt     PetscInt c_start, c_end;
1212b730f8bSJeremy L Thompson     PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
122de1229c5Srezgarshakeri     DMPolytopeType cell_type;
1232b730f8bSJeremy L Thompson     PetscCall(DMPlexGetCellType(dm, c_start, &cell_type));
124de1229c5Srezgarshakeri     CeedElemTopology elem_topo = ElemTopologyP2C(cell_type);
12553a0f73bSJed Brown     PetscMPIInt      comm_size;
1262b730f8bSJeremy L Thompson     PetscCall(MPI_Comm_size(rp->comm, &comm_size));
1272b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(rp->comm,
128990fdeb6SJeremy L Thompson                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
12953a0f73bSJed Brown                           "  MPI:\n"
13053a0f73bSJed Brown                           "    Hostname                                : %s\n"
13153a0f73bSJed Brown                           "    Total ranks                             : %d\n"
13253a0f73bSJed Brown                           "    Ranks per compute node                  : %d\n"
1331e284482Svaleriabarra                           "  PETSc:\n"
1341e284482Svaleriabarra                           "    PETSc Vec Type                          : %s\n"
1351e284482Svaleriabarra                           "  libCEED:\n"
1361e284482Svaleriabarra                           "    libCEED Backend                         : %s\n"
1371e284482Svaleriabarra                           "    libCEED Backend MemType                 : %s\n"
1381e284482Svaleriabarra                           "  Mesh:\n"
1394d00b080SJeremy L Thompson                           "    Solution Order (P)                      : %" PetscInt_FMT "\n"
1404d00b080SJeremy L Thompson                           "    Quadrature  Order (Q)                   : %" PetscInt_FMT "\n"
1414d00b080SJeremy L Thompson                           "    Additional quadrature points (q_extra)  : %" PetscInt_FMT "\n"
14208140895SJed Brown                           "    Global nodes                            : %" PetscInt_FMT "\n"
14308140895SJed Brown                           "    Local Elements                          : %" PetscInt_FMT "\n"
14451ad7d5bSrezgarshakeri                           "    Element topology                        : %s\n"
14508140895SJed Brown                           "    Owned nodes                             : %" PetscInt_FMT "\n"
14608140895SJed Brown                           "    DoF per node                            : %" PetscInt_FMT "\n",
1472b730f8bSJeremy L Thompson                           rp->bp_choice + 1, rp->hostname, comm_size, rp->ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P,
1482b730f8bSJeremy 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,
1492b730f8bSJeremy L Thompson                           rp->num_comp_u));
1501e284482Svaleriabarra   }
1511e284482Svaleriabarra 
1521e284482Svaleriabarra   // Create RHS vector
1532b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc, &rhs_loc));
1542b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs_loc));
1559b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
156179e5961SZach Atkins   PetscCall(VecP2C(rhs_loc, &mem_type, rhs_ceed));
1571e284482Svaleriabarra 
1582b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &ceed_data));
1592b730f8bSJeremy 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],
1604dbe2ad5SJeremy L Thompson                                  ceed_data, true, true, rhs_ceed, &target));
1611e284482Svaleriabarra 
1621e284482Svaleriabarra   // Gather RHS
163179e5961SZach Atkins   PetscCall(VecC2P(rhs_ceed, mem_type, rhs_loc));
1642b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs));
1652b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs));
1669b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
1671e284482Svaleriabarra 
1686a6c615bSJeremy L Thompson   // Create the error QFunction
1692b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error, bp_options[rp->bp_choice].error_loc, &qf_error);
1709b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP);
1719b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE);
1722b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE);
17338f32c05Srezgarshakeri   CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_INTERP);
1741e284482Svaleriabarra 
1751e284482Svaleriabarra   // Create the error operator
1762b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
1772b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
178356036faSJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_NONE, target);
179356036faSJeremy L Thompson   CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data);
1802b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
1811e284482Svaleriabarra 
1826c88e6a2Srezgarshakeri   // Set up apply operator context
1832b730f8bSJeremy L Thompson   PetscCall(SetupApplyOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_apply_ctx));
1842b730f8bSJeremy L Thompson   PetscCall(KSPCreate(rp->comm, &ksp));
1851e284482Svaleriabarra   {
1861e284482Svaleriabarra     PC pc;
1872b730f8bSJeremy L Thompson     PetscCall(KSPGetPC(ksp, &pc));
188*3c11f1fcSJeremy L Thompson     if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2 || rp->bp_choice == CEED_BP13 || rp->bp_choice == CEED_BP24) {
1892b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCJACOBI));
190*3c11f1fcSJeremy L Thompson       if (rp->simplex || rp->bp_choice == CEED_BP13 || rp->bp_choice == CEED_BP24) {
1912b730f8bSJeremy L Thompson         PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
192345f5e33Srezgarshakeri       } else {
1932b730f8bSJeremy L Thompson         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
194345f5e33Srezgarshakeri       }
1951e284482Svaleriabarra     } else {
1962b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCNONE));
1971e284482Svaleriabarra     }
1982b730f8bSJeremy L Thompson     PetscCall(KSPSetType(ksp, KSPCG));
1992b730f8bSJeremy L Thompson     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
2002b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
2011e284482Svaleriabarra   }
2022b730f8bSJeremy L Thompson   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
2031e284482Svaleriabarra 
204da9108adSvaleriabarra   // First run's performance log is not considered for benchmarking purposes
2052b730f8bSJeremy L Thompson   PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
2061e284482Svaleriabarra   my_rt_start = MPI_Wtime();
2072b730f8bSJeremy L Thompson   PetscCall(KSPSolve(ksp, rhs, X));
2081e284482Svaleriabarra   my_rt = MPI_Wtime() - my_rt_start;
2092b730f8bSJeremy L Thompson   PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm));
2101e284482Svaleriabarra   // Set maxits based on first iteration timing
2111e284482Svaleriabarra   if (my_rt > 0.02) {
2122b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[0]));
2131e284482Svaleriabarra   } else {
2142b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[1]));
2151e284482Svaleriabarra   }
2162b730f8bSJeremy L Thompson   PetscCall(KSPSetFromOptions(ksp));
2171e284482Svaleriabarra 
2181e284482Svaleriabarra   // Timed solve
2192b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(X));
2202b730f8bSJeremy L Thompson   PetscCall(PetscBarrier((PetscObject)ksp));
2211e284482Svaleriabarra 
2221e284482Svaleriabarra   // -- Performance logging
2232b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(rp->solve_stage));
2241e284482Svaleriabarra 
2251e284482Svaleriabarra   // -- Solve
2261e284482Svaleriabarra   my_rt_start = MPI_Wtime();
2272b730f8bSJeremy L Thompson   PetscCall(KSPSolve(ksp, rhs, X));
2281e284482Svaleriabarra   my_rt = MPI_Wtime() - my_rt_start;
2291e284482Svaleriabarra 
2301e284482Svaleriabarra   // -- Performance logging
2312b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
2321e284482Svaleriabarra 
2331e284482Svaleriabarra   // Output results
2341e284482Svaleriabarra   {
2359b072555Sjeremylt     KSPType            ksp_type;
2361e284482Svaleriabarra     KSPConvergedReason reason;
2371e284482Svaleriabarra     PetscReal          rnorm;
2381e284482Svaleriabarra     PetscInt           its;
2392b730f8bSJeremy L Thompson     PetscCall(KSPGetType(ksp, &ksp_type));
2402b730f8bSJeremy L Thompson     PetscCall(KSPGetConvergedReason(ksp, &reason));
2412b730f8bSJeremy L Thompson     PetscCall(KSPGetIterationNumber(ksp, &its));
2422b730f8bSJeremy L Thompson     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
2431e284482Svaleriabarra     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
2442b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(rp->comm,
2451e284482Svaleriabarra                             "  KSP:\n"
2461e284482Svaleriabarra                             "    KSP Type                                : %s\n"
2471e284482Svaleriabarra                             "    KSP Convergence                         : %s\n"
248a9b2c5ddSrezgarshakeri                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
2491e284482Svaleriabarra                             "    Final rnorm                             : %e\n",
2502b730f8bSJeremy L Thompson                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
2511e284482Svaleriabarra     }
2521e284482Svaleriabarra     if (!rp->test_mode) {
2532b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(rp->comm, "  Performance:\n"));
2541e284482Svaleriabarra     }
2551e284482Svaleriabarra     {
2566c88e6a2Srezgarshakeri       // Set up error operator context
2572b730f8bSJeremy L Thompson       PetscCall(SetupErrorOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx));
25838f32c05Srezgarshakeri       PetscScalar l2_error;
2592b730f8bSJeremy L Thompson       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
260*3c11f1fcSJeremy L Thompson       // Tighter tol for BP1, BP2
261*3c11f1fcSJeremy L Thompson       // Looser tol for BP3, BP4, BP5, and BP6 with extra for vector valued problems
262*3c11f1fcSJeremy L Thompson       // BP1+3 and BP2+4 follow the pattern for BP3 and BP4
263*3c11f1fcSJeremy L Thompson       PetscReal tol = rp->bp_choice < CEED_BP3 ? 5e-4 : (5e-2 + (rp->bp_choice % 2 == 1 ? 5e-3 : 0));
26438f32c05Srezgarshakeri       if (!rp->test_mode || l2_error > tol) {
2652b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm));
2662b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm));
2672b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(rp->comm,
26838f32c05Srezgarshakeri                               "    L2 Error                                : %e\n"
2691e284482Svaleriabarra                               "    CG Solve Time                           : %g (%g) sec\n",
2702b730f8bSJeremy L Thompson                               (double)l2_error, rt_max, rt_min));
2711e284482Svaleriabarra       }
2721e284482Svaleriabarra     }
2734c583f1fSvaleriabarra     if (!rp->test_mode) {
2742b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(rp->comm, "    DoFs/Sec in CG                          : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
2752b730f8bSJeremy L Thompson                             1e-6 * g_size * its / rt_min));
2761e284482Svaleriabarra     }
2771e284482Svaleriabarra   }
2781e284482Svaleriabarra 
2791e284482Svaleriabarra   if (rp->write_solution) {
2809b072555Sjeremylt     PetscViewer vtk_viewer_soln;
2811e284482Svaleriabarra 
2822b730f8bSJeremy L Thompson     PetscCall(PetscViewerCreate(rp->comm, &vtk_viewer_soln));
2832b730f8bSJeremy L Thompson     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
2842b730f8bSJeremy L Thompson     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
2852b730f8bSJeremy L Thompson     PetscCall(VecView(X, vtk_viewer_soln));
2862b730f8bSJeremy L Thompson     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
2871e284482Svaleriabarra   }
2881e284482Svaleriabarra 
2891e284482Svaleriabarra   // Cleanup
2902b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&X));
2912b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&X_loc));
2922b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
2932b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
2942b730f8bSJeremy L Thompson   PetscCall(MatDestroy(&mat_O));
2952b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_apply_ctx));
2962b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_error_ctx));
2972b730f8bSJeremy L Thompson   PetscCall(CeedDataDestroy(0, ceed_data));
2981e284482Svaleriabarra 
2992b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs));
3002b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs_loc));
3012b730f8bSJeremy L Thompson   PetscCall(KSPDestroy(&ksp));
3021e284482Svaleriabarra   CeedVectorDestroy(&target);
3039b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
3049b072555Sjeremylt   CeedOperatorDestroy(&op_error);
3051e284482Svaleriabarra   CeedDestroy(&ceed);
306ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3071e284482Svaleriabarra }
3081e284482Svaleriabarra 
3092b730f8bSJeremy L Thompson static PetscErrorCode Run(RunParams rp, PetscInt num_resources, char *const *ceed_resources, PetscInt num_bp_choices, const BPType *bp_choices) {
3105f284d84SJed Brown   DM dm;
3115f284d84SJed Brown 
3125f284d84SJed Brown   PetscFunctionBeginUser;
3135f284d84SJed Brown   // Setup DM
3142b730f8bSJeremy L Thompson   PetscCall(CreateDistributedDM(rp, &dm));
3155f284d84SJed Brown 
3169b072555Sjeremylt   for (PetscInt b = 0; b < num_bp_choices; b++) {
3175f284d84SJed Brown     DM       dm_deg;
3189b072555Sjeremylt     VecType  vec_type;
3199b072555Sjeremylt     PetscInt q_extra = rp->q_extra;
3209b072555Sjeremylt     rp->bp_choice    = bp_choices[b];
3219b072555Sjeremylt     rp->num_comp_u   = bp_options[rp->bp_choice].num_comp_u;
3229b072555Sjeremylt     rp->q_extra      = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra;
3232b730f8bSJeremy L Thompson     PetscCall(DMClone(dm, &dm_deg));
3242b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
3252b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm_deg, vec_type));
3265f284d84SJed Brown     // Create DM
327e83e87a5Sjeremylt     PetscInt dim;
3282b730f8bSJeremy L Thompson     PetscCall(DMGetDimension(dm_deg, &dim));
3292b730f8bSJeremy L Thompson     PetscCall(SetupDMByDegree(dm_deg, rp->degree, rp->q_extra, rp->num_comp_u, dim, bp_options[rp->bp_choice].enforce_bc));
3305f284d84SJed Brown     for (PetscInt r = 0; r < num_resources; r++) {
3312b730f8bSJeremy L Thompson       PetscCall(RunWithDM(rp, dm_deg, ceed_resources[r]));
3325f284d84SJed Brown     }
3332b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&dm_deg));
3349b072555Sjeremylt     rp->q_extra = q_extra;
3355f284d84SJed Brown   }
3365f284d84SJed Brown 
3372b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm));
338ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3395f284d84SJed Brown }
3405f284d84SJed Brown 
3411e284482Svaleriabarra int main(int argc, char **argv) {
3424d00b080SJeremy L Thompson   PetscMPIInt comm_size;
3431e284482Svaleriabarra   RunParams   rp;
34453a0f73bSJed Brown   MPI_Comm    comm;
34553a0f73bSJed Brown   char        filename[PETSC_MAX_PATH_LEN];
3469b072555Sjeremylt   char       *ceed_resources[30];
3479b072555Sjeremylt   PetscInt    num_ceed_resources = 30;
34853a0f73bSJed Brown   char        hostname[PETSC_MAX_PATH_LEN];
3491e284482Svaleriabarra 
3509b072555Sjeremylt   PetscInt    dim = 3, mesh_elem[3] = {3, 3, 3};
351a8d73ac5SSebastian Grimberg   PetscInt    num_degrees = 30, degree[30] = {0}, num_local_nodes = 2, local_nodes[2] = {0};
3529b072555Sjeremylt   PetscMPIInt ranks_per_node;
353c36f77d8SJed Brown   PetscBool   degree_set;
3549b072555Sjeremylt   BPType      bp_choices[10];
3559b072555Sjeremylt   PetscInt    num_bp_choices = 10;
3560c59ef15SJeremy L Thompson 
3571e284482Svaleriabarra   // Initialize PETSc
3582b730f8bSJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3590c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
3602b730f8bSJeremy L Thompson   PetscCall(MPI_Comm_size(comm, &comm_size));
36153a0f73bSJed Brown #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
36253a0f73bSJed Brown   {
36353a0f73bSJed Brown     MPI_Comm splitcomm;
3644d00b080SJeremy L Thompson 
3652b730f8bSJeremy L Thompson     PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm));
3662b730f8bSJeremy L Thompson     PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node));
3672b730f8bSJeremy L Thompson     PetscCall(MPI_Comm_free(&splitcomm));
36853a0f73bSJed Brown   }
36953a0f73bSJed Brown #else
3709b072555Sjeremylt   ranks_per_node = -1;  // Unknown
37153a0f73bSJed Brown #endif
372cb32e2e7SValeria Barra 
373c36f77d8SJed Brown   // Setup all parameters needed in Run()
3742b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &rp));
375c36f77d8SJed Brown   rp->comm = comm;
376c36f77d8SJed Brown 
37732d2ee49SValeria Barra   // Read command line options
37867490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
379565a3730SJed Brown   {
380565a3730SJed Brown     PetscBool set;
3812b730f8bSJeremy L Thompson     PetscCall(PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set));
382565a3730SJed Brown     if (!set) {
3839b072555Sjeremylt       bp_choices[0]  = CEED_BP1;
3849b072555Sjeremylt       num_bp_choices = 1;
385565a3730SJed Brown     }
386565a3730SJed Brown   }
387c36f77d8SJed Brown   rp->test_mode = PETSC_FALSE;
3882b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, rp->test_mode, &rp->test_mode, NULL));
389c36f77d8SJed Brown   rp->write_solution = PETSC_FALSE;
3902b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, rp->write_solution, &rp->write_solution, NULL));
391de1229c5Srezgarshakeri   rp->simplex = PETSC_FALSE;
3922b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, rp->simplex, &rp->simplex, NULL));
393eb6a3b92Srezgarshakeri   if ((bp_choices[0] == CEED_BP5 || bp_choices[0] == CEED_BP6) && (rp->simplex)) {
3942b730f8bSJeremy L Thompson     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex");
395eb6a3b92Srezgarshakeri   }
396c36f77d8SJed Brown   degree[0] = rp->test_mode ? 3 : 2;
3972b730f8bSJeremy L Thompson   PetscCall(PetscOptionsIntArray("-degree", "Polynomial degree of tensor product basis", NULL, degree, &num_degrees, &degree_set));
3982b730f8bSJeremy L Thompson   if (!degree_set) num_degrees = 1;
3999b072555Sjeremylt   rp->q_extra = PETSC_DECIDE;
4002b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points (-1 for auto)", NULL, rp->q_extra, &rp->q_extra, NULL));
401565a3730SJed Brown   {
402565a3730SJed Brown     PetscBool set;
4032b730f8bSJeremy L Thompson     PetscCall(PetscOptionsStringArray("-ceed", "CEED resource specifier (comma-separated list)", NULL, ceed_resources, &num_ceed_resources, &set));
404565a3730SJed Brown     if (!set) {
4052b730f8bSJeremy L Thompson       PetscCall(PetscStrallocpy("/cpu/self", &ceed_resources[0]));
4069b072555Sjeremylt       num_ceed_resources = 1;
407565a3730SJed Brown     }
408565a3730SJed Brown   }
4092b730f8bSJeremy L Thompson   PetscCall(PetscGetHostName(hostname, sizeof hostname));
4102b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL));
411c36f77d8SJed Brown   rp->read_mesh = PETSC_FALSE;
4122b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &rp->read_mesh));
413c36f77d8SJed Brown   rp->filename = filename;
414c36f77d8SJed Brown   if (!rp->read_mesh) {
415cb32e2e7SValeria Barra     PetscInt tmp = dim;
4162b730f8bSJeremy L Thompson     PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL));
417cb32e2e7SValeria Barra   }
4189b072555Sjeremylt   local_nodes[0] = 1000;
4192b730f8bSJeremy L Thompson   PetscCall(PetscOptionsIntArray("-local_nodes",
42053a0f73bSJed Brown                                  "Target number of locally owned nodes per "
42153a0f73bSJed Brown                                  "process (single value or min,max)",
4222b730f8bSJeremy L Thompson                                  NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes));
4232b730f8bSJeremy L Thompson   if (num_local_nodes < 2) local_nodes[1] = 2 * local_nodes[0];
424c36f77d8SJed Brown   {
425c36f77d8SJed Brown     PetscInt two           = 2;
426c36f77d8SJed Brown     rp->ksp_max_it_clip[0] = 5;
427c36f77d8SJed Brown     rp->ksp_max_it_clip[1] = 20;
4282b730f8bSJeremy 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,
4292b730f8bSJeremy L Thompson                                    NULL));
430c36f77d8SJed Brown   }
431da9108adSvaleriabarra   if (!degree_set) {
4329b072555Sjeremylt     PetscInt max_degree = 8;
4332b730f8bSJeremy L Thompson     PetscCall(PetscOptionsInt("-max_degree", "Range of degrees [1, max_degree] to run with", NULL, max_degree, &max_degree, NULL));
4342b730f8bSJeremy L Thompson     for (PetscInt i = 0; i < max_degree; i++) degree[i] = i + 1;
4359b072555Sjeremylt     num_degrees = max_degree;
43653a0f73bSJed Brown   }
43753a0f73bSJed Brown   {
43853a0f73bSJed Brown     PetscBool flg;
4399b072555Sjeremylt     PetscInt  p = ranks_per_node;
4402b730f8bSJeremy L Thompson     PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg));
4419b072555Sjeremylt     if (flg) ranks_per_node = p;
44253a0f73bSJed Brown   }
4439396343dSjeremylt 
44467490bc6SJeremy L Thompson   PetscOptionsEnd();
445cb32e2e7SValeria Barra 
4461e284482Svaleriabarra   // Register PETSc logging stage
4472b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("Solve Stage", &rp->solve_stage));
44809a940d7Sjeremylt 
44953a0f73bSJed Brown   rp->hostname       = hostname;
4501e284482Svaleriabarra   rp->dim            = dim;
4519b072555Sjeremylt   rp->mesh_elem      = mesh_elem;
4529b072555Sjeremylt   rp->ranks_per_node = ranks_per_node;
453cb32e2e7SValeria Barra 
45453a0f73bSJed Brown   for (PetscInt d = 0; d < num_degrees; d++) {
45553a0f73bSJed Brown     PetscInt deg = degree[d];
4569b072555Sjeremylt     for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) {
45753a0f73bSJed Brown       rp->degree      = deg;
4589b072555Sjeremylt       rp->local_nodes = n;
4592b730f8bSJeremy L Thompson       PetscCall(Run(rp, num_ceed_resources, ceed_resources, num_bp_choices, bp_choices));
460565a3730SJed Brown     }
461565a3730SJed Brown   }
4621e284482Svaleriabarra   // Clear memory
4632b730f8bSJeremy L Thompson   PetscCall(PetscFree(rp));
4649b072555Sjeremylt   for (PetscInt i = 0; i < num_ceed_resources; i++) {
4652b730f8bSJeremy L Thompson     PetscCall(PetscFree(ceed_resources[i]));
466565a3730SJed Brown   }
4670c59ef15SJeremy L Thompson   return PetscFinalize();
4680c59ef15SJeremy L Thompson }
469