xref: /libCEED/examples/petsc/bps.c (revision 9b072555b57804a6f4e0fc2b1ad83be89838f0e5)
1819eb1b3Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2819eb1b3Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3819eb1b3Sjeremylt // reserved. See files LICENSE and NOTICE for details.
4819eb1b3Sjeremylt //
5819eb1b3Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
6819eb1b3Sjeremylt // libraries and APIs for efficient high-order finite element and spectral
7819eb1b3Sjeremylt // element discretizations for exascale applications. For more information and
8819eb1b3Sjeremylt // source code availability see http://github.com/ceed.
9819eb1b3Sjeremylt //
10819eb1b3Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11819eb1b3Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
12819eb1b3Sjeremylt // of Science and the National Nuclear Security Administration) responsible for
13819eb1b3Sjeremylt // the planning and preparation of a capable exascale ecosystem, including
14819eb1b3Sjeremylt // software, applications, hardware, advanced system engineering and early
15819eb1b3Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
16819eb1b3Sjeremylt 
170c59ef15SJeremy L Thompson //                        libCEED + PETSc Example: CEED BPs
180c59ef15SJeremy L Thompson //
190c59ef15SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the
200c59ef15SJeremy L Thompson // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
210c59ef15SJeremy L Thompson //
22cb32e2e7SValeria Barra // The code uses higher level communication protocols in DMPlex.
230c59ef15SJeremy L Thompson //
240c59ef15SJeremy L Thompson // Build with:
250c59ef15SJeremy L Thompson //
260c59ef15SJeremy L Thompson //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
270c59ef15SJeremy L Thompson //
280c59ef15SJeremy L Thompson // Sample runs:
290c59ef15SJeremy L Thompson //
3032d2ee49SValeria Barra //     ./bps -problem bp1 -degree 3
3128688798Sjeremylt //     ./bps -problem bp2 -degree 3
3228688798Sjeremylt //     ./bps -problem bp3 -degree 3
3328688798Sjeremylt //     ./bps -problem bp4 -degree 3
3428688798Sjeremylt //     ./bps -problem bp5 -degree 3 -ceed /cpu/self
3528688798Sjeremylt //     ./bps -problem bp6 -degree 3 -ceed /gpu/cuda
360c59ef15SJeremy L Thompson //
372fbc6e41SJeremy L Thompson //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 15,15
380c59ef15SJeremy L Thompson 
390c59ef15SJeremy L Thompson /// @file
40cb32e2e7SValeria Barra /// CEED BPs example using PETSc with DMPlex
41cb32e2e7SValeria Barra /// See bpsraw.c for a "raw" implementation using a structured grid.
42cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc with DMPlex\n";
430c59ef15SJeremy L Thompson 
443d576824SJeremy L Thompson #include <ceed.h>
453d576824SJeremy L Thompson #include <petscdmplex.h>
463d576824SJeremy L Thompson #include <petscksp.h>
470c59ef15SJeremy L Thompson #include <stdbool.h>
480c59ef15SJeremy L Thompson #include <string.h>
49e83e87a5Sjeremylt #include "bps.h"
500c59ef15SJeremy L Thompson 
51d5b2ba77SJed Brown // -----------------------------------------------------------------------------
52d5b2ba77SJed Brown // Utilities
53d5b2ba77SJed Brown // -----------------------------------------------------------------------------
54d5b2ba77SJed Brown 
55d5b2ba77SJed Brown // Utility function, compute three factors of an integer
56d5b2ba77SJed Brown static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
57*9b072555Sjeremylt   for (PetscInt d=0, size_left=size; d<3; d++) {
58*9b072555Sjeremylt     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d)));
59*9b072555Sjeremylt     while (try * (size_left / try) != size_left) try++;
60d5b2ba77SJed Brown     m[reverse ? 2-d : d] = try;
61*9b072555Sjeremylt     size_left /= try;
62d5b2ba77SJed Brown   }
63d5b2ba77SJed Brown }
64d5b2ba77SJed Brown 
65d5b2ba77SJed Brown static int Max3(const PetscInt a[3]) {
66d5b2ba77SJed Brown   return PetscMax(a[0], PetscMax(a[1], a[2]));
67d5b2ba77SJed Brown }
68d5b2ba77SJed Brown 
69d5b2ba77SJed Brown static int Min3(const PetscInt a[3]) {
70d5b2ba77SJed Brown   return PetscMin(a[0], PetscMin(a[1], a[2]));
71d5b2ba77SJed Brown }
72d5b2ba77SJed Brown 
731e284482Svaleriabarra // -----------------------------------------------------------------------------
741e284482Svaleriabarra // Parameter structure for running problems
751e284482Svaleriabarra // -----------------------------------------------------------------------------
761e284482Svaleriabarra typedef struct RunParams_ *RunParams;
771e284482Svaleriabarra struct RunParams_ {
780c59ef15SJeremy L Thompson   MPI_Comm comm;
79*9b072555Sjeremylt   PetscBool test_mode, read_mesh, user_l_nodes, write_solution;
805f284d84SJed Brown   char *filename, *hostname;
81*9b072555Sjeremylt   PetscInt local_nodes, degree, q_extra, dim, num_comp_u, *mesh_elem;
82c36f77d8SJed Brown   PetscInt ksp_max_it_clip[2];
83*9b072555Sjeremylt   PetscMPIInt ranks_per_node;
84*9b072555Sjeremylt   BPType bp_choice;
85*9b072555Sjeremylt   PetscLogStage solve_stage;
861e284482Svaleriabarra };
871e284482Svaleriabarra 
881e284482Svaleriabarra // -----------------------------------------------------------------------------
891e284482Svaleriabarra // Main body of program, called in a loop for performance benchmarking purposes
901e284482Svaleriabarra // -----------------------------------------------------------------------------
915f284d84SJed Brown static PetscErrorCode RunWithDM(RunParams rp, DM dm,
92*9b072555Sjeremylt                                 const char *ceed_resource) {
935f284d84SJed Brown   PetscErrorCode ierr;
941e284482Svaleriabarra   double my_rt_start, my_rt, rt_min, rt_max;
95*9b072555Sjeremylt   PetscInt xl_size, l_size, g_size;
961e284482Svaleriabarra   PetscScalar *r;
97*9b072555Sjeremylt   Vec X, X_loc, rhs, rhs_loc;
98*9b072555Sjeremylt   Mat mat_O;
99819eb1b3Sjeremylt   KSP ksp;
100*9b072555Sjeremylt   UserO user_O;
1010c59ef15SJeremy L Thompson   Ceed ceed;
102*9b072555Sjeremylt   CeedData ceed_data;
103*9b072555Sjeremylt   CeedQFunction qf_error;
104*9b072555Sjeremylt   CeedOperator op_error;
105*9b072555Sjeremylt   CeedVector rhs_ceed, target;
106*9b072555Sjeremylt   VecType vec_type;
107*9b072555Sjeremylt   PetscMemType mem_type;
1081e284482Svaleriabarra 
1091e284482Svaleriabarra   PetscFunctionBeginUser;
1101e284482Svaleriabarra   // Set up libCEED
111*9b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
112*9b072555Sjeremylt   CeedMemType mem_type_backend;
113*9b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1141e284482Svaleriabarra 
115*9b072555Sjeremylt   ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
116*9b072555Sjeremylt   if (!vec_type) { // Not yet set by user -dm_vec_type
117*9b072555Sjeremylt     switch (mem_type_backend) {
118*9b072555Sjeremylt     case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
119b68a8d79SJed Brown     case CEED_MEM_DEVICE: {
120b68a8d79SJed Brown       const char *resolved;
121b68a8d79SJed Brown       CeedGetResource(ceed, &resolved);
122*9b072555Sjeremylt       if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
123b68a8d79SJed Brown       else if (strstr(resolved, "/gpu/hip/occa"))
124*9b072555Sjeremylt         vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
125*9b072555Sjeremylt       else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
126*9b072555Sjeremylt       else vec_type = VECSTANDARD;
1271e284482Svaleriabarra     }
128b68a8d79SJed Brown     }
129*9b072555Sjeremylt     ierr = DMSetVecType(dm, vec_type); CHKERRQ(ierr);
130b68a8d79SJed Brown   }
131b68a8d79SJed Brown 
132b68a8d79SJed Brown   // Create global and local solution vectors
1331e284482Svaleriabarra   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
134*9b072555Sjeremylt   ierr = VecGetLocalSize(X, &l_size); CHKERRQ(ierr);
135*9b072555Sjeremylt   ierr = VecGetSize(X, &g_size); CHKERRQ(ierr);
136*9b072555Sjeremylt   ierr = DMCreateLocalVector(dm, &X_loc); CHKERRQ(ierr);
137*9b072555Sjeremylt   ierr = VecGetSize(X_loc, &xl_size); CHKERRQ(ierr);
1381e284482Svaleriabarra   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
1391e284482Svaleriabarra 
1401e284482Svaleriabarra   // Operator
141*9b072555Sjeremylt   ierr = PetscMalloc1(1, &user_O); CHKERRQ(ierr);
142*9b072555Sjeremylt   ierr = MatCreateShell(rp->comm, l_size, l_size, g_size, g_size,
143*9b072555Sjeremylt                         user_O, &mat_O); CHKERRQ(ierr);
144*9b072555Sjeremylt   ierr = MatShellSetOperation(mat_O, MATOP_MULT,
1451e284482Svaleriabarra                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
146*9b072555Sjeremylt   ierr = MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL,
1471e284482Svaleriabarra                               (void(*)(void))MatGetDiag); CHKERRQ(ierr);
148*9b072555Sjeremylt   ierr = MatShellSetVecType(mat_O, vec_type); CHKERRQ(ierr);
1491e284482Svaleriabarra 
1501e284482Svaleriabarra   // Print summary
1511e284482Svaleriabarra   if (!rp->test_mode) {
152*9b072555Sjeremylt     PetscInt P = rp->degree + 1, Q = P + rp->q_extra;
1531e284482Svaleriabarra 
154*9b072555Sjeremylt     const char *used_resource;
155*9b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
1561e284482Svaleriabarra 
157*9b072555Sjeremylt     VecType vec_type;
158*9b072555Sjeremylt     ierr = VecGetType(X, &vec_type); CHKERRQ(ierr);
1591e284482Svaleriabarra 
160*9b072555Sjeremylt     PetscInt c_start, c_end;
161*9b072555Sjeremylt     ierr = DMPlexGetHeightStratum(dm, 0, &c_start, &c_end); CHKERRQ(ierr);
16253a0f73bSJed Brown     PetscMPIInt comm_size;
16353a0f73bSJed Brown     ierr = MPI_Comm_size(rp->comm, &comm_size); CHKERRQ(ierr);
1641e284482Svaleriabarra     ierr = PetscPrintf(rp->comm,
1651e284482Svaleriabarra                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
16653a0f73bSJed Brown                        "  MPI:\n"
16753a0f73bSJed Brown                        "    Hostname                           : %s\n"
16853a0f73bSJed Brown                        "    Total ranks                        : %d\n"
16953a0f73bSJed Brown                        "    Ranks per compute node             : %d\n"
1701e284482Svaleriabarra                        "  PETSc:\n"
1711e284482Svaleriabarra                        "    PETSc Vec Type                     : %s\n"
1721e284482Svaleriabarra                        "  libCEED:\n"
1731e284482Svaleriabarra                        "    libCEED Backend                    : %s\n"
1741e284482Svaleriabarra                        "    libCEED Backend MemType            : %s\n"
1751e284482Svaleriabarra                        "  Mesh:\n"
1761e284482Svaleriabarra                        "    Number of 1D Basis Nodes (P)       : %d\n"
1771e284482Svaleriabarra                        "    Number of 1D Quadrature Points (Q) : %d\n"
1781e284482Svaleriabarra                        "    Global nodes                       : %D\n"
1791e284482Svaleriabarra                        "    Local Elements                     : %D\n"
1801e284482Svaleriabarra                        "    Owned nodes                        : %D\n"
1811e284482Svaleriabarra                        "    DoF per node                       : %D\n",
182*9b072555Sjeremylt                        rp->bp_choice+1, rp->hostname, comm_size,
183*9b072555Sjeremylt                        rp->ranks_per_node, vec_type, used_resource,
184*9b072555Sjeremylt                        CeedMemTypes[mem_type_backend],
185*9b072555Sjeremylt                        P, Q, g_size/rp->num_comp_u, c_end - c_start, l_size/rp->num_comp_u,
186*9b072555Sjeremylt                        rp->num_comp_u);
1871e284482Svaleriabarra     CHKERRQ(ierr);
1881e284482Svaleriabarra   }
1891e284482Svaleriabarra 
1901e284482Svaleriabarra   // Create RHS vector
191*9b072555Sjeremylt   ierr = VecDuplicate(X_loc, &rhs_loc); CHKERRQ(ierr);
192*9b072555Sjeremylt   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
193*9b072555Sjeremylt   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
194*9b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
195*9b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
1961e284482Svaleriabarra 
197*9b072555Sjeremylt   ierr = PetscMalloc1(1, &ceed_data); CHKERRQ(ierr);
198*9b072555Sjeremylt   ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->q_extra,
199*9b072555Sjeremylt                               rp->dim, rp->num_comp_u, g_size, xl_size, bp_options[rp->bp_choice],
200*9b072555Sjeremylt                               ceed_data, true, rhs_ceed, &target); CHKERRQ(ierr);
2011e284482Svaleriabarra 
2021e284482Svaleriabarra   // Gather RHS
203*9b072555Sjeremylt   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
204*9b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
2051e284482Svaleriabarra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
206*9b072555Sjeremylt   ierr = DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr);
207*9b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
2081e284482Svaleriabarra 
2096a6c615bSJeremy L Thompson   // Create the error QFunction
210*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error,
211*9b072555Sjeremylt                               bp_options[rp->bp_choice].error_loc, &qf_error);
212*9b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP);
213*9b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE);
214*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_NONE);
2151e284482Svaleriabarra 
2161e284482Svaleriabarra   // Create the error operator
217*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
218*9b072555Sjeremylt                      &op_error);
219*9b072555Sjeremylt   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u,
220*9b072555Sjeremylt                        ceed_data->basis_u, CEED_VECTOR_ACTIVE);
221*9b072555Sjeremylt   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i,
2221e284482Svaleriabarra                        CEED_BASIS_COLLOCATED, target);
223*9b072555Sjeremylt   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u_i,
2241e284482Svaleriabarra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
2251e284482Svaleriabarra 
2261e284482Svaleriabarra   // Set up Mat
227*9b072555Sjeremylt   user_O->comm = rp->comm;
228*9b072555Sjeremylt   user_O->dm = dm;
229*9b072555Sjeremylt   user_O->X_loc = X_loc;
230*9b072555Sjeremylt   ierr = VecDuplicate(X_loc, &user_O->Y_loc); CHKERRQ(ierr);
231*9b072555Sjeremylt   user_O->x_ceed = ceed_data->x_ceed;
232*9b072555Sjeremylt   user_O->y_ceed = ceed_data->y_ceed;
233*9b072555Sjeremylt   user_O->op = ceed_data->op_apply;
234*9b072555Sjeremylt   user_O->ceed = ceed;
2351e284482Svaleriabarra 
2361e284482Svaleriabarra   ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr);
2371e284482Svaleriabarra   {
2381e284482Svaleriabarra     PC pc;
2391e284482Svaleriabarra     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
240*9b072555Sjeremylt     if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) {
2411e284482Svaleriabarra       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
2421e284482Svaleriabarra       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
2431e284482Svaleriabarra     } else {
2441e284482Svaleriabarra       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
2451e284482Svaleriabarra     }
2461e284482Svaleriabarra     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
2471e284482Svaleriabarra     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
2481e284482Svaleriabarra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
2491e284482Svaleriabarra                             PETSC_DEFAULT); CHKERRQ(ierr);
2501e284482Svaleriabarra   }
251*9b072555Sjeremylt   ierr = KSPSetOperators(ksp, mat_O, mat_O); CHKERRQ(ierr);
2521e284482Svaleriabarra 
253da9108adSvaleriabarra   // First run's performance log is not considered for benchmarking purposes
2541e284482Svaleriabarra   ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
2551e284482Svaleriabarra   CHKERRQ(ierr);
2561e284482Svaleriabarra   my_rt_start = MPI_Wtime();
2571e284482Svaleriabarra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
2581e284482Svaleriabarra   my_rt = MPI_Wtime() - my_rt_start;
2591e284482Svaleriabarra   ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
2601e284482Svaleriabarra   CHKERRQ(ierr);
2611e284482Svaleriabarra   // Set maxits based on first iteration timing
2621e284482Svaleriabarra   if (my_rt > 0.02) {
263565a3730SJed Brown     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
264565a3730SJed Brown                             rp->ksp_max_it_clip[0]);
2651e284482Svaleriabarra     CHKERRQ(ierr);
2661e284482Svaleriabarra   } else {
267565a3730SJed Brown     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
268565a3730SJed Brown                             rp->ksp_max_it_clip[1]);
2691e284482Svaleriabarra     CHKERRQ(ierr);
2701e284482Svaleriabarra   }
2711e284482Svaleriabarra   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
2721e284482Svaleriabarra 
2731e284482Svaleriabarra   // Timed solve
2741e284482Svaleriabarra   ierr = VecZeroEntries(X); CHKERRQ(ierr);
2751e284482Svaleriabarra   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
2761e284482Svaleriabarra 
2771e284482Svaleriabarra   // -- Performance logging
278*9b072555Sjeremylt   ierr = PetscLogStagePush(rp->solve_stage); CHKERRQ(ierr);
2791e284482Svaleriabarra 
2801e284482Svaleriabarra   // -- Solve
2811e284482Svaleriabarra   my_rt_start = MPI_Wtime();
2821e284482Svaleriabarra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
2831e284482Svaleriabarra   my_rt = MPI_Wtime() - my_rt_start;
2841e284482Svaleriabarra 
2851e284482Svaleriabarra   // -- Performance logging
2861e284482Svaleriabarra   ierr = PetscLogStagePop();
2871e284482Svaleriabarra 
2881e284482Svaleriabarra   // Output results
2891e284482Svaleriabarra   {
290*9b072555Sjeremylt     KSPType ksp_type;
2911e284482Svaleriabarra     KSPConvergedReason reason;
2921e284482Svaleriabarra     PetscReal rnorm;
2931e284482Svaleriabarra     PetscInt its;
294*9b072555Sjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
2951e284482Svaleriabarra     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
2961e284482Svaleriabarra     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
2971e284482Svaleriabarra     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
2981e284482Svaleriabarra     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
2991e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,
3001e284482Svaleriabarra                          "  KSP:\n"
3011e284482Svaleriabarra                          "    KSP Type                           : %s\n"
3021e284482Svaleriabarra                          "    KSP Convergence                    : %s\n"
3031e284482Svaleriabarra                          "    Total KSP Iterations               : %D\n"
3041e284482Svaleriabarra                          "    Final rnorm                        : %e\n",
305*9b072555Sjeremylt                          ksp_type, KSPConvergedReasons[reason], its,
3061e284482Svaleriabarra                          (double)rnorm); CHKERRQ(ierr);
3071e284482Svaleriabarra     }
3081e284482Svaleriabarra     if (!rp->test_mode) {
3091e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,"  Performance:\n"); CHKERRQ(ierr);
3101e284482Svaleriabarra     }
3111e284482Svaleriabarra     {
312*9b072555Sjeremylt       PetscReal max_error;
313*9b072555Sjeremylt       ierr = ComputeErrorMax(user_O, op_error, X, target, &max_error);
3141e284482Svaleriabarra       CHKERRQ(ierr);
3151e284482Svaleriabarra       PetscReal tol = 5e-2;
316*9b072555Sjeremylt       if (!rp->test_mode || max_error > tol) {
3171e284482Svaleriabarra         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
3181e284482Svaleriabarra         CHKERRQ(ierr);
3191e284482Svaleriabarra         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm);
3201e284482Svaleriabarra         CHKERRQ(ierr);
3211e284482Svaleriabarra         ierr = PetscPrintf(rp->comm,
3221e284482Svaleriabarra                            "    Pointwise Error (max)              : %e\n"
3231e284482Svaleriabarra                            "    CG Solve Time                      : %g (%g) sec\n",
324*9b072555Sjeremylt                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
3251e284482Svaleriabarra       }
3261e284482Svaleriabarra     }
3274c583f1fSvaleriabarra     if (!rp->test_mode) {
3281e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,
3291e284482Svaleriabarra                          "    DoFs/Sec in CG                     : %g (%g) million\n",
330*9b072555Sjeremylt                          1e-6*g_size*its/rt_max,
331*9b072555Sjeremylt                          1e-6*g_size*its/rt_min); CHKERRQ(ierr);
3321e284482Svaleriabarra     }
3331e284482Svaleriabarra   }
3341e284482Svaleriabarra 
3351e284482Svaleriabarra   if (rp->write_solution) {
336*9b072555Sjeremylt     PetscViewer vtk_viewer_soln;
3371e284482Svaleriabarra 
338*9b072555Sjeremylt     ierr = PetscViewerCreate(rp->comm, &vtk_viewer_soln); CHKERRQ(ierr);
339*9b072555Sjeremylt     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
340*9b072555Sjeremylt     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
341*9b072555Sjeremylt     ierr = VecView(X, vtk_viewer_soln); CHKERRQ(ierr);
342*9b072555Sjeremylt     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
3431e284482Svaleriabarra   }
3441e284482Svaleriabarra 
3451e284482Svaleriabarra   // Cleanup
3461e284482Svaleriabarra   ierr = VecDestroy(&X); CHKERRQ(ierr);
347*9b072555Sjeremylt   ierr = VecDestroy(&X_loc); CHKERRQ(ierr);
348*9b072555Sjeremylt   ierr = VecDestroy(&user_O->Y_loc); CHKERRQ(ierr);
349*9b072555Sjeremylt   ierr = MatDestroy(&mat_O); CHKERRQ(ierr);
350*9b072555Sjeremylt   ierr = PetscFree(user_O); CHKERRQ(ierr);
351*9b072555Sjeremylt   ierr = CeedDataDestroy(0, ceed_data); CHKERRQ(ierr);
3521e284482Svaleriabarra 
3531e284482Svaleriabarra   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
354*9b072555Sjeremylt   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
3551e284482Svaleriabarra   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
3561e284482Svaleriabarra   CeedVectorDestroy(&target);
357*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
358*9b072555Sjeremylt   CeedOperatorDestroy(&op_error);
3591e284482Svaleriabarra   CeedDestroy(&ceed);
3601e284482Svaleriabarra   PetscFunctionReturn(0);
3611e284482Svaleriabarra }
3621e284482Svaleriabarra 
363e83e87a5Sjeremylt static PetscErrorCode Run(RunParams rp, PetscInt num_resources,
364*9b072555Sjeremylt                           char *const *ceed_resources, PetscInt num_bp_choices,
365*9b072555Sjeremylt                           const BPType *bp_choices) {
3665f284d84SJed Brown   PetscInt ierr;
3675f284d84SJed Brown   DM dm;
3685f284d84SJed Brown 
3695f284d84SJed Brown   PetscFunctionBeginUser;
3705f284d84SJed Brown   // Setup DM
3715f284d84SJed Brown   if (rp->read_mesh) {
3725f284d84SJed Brown     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, PETSC_TRUE, &dm);
3735f284d84SJed Brown     CHKERRQ(ierr);
3745f284d84SJed Brown   } else {
375*9b072555Sjeremylt     if (rp->user_l_nodes) {
3765f284d84SJed Brown       // Find a nicely composite number of elements no less than global nodes
3775f284d84SJed Brown       PetscMPIInt size;
3785f284d84SJed Brown       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
379*9b072555Sjeremylt       for (PetscInt g_elem =
380*9b072555Sjeremylt              PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim));
3815f284d84SJed Brown            ;
382*9b072555Sjeremylt            g_elem++) {
383*9b072555Sjeremylt         Split3(g_elem, rp->mesh_elem, true);
384*9b072555Sjeremylt         if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break;
3855f284d84SJed Brown       }
3865f284d84SJed Brown     }
387*9b072555Sjeremylt     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, PETSC_FALSE,
388*9b072555Sjeremylt                                rp->mesh_elem,
3895f284d84SJed Brown                                NULL, NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
3905f284d84SJed Brown   }
3915f284d84SJed Brown 
3925f284d84SJed Brown   {
393*9b072555Sjeremylt     DM dm_dist = NULL;
3945f284d84SJed Brown     PetscPartitioner part;
3955f284d84SJed Brown 
3965f284d84SJed Brown     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
3975f284d84SJed Brown     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
398*9b072555Sjeremylt     ierr = DMPlexDistribute(dm, 0, NULL, &dm_dist); CHKERRQ(ierr);
399*9b072555Sjeremylt     if (dm_dist) {
4005f284d84SJed Brown       ierr = DMDestroy(&dm); CHKERRQ(ierr);
401*9b072555Sjeremylt       dm  = dm_dist;
4025f284d84SJed Brown     }
4035f284d84SJed Brown   }
404b68a8d79SJed Brown   // Disable default VECSTANDARD *after* distribution (which creates a Vec)
405b68a8d79SJed Brown   ierr = DMSetVecType(dm, NULL); CHKERRQ(ierr);
4065f284d84SJed Brown 
407*9b072555Sjeremylt   for (PetscInt b = 0; b < num_bp_choices; b++) {
4085f284d84SJed Brown     DM dm_deg;
409*9b072555Sjeremylt     VecType vec_type;
410*9b072555Sjeremylt     PetscInt q_extra = rp->q_extra;
411*9b072555Sjeremylt     rp->bp_choice = bp_choices[b];
412*9b072555Sjeremylt     rp->num_comp_u = bp_options[rp->bp_choice].num_comp_u;
413*9b072555Sjeremylt     rp->q_extra = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra;
4145f284d84SJed Brown     ierr = DMClone(dm, &dm_deg); CHKERRQ(ierr);
415*9b072555Sjeremylt     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
416*9b072555Sjeremylt     ierr = DMSetVecType(dm_deg, vec_type); CHKERRQ(ierr);
4175f284d84SJed Brown     // Create DM
418e83e87a5Sjeremylt     PetscInt dim;
419e83e87a5Sjeremylt     ierr = DMGetDimension(dm_deg, &dim); CHKERRQ(ierr);
420*9b072555Sjeremylt     ierr = SetupDMByDegree(dm_deg, rp->degree, rp->num_comp_u, dim,
421*9b072555Sjeremylt                            bp_options[rp->bp_choice].enforce_bc,
422*9b072555Sjeremylt                            bp_options[rp->bp_choice].bc_func); CHKERRQ(ierr);
4235f284d84SJed Brown     for (PetscInt r = 0; r < num_resources; r++) {
424*9b072555Sjeremylt       ierr = RunWithDM(rp, dm_deg, ceed_resources[r]); CHKERRQ(ierr);
4255f284d84SJed Brown     }
4265f284d84SJed Brown     ierr = DMDestroy(&dm_deg); CHKERRQ(ierr);
427*9b072555Sjeremylt     rp->q_extra = q_extra;
4285f284d84SJed Brown   }
4295f284d84SJed Brown 
4305f284d84SJed Brown   ierr = DMDestroy(&dm); CHKERRQ(ierr);
4315f284d84SJed Brown   PetscFunctionReturn(0);
4325f284d84SJed Brown }
4335f284d84SJed Brown 
4341e284482Svaleriabarra int main(int argc, char **argv) {
435*9b072555Sjeremylt   PetscInt ierr, comm_size;
4361e284482Svaleriabarra   RunParams rp;
43753a0f73bSJed Brown   MPI_Comm comm;
43853a0f73bSJed Brown   char filename[PETSC_MAX_PATH_LEN];
439*9b072555Sjeremylt   char *ceed_resources[30];
440*9b072555Sjeremylt   PetscInt num_ceed_resources = 30;
44153a0f73bSJed Brown   char hostname[PETSC_MAX_PATH_LEN];
4421e284482Svaleriabarra 
443*9b072555Sjeremylt   PetscInt dim = 3, mesh_elem[3] = {3, 3, 3};
444*9b072555Sjeremylt   PetscInt num_degrees = 30, degree[30] = {}, num_local_nodes = 2,
445*9b072555Sjeremylt                                           local_nodes[2] = {};
446*9b072555Sjeremylt   PetscMPIInt ranks_per_node;
447c36f77d8SJed Brown   PetscBool degree_set;
448*9b072555Sjeremylt   BPType bp_choices[10];
449*9b072555Sjeremylt   PetscInt num_bp_choices = 10;
4500c59ef15SJeremy L Thompson 
4511e284482Svaleriabarra   // Initialize PETSc
4520c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
4530c59ef15SJeremy L Thompson   if (ierr) return ierr;
4540c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
455*9b072555Sjeremylt   ierr = MPI_Comm_size(comm, &comm_size);
4561e284482Svaleriabarra   if (ierr != MPI_SUCCESS) return ierr;
45753a0f73bSJed Brown   #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
45853a0f73bSJed Brown   {
45953a0f73bSJed Brown     MPI_Comm splitcomm;
4601e284482Svaleriabarra     ierr = MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL,
4611e284482Svaleriabarra                                &splitcomm);
46253a0f73bSJed Brown     CHKERRQ(ierr);
463*9b072555Sjeremylt     ierr = MPI_Comm_size(splitcomm, &ranks_per_node); CHKERRQ(ierr);
46453a0f73bSJed Brown     ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr);
46553a0f73bSJed Brown   }
46653a0f73bSJed Brown   #else
467*9b072555Sjeremylt   ranks_per_node = -1; // Unknown
46853a0f73bSJed Brown   #endif
469cb32e2e7SValeria Barra 
470c36f77d8SJed Brown   // Setup all parameters needed in Run()
471c36f77d8SJed Brown   ierr = PetscMalloc1(1, &rp); CHKERRQ(ierr);
472c36f77d8SJed Brown   rp->comm = comm;
473c36f77d8SJed Brown 
47432d2ee49SValeria Barra   // Read command line options
4751e284482Svaleriabarra   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
4761e284482Svaleriabarra   CHKERRQ(ierr);
477565a3730SJed Brown   {
478565a3730SJed Brown     PetscBool set;
479565a3730SJed Brown     ierr = PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve",
480565a3730SJed Brown                                  NULL,
481*9b072555Sjeremylt                                  bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set);
482565a3730SJed Brown     CHKERRQ(ierr);
483565a3730SJed Brown     if (!set) {
484*9b072555Sjeremylt       bp_choices[0] = CEED_BP1;
485*9b072555Sjeremylt       num_bp_choices = 1;
486565a3730SJed Brown     }
487565a3730SJed Brown   }
488c36f77d8SJed Brown   rp->test_mode = PETSC_FALSE;
4890c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
4900c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
491c36f77d8SJed Brown                           NULL, rp->test_mode, &rp->test_mode, NULL); CHKERRQ(ierr);
492c36f77d8SJed Brown   rp->write_solution = PETSC_FALSE;
49353a0f73bSJed Brown   ierr = PetscOptionsBool("-write_solution", "Write solution for visualization",
494c36f77d8SJed Brown                           NULL, rp->write_solution, &rp->write_solution, NULL);
4956bd9afcaSjeremylt   CHKERRQ(ierr);
496c36f77d8SJed Brown   degree[0] = rp->test_mode ? 3 : 2;
49753a0f73bSJed Brown   ierr = PetscOptionsIntArray("-degree",
49853a0f73bSJed Brown                               "Polynomial degree of tensor product basis", NULL,
49953a0f73bSJed Brown                               degree, &num_degrees, &degree_set); CHKERRQ(ierr);
50053a0f73bSJed Brown   if (!degree_set)
50153a0f73bSJed Brown     num_degrees = 1;
502*9b072555Sjeremylt   rp->q_extra = PETSC_DECIDE;
503*9b072555Sjeremylt   ierr = PetscOptionsInt("-q_extra",
5045f284d84SJed Brown                          "Number of extra quadrature points (-1 for auto)", NULL,
505*9b072555Sjeremylt                          rp->q_extra, &rp->q_extra, NULL); CHKERRQ(ierr);
506565a3730SJed Brown   {
507565a3730SJed Brown     PetscBool set;
508565a3730SJed Brown     ierr = PetscOptionsStringArray("-ceed",
509565a3730SJed Brown                                    "CEED resource specifier (comma-separated list)", NULL,
510*9b072555Sjeremylt                                    ceed_resources, &num_ceed_resources, &set); CHKERRQ(ierr);
511565a3730SJed Brown     if (!set) {
512*9b072555Sjeremylt       ierr = PetscStrallocpy( "/cpu/self", &ceed_resources[0]); CHKERRQ(ierr);
513*9b072555Sjeremylt       num_ceed_resources = 1;
514565a3730SJed Brown     }
515565a3730SJed Brown   }
51653a0f73bSJed Brown   ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
51753a0f73bSJed Brown   ierr = PetscOptionsString("-hostname", "Hostname for output", NULL, hostname,
51853a0f73bSJed Brown                             hostname, sizeof(hostname), NULL); CHKERRQ(ierr);
519c36f77d8SJed Brown   rp->read_mesh = PETSC_FALSE;
52053a0f73bSJed Brown   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, filename,
521c36f77d8SJed Brown                             filename, sizeof(filename), &rp->read_mesh);
5220c59ef15SJeremy L Thompson   CHKERRQ(ierr);
523c36f77d8SJed Brown   rp->filename = filename;
524c36f77d8SJed Brown   if (!rp->read_mesh) {
525cb32e2e7SValeria Barra     PetscInt tmp = dim;
526cb32e2e7SValeria Barra     ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL,
527*9b072555Sjeremylt                                 mesh_elem, &tmp, NULL); CHKERRQ(ierr);
528cb32e2e7SValeria Barra   }
529*9b072555Sjeremylt   local_nodes[0] = 1000;
53053a0f73bSJed Brown   ierr = PetscOptionsIntArray("-local_nodes",
53153a0f73bSJed Brown                               "Target number of locally owned nodes per "
53253a0f73bSJed Brown                               "process (single value or min,max)",
533*9b072555Sjeremylt                               NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes);
534153bcb04Svaleriabarra   CHKERRQ(ierr);
535*9b072555Sjeremylt   if (num_local_nodes < 2)
536*9b072555Sjeremylt     local_nodes[1] = 2 * local_nodes[0];
537c36f77d8SJed Brown   {
538c36f77d8SJed Brown     PetscInt two = 2;
539c36f77d8SJed Brown     rp->ksp_max_it_clip[0] = 5;
540c36f77d8SJed Brown     rp->ksp_max_it_clip[1] = 20;
541c36f77d8SJed Brown     ierr = PetscOptionsIntArray("-ksp_max_it_clip",
542c36f77d8SJed Brown                                 "Min and max number of iterations to use during benchmarking",
543c36f77d8SJed Brown                                 NULL, rp->ksp_max_it_clip, &two, NULL); CHKERRQ(ierr);
544c36f77d8SJed Brown   }
545da9108adSvaleriabarra   if (!degree_set) {
546*9b072555Sjeremylt     PetscInt max_degree = 8;
5475f284d84SJed Brown     ierr = PetscOptionsInt("-max_degree",
548*9b072555Sjeremylt                            "Range of degrees [1, max_degree] to run with",
549*9b072555Sjeremylt                            NULL, max_degree, &max_degree, NULL);
5501e284482Svaleriabarra     CHKERRQ(ierr);
551*9b072555Sjeremylt     for (PetscInt i = 0; i < max_degree; i++)
55253a0f73bSJed Brown       degree[i] = i + 1;
553*9b072555Sjeremylt     num_degrees = max_degree;
55453a0f73bSJed Brown   }
55553a0f73bSJed Brown   {
55653a0f73bSJed Brown     PetscBool flg;
557*9b072555Sjeremylt     PetscInt p = ranks_per_node;
55853a0f73bSJed Brown     ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL,
55953a0f73bSJed Brown                            p, &p, &flg);
56053a0f73bSJed Brown     CHKERRQ(ierr);
561*9b072555Sjeremylt     if (flg) ranks_per_node = p;
56253a0f73bSJed Brown   }
5639396343dSjeremylt 
56453a0f73bSJed Brown   ierr = PetscOptionsEnd();
56553a0f73bSJed Brown   CHKERRQ(ierr);
566cb32e2e7SValeria Barra 
5671e284482Svaleriabarra   // Register PETSc logging stage
568*9b072555Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &rp->solve_stage);
56953a0f73bSJed Brown   CHKERRQ(ierr);
57009a940d7Sjeremylt 
57153a0f73bSJed Brown   rp->hostname = hostname;
5721e284482Svaleriabarra   rp->dim = dim;
573*9b072555Sjeremylt   rp->mesh_elem = mesh_elem;
574*9b072555Sjeremylt   rp->ranks_per_node = ranks_per_node;
575cb32e2e7SValeria Barra 
57653a0f73bSJed Brown   for (PetscInt d = 0; d < num_degrees; d++) {
57753a0f73bSJed Brown     PetscInt deg = degree[d];
578*9b072555Sjeremylt     for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) {
57953a0f73bSJed Brown       rp->degree = deg;
580*9b072555Sjeremylt       rp->local_nodes = n;
581*9b072555Sjeremylt       ierr = Run(rp, num_ceed_resources, ceed_resources,
582*9b072555Sjeremylt                  num_bp_choices, bp_choices); CHKERRQ(ierr);
583565a3730SJed Brown     }
584565a3730SJed Brown   }
5851e284482Svaleriabarra   // Clear memory
5861e284482Svaleriabarra   ierr = PetscFree(rp); CHKERRQ(ierr);
587*9b072555Sjeremylt   for (PetscInt i=0; i<num_ceed_resources; i++) {
588*9b072555Sjeremylt     ierr = PetscFree(ceed_resources[i]); CHKERRQ(ierr);
589565a3730SJed Brown   }
5900c59ef15SJeremy L Thompson   return PetscFinalize();
5910c59ef15SJeremy L Thompson }
592