xref: /libCEED/examples/petsc/bps.c (revision 6a6c615b31508fbb9571bc7a279f860841ca2097)
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
3132d2ee49SValeria Barra //     ./bps -problem bp2 -ceed /cpu/self -degree 3
3232d2ee49SValeria Barra //     ./bps -problem bp3 -ceed /gpu/occa -degree 3
3332d2ee49SValeria Barra //     ./bps -problem bp4 -ceed /cpu/occa -degree 3
3432d2ee49SValeria Barra //     ./bps -problem bp5 -ceed /omp/occa -degree 3
3532d2ee49SValeria Barra //     ./bps -problem bp6 -ceed /ocl/occa -degree 3
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 
440c59ef15SJeremy L Thompson #include <stdbool.h>
450c59ef15SJeremy L Thompson #include <string.h>
464d537eeaSYohann #include <petscksp.h>
47cb32e2e7SValeria Barra #include <petscdmplex.h>
484d537eeaSYohann #include <ceed.h>
49cb32e2e7SValeria Barra #include "setup.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) {
57d5b2ba77SJed Brown   for (PetscInt d=0,sizeleft=size; d<3; d++) {
58d5b2ba77SJed Brown     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d)));
59d5b2ba77SJed Brown     while (try * (sizeleft / try) != sizeleft) try++;
60d5b2ba77SJed Brown     m[reverse ? 2-d : d] = try;
61d5b2ba77SJed Brown     sizeleft /= 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;
79da9108adSvaleriabarra   PetscBool test_mode, read_mesh, userlnodes, setmemtyperequest,
801e284482Svaleriabarra             petschavecuda, write_solution;
815f284d84SJed Brown   char *filename, *hostname;
821e284482Svaleriabarra   PetscInt localnodes, degree, qextra, dim, ncompu, *melem;
83c36f77d8SJed Brown   PetscInt ksp_max_it_clip[2];
8453a0f73bSJed Brown   PetscMPIInt rankspernode;
851e284482Svaleriabarra   bpType bpchoice;
861e284482Svaleriabarra   CeedMemType memtyperequested;
8709a940d7Sjeremylt   PetscLogStage solvestage;
881e284482Svaleriabarra };
891e284482Svaleriabarra 
901e284482Svaleriabarra // -----------------------------------------------------------------------------
911e284482Svaleriabarra // Main body of program, called in a loop for performance benchmarking purposes
921e284482Svaleriabarra // -----------------------------------------------------------------------------
935f284d84SJed Brown static PetscErrorCode RunWithDM(RunParams rp, DM dm,
945f284d84SJed Brown                                 const char *ceedresource) {
955f284d84SJed Brown   PetscErrorCode ierr;
961e284482Svaleriabarra   double my_rt_start, my_rt, rt_min, rt_max;
971e284482Svaleriabarra   PetscInt xlsize, lsize, gsize;
981e284482Svaleriabarra   PetscScalar *r;
99819eb1b3Sjeremylt   Vec X, Xloc, rhs, rhsloc;
100cb32e2e7SValeria Barra   Mat matO;
101819eb1b3Sjeremylt   KSP ksp;
102cb32e2e7SValeria Barra   UserO userO;
1030c59ef15SJeremy L Thompson   Ceed ceed;
104cb32e2e7SValeria Barra   CeedData ceeddata;
105c70bd2a0Sjeremylt   CeedQFunction qferror;
106c70bd2a0Sjeremylt   CeedOperator operror;
107cb32e2e7SValeria Barra   CeedVector rhsceed, target;
1081e284482Svaleriabarra 
1091e284482Svaleriabarra   PetscFunctionBeginUser;
1101e284482Svaleriabarra   // Set up libCEED
1115f284d84SJed Brown   CeedInit(ceedresource, &ceed);
1121e284482Svaleriabarra   CeedMemType memtypebackend;
1131e284482Svaleriabarra   CeedGetPreferredMemType(ceed, &memtypebackend);
1141e284482Svaleriabarra 
1151e284482Svaleriabarra   // Check memtype compatibility
1161e284482Svaleriabarra   if (!rp->setmemtyperequest)
1171e284482Svaleriabarra     rp->memtyperequested = memtypebackend;
1181e284482Svaleriabarra   else if (!rp->petschavecuda && rp->memtyperequested == CEED_MEM_DEVICE)
1191e284482Svaleriabarra     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
1201e284482Svaleriabarra              "PETSc was not built with CUDA. "
1211e284482Svaleriabarra              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
1221e284482Svaleriabarra 
1231e284482Svaleriabarra   // Create vectors
1241e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_DEVICE) {
1251e284482Svaleriabarra     ierr = DMSetVecType(dm, VECCUDA); CHKERRQ(ierr);
1261e284482Svaleriabarra   }
1271e284482Svaleriabarra   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
1281e284482Svaleriabarra   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
1291e284482Svaleriabarra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
1301e284482Svaleriabarra   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
1311e284482Svaleriabarra   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
1321e284482Svaleriabarra   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
1331e284482Svaleriabarra 
1341e284482Svaleriabarra   // Operator
1351e284482Svaleriabarra   ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr);
1361e284482Svaleriabarra   ierr = MatCreateShell(rp->comm, lsize, lsize, gsize, gsize,
1371e284482Svaleriabarra                         userO, &matO); CHKERRQ(ierr);
1381e284482Svaleriabarra   ierr = MatShellSetOperation(matO, MATOP_MULT,
1391e284482Svaleriabarra                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
1401e284482Svaleriabarra   ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL,
1411e284482Svaleriabarra                               (void(*)(void))MatGetDiag); CHKERRQ(ierr);
1421e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_DEVICE) {
1431e284482Svaleriabarra     ierr = MatShellSetVecType(matO, VECCUDA); CHKERRQ(ierr);
1441e284482Svaleriabarra   }
1451e284482Svaleriabarra 
1461e284482Svaleriabarra   // Print summary
1471e284482Svaleriabarra   if (!rp->test_mode) {
1481e284482Svaleriabarra     PetscInt P = rp->degree + 1, Q = P + rp->qextra;
1491e284482Svaleriabarra 
1501e284482Svaleriabarra     const char *usedresource;
1511e284482Svaleriabarra     CeedGetResource(ceed, &usedresource);
1521e284482Svaleriabarra 
1531e284482Svaleriabarra     VecType vectype;
1541e284482Svaleriabarra     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
1551e284482Svaleriabarra 
1561e284482Svaleriabarra     PetscInt cStart, cEnd;
1571e284482Svaleriabarra     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
15853a0f73bSJed Brown     PetscMPIInt comm_size;
15953a0f73bSJed Brown     ierr = MPI_Comm_size(rp->comm, &comm_size); CHKERRQ(ierr);
1601e284482Svaleriabarra     ierr = PetscPrintf(rp->comm,
1611e284482Svaleriabarra                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
16253a0f73bSJed Brown                        "  MPI:\n"
16353a0f73bSJed Brown                        "    Hostname                           : %s\n"
16453a0f73bSJed Brown                        "    Total ranks                        : %d\n"
16553a0f73bSJed Brown                        "    Ranks per compute node             : %d\n"
1661e284482Svaleriabarra                        "  PETSc:\n"
1671e284482Svaleriabarra                        "    PETSc Vec Type                     : %s\n"
1681e284482Svaleriabarra                        "  libCEED:\n"
1691e284482Svaleriabarra                        "    libCEED Backend                    : %s\n"
1701e284482Svaleriabarra                        "    libCEED Backend MemType            : %s\n"
1711e284482Svaleriabarra                        "    libCEED User Requested MemType     : %s\n"
1721e284482Svaleriabarra                        "  Mesh:\n"
1731e284482Svaleriabarra                        "    Number of 1D Basis Nodes (P)       : %d\n"
1741e284482Svaleriabarra                        "    Number of 1D Quadrature Points (Q) : %d\n"
1751e284482Svaleriabarra                        "    Global nodes                       : %D\n"
1761e284482Svaleriabarra                        "    Local Elements                     : %D\n"
1771e284482Svaleriabarra                        "    Owned nodes                        : %D\n"
1781e284482Svaleriabarra                        "    DoF per node                       : %D\n",
17953a0f73bSJed Brown                        rp->bpchoice+1,
18053a0f73bSJed Brown                        rp->hostname,
18153a0f73bSJed Brown                        comm_size, rp->rankspernode,
18253a0f73bSJed Brown                        vectype, usedresource,
1831e284482Svaleriabarra                        CeedMemTypes[memtypebackend],
1841e284482Svaleriabarra                        (rp->setmemtyperequest) ?
1851e284482Svaleriabarra                        CeedMemTypes[rp->memtyperequested] : "none",
1861e284482Svaleriabarra                        P, Q, gsize/rp->ncompu, cEnd - cStart, lsize/rp->ncompu,
1871e284482Svaleriabarra                        rp->ncompu);
1881e284482Svaleriabarra     CHKERRQ(ierr);
1891e284482Svaleriabarra   }
1901e284482Svaleriabarra 
1911e284482Svaleriabarra   // Create RHS vector
1921e284482Svaleriabarra   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
1931e284482Svaleriabarra   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
1941e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_HOST) {
1951e284482Svaleriabarra     ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
1961e284482Svaleriabarra   } else {
1971e284482Svaleriabarra     ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr);
1981e284482Svaleriabarra   }
1991e284482Svaleriabarra   CeedVectorCreate(ceed, xlsize, &rhsceed);
2001e284482Svaleriabarra   CeedVectorSetArray(rhsceed, rp->memtyperequested, CEED_USE_POINTER, r);
2011e284482Svaleriabarra 
2021e284482Svaleriabarra   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
2031e284482Svaleriabarra   ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->qextra,
2041e284482Svaleriabarra                               rp->ncompu, gsize, xlsize, rp->bpchoice, ceeddata,
2051e284482Svaleriabarra                               true, rhsceed, &target); CHKERRQ(ierr);
2061e284482Svaleriabarra 
2071e284482Svaleriabarra   // Gather RHS
208*6a6c615bSJeremy L Thompson   CeedVectorTakeArray(rhsceed, rp->memtyperequested, NULL);
2091e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_HOST) {
2101e284482Svaleriabarra     ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
2111e284482Svaleriabarra   } else {
2121e284482Svaleriabarra     ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr);
2131e284482Svaleriabarra   }
2141e284482Svaleriabarra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
2151e284482Svaleriabarra   ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
2161e284482Svaleriabarra   ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
2171e284482Svaleriabarra   CeedVectorDestroy(&rhsceed);
2181e284482Svaleriabarra 
219*6a6c615bSJeremy L Thompson   // Create the error QFunction
2201e284482Svaleriabarra   CeedQFunctionCreateInterior(ceed, 1, bpOptions[rp->bpchoice].error,
2211e284482Svaleriabarra                               bpOptions[rp->bpchoice].errorfname, &qferror);
2221e284482Svaleriabarra   CeedQFunctionAddInput(qferror, "u", rp->ncompu, CEED_EVAL_INTERP);
2231e284482Svaleriabarra   CeedQFunctionAddInput(qferror, "true_soln", rp->ncompu, CEED_EVAL_NONE);
2241e284482Svaleriabarra   CeedQFunctionAddOutput(qferror, "error", rp->ncompu, CEED_EVAL_NONE);
2251e284482Svaleriabarra 
2261e284482Svaleriabarra   // Create the error operator
2271e284482Svaleriabarra   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
2281e284482Svaleriabarra                      &operror);
2291e284482Svaleriabarra   CeedOperatorSetField(operror, "u", ceeddata->Erestrictu,
2301e284482Svaleriabarra                        ceeddata->basisu, CEED_VECTOR_ACTIVE);
2311e284482Svaleriabarra   CeedOperatorSetField(operror, "true_soln", ceeddata->Erestrictui,
2321e284482Svaleriabarra                        CEED_BASIS_COLLOCATED, target);
2331e284482Svaleriabarra   CeedOperatorSetField(operror, "error", ceeddata->Erestrictui,
2341e284482Svaleriabarra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
2351e284482Svaleriabarra 
2361e284482Svaleriabarra   // Set up Mat
2371e284482Svaleriabarra   userO->comm = rp->comm;
2381e284482Svaleriabarra   userO->dm = dm;
2391e284482Svaleriabarra   userO->Xloc = Xloc;
2401e284482Svaleriabarra   ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr);
2411e284482Svaleriabarra   userO->xceed = ceeddata->xceed;
2421e284482Svaleriabarra   userO->yceed = ceeddata->yceed;
2431e284482Svaleriabarra   userO->op = ceeddata->opapply;
2441e284482Svaleriabarra   userO->ceed = ceed;
2451e284482Svaleriabarra   userO->memtype = rp->memtyperequested;
2461e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_HOST) {
2471e284482Svaleriabarra     userO->VecGetArray = VecGetArray;
2481e284482Svaleriabarra     userO->VecGetArrayRead = VecGetArrayRead;
2491e284482Svaleriabarra     userO->VecRestoreArray = VecRestoreArray;
2501e284482Svaleriabarra     userO->VecRestoreArrayRead = VecRestoreArrayRead;
2511e284482Svaleriabarra   } else {
2521e284482Svaleriabarra     userO->VecGetArray = VecCUDAGetArray;
2531e284482Svaleriabarra     userO->VecGetArrayRead = VecCUDAGetArrayRead;
2541e284482Svaleriabarra     userO->VecRestoreArray = VecCUDARestoreArray;
2551e284482Svaleriabarra     userO->VecRestoreArrayRead = VecCUDARestoreArrayRead;
2561e284482Svaleriabarra   }
2571e284482Svaleriabarra 
2581e284482Svaleriabarra   ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr);
2591e284482Svaleriabarra   {
2601e284482Svaleriabarra     PC pc;
2611e284482Svaleriabarra     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
2621e284482Svaleriabarra     if (rp->bpchoice == CEED_BP1 || rp->bpchoice == CEED_BP2) {
2631e284482Svaleriabarra       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
2641e284482Svaleriabarra       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
2651e284482Svaleriabarra     } else {
2661e284482Svaleriabarra       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
2671e284482Svaleriabarra     }
2681e284482Svaleriabarra     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
2691e284482Svaleriabarra     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
2701e284482Svaleriabarra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
2711e284482Svaleriabarra                             PETSC_DEFAULT); CHKERRQ(ierr);
2721e284482Svaleriabarra   }
2731e284482Svaleriabarra   ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr);
2741e284482Svaleriabarra 
275da9108adSvaleriabarra   // First run's performance log is not considered for benchmarking purposes
2761e284482Svaleriabarra   ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
2771e284482Svaleriabarra   CHKERRQ(ierr);
2781e284482Svaleriabarra   my_rt_start = MPI_Wtime();
2791e284482Svaleriabarra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
2801e284482Svaleriabarra   my_rt = MPI_Wtime() - my_rt_start;
2811e284482Svaleriabarra   ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
2821e284482Svaleriabarra   CHKERRQ(ierr);
2831e284482Svaleriabarra   // Set maxits based on first iteration timing
2841e284482Svaleriabarra   if (my_rt > 0.02) {
285565a3730SJed Brown     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
286565a3730SJed Brown                             rp->ksp_max_it_clip[0]);
2871e284482Svaleriabarra     CHKERRQ(ierr);
2881e284482Svaleriabarra   } else {
289565a3730SJed Brown     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
290565a3730SJed Brown                             rp->ksp_max_it_clip[1]);
2911e284482Svaleriabarra     CHKERRQ(ierr);
2921e284482Svaleriabarra   }
2931e284482Svaleriabarra   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
2941e284482Svaleriabarra 
2951e284482Svaleriabarra   // Timed solve
2961e284482Svaleriabarra   ierr = VecZeroEntries(X); CHKERRQ(ierr);
2971e284482Svaleriabarra   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
2981e284482Svaleriabarra 
2991e284482Svaleriabarra   // -- Performance logging
3001e284482Svaleriabarra   ierr = PetscLogStagePush(rp->solvestage); CHKERRQ(ierr);
3011e284482Svaleriabarra 
3021e284482Svaleriabarra   // -- Solve
3031e284482Svaleriabarra   my_rt_start = MPI_Wtime();
3041e284482Svaleriabarra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
3051e284482Svaleriabarra   my_rt = MPI_Wtime() - my_rt_start;
3061e284482Svaleriabarra 
3071e284482Svaleriabarra   // -- Performance logging
3081e284482Svaleriabarra   ierr = PetscLogStagePop();
3091e284482Svaleriabarra 
3101e284482Svaleriabarra   // Output results
3111e284482Svaleriabarra   {
3121e284482Svaleriabarra     KSPType ksptype;
3131e284482Svaleriabarra     KSPConvergedReason reason;
3141e284482Svaleriabarra     PetscReal rnorm;
3151e284482Svaleriabarra     PetscInt its;
3161e284482Svaleriabarra     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
3171e284482Svaleriabarra     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
3181e284482Svaleriabarra     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
3191e284482Svaleriabarra     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
3201e284482Svaleriabarra     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
3211e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,
3221e284482Svaleriabarra                          "  KSP:\n"
3231e284482Svaleriabarra                          "    KSP Type                           : %s\n"
3241e284482Svaleriabarra                          "    KSP Convergence                    : %s\n"
3251e284482Svaleriabarra                          "    Total KSP Iterations               : %D\n"
3261e284482Svaleriabarra                          "    Final rnorm                        : %e\n",
3271e284482Svaleriabarra                          ksptype, KSPConvergedReasons[reason], its,
3281e284482Svaleriabarra                          (double)rnorm); CHKERRQ(ierr);
3291e284482Svaleriabarra     }
3301e284482Svaleriabarra     if (!rp->test_mode) {
3311e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,"  Performance:\n"); CHKERRQ(ierr);
3321e284482Svaleriabarra     }
3331e284482Svaleriabarra     {
3341e284482Svaleriabarra       PetscReal maxerror;
3351e284482Svaleriabarra       ierr = ComputeErrorMax(userO, operror, X, target, &maxerror);
3361e284482Svaleriabarra       CHKERRQ(ierr);
3371e284482Svaleriabarra       PetscReal tol = 5e-2;
3381e284482Svaleriabarra       if (!rp->test_mode || maxerror > tol) {
3391e284482Svaleriabarra         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
3401e284482Svaleriabarra         CHKERRQ(ierr);
3411e284482Svaleriabarra         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm);
3421e284482Svaleriabarra         CHKERRQ(ierr);
3431e284482Svaleriabarra         ierr = PetscPrintf(rp->comm,
3441e284482Svaleriabarra                            "    Pointwise Error (max)              : %e\n"
3451e284482Svaleriabarra                            "    CG Solve Time                      : %g (%g) sec\n",
3461e284482Svaleriabarra                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
3471e284482Svaleriabarra       }
3481e284482Svaleriabarra     }
3494c583f1fSvaleriabarra     if (!rp->test_mode) {
3501e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,
3511e284482Svaleriabarra                          "    DoFs/Sec in CG                     : %g (%g) million\n",
3521e284482Svaleriabarra                          1e-6*gsize*its/rt_max,
3531e284482Svaleriabarra                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
3541e284482Svaleriabarra     }
3551e284482Svaleriabarra   }
3561e284482Svaleriabarra 
3571e284482Svaleriabarra   if (rp->write_solution) {
3581e284482Svaleriabarra     PetscViewer vtkviewersoln;
3591e284482Svaleriabarra 
3601e284482Svaleriabarra     ierr = PetscViewerCreate(rp->comm, &vtkviewersoln); CHKERRQ(ierr);
3611e284482Svaleriabarra     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
3621e284482Svaleriabarra     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
3631e284482Svaleriabarra     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
3641e284482Svaleriabarra     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
3651e284482Svaleriabarra   }
3661e284482Svaleriabarra 
3671e284482Svaleriabarra   // Cleanup
3681e284482Svaleriabarra   ierr = VecDestroy(&X); CHKERRQ(ierr);
3691e284482Svaleriabarra   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
3701e284482Svaleriabarra   ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr);
3711e284482Svaleriabarra   ierr = MatDestroy(&matO); CHKERRQ(ierr);
3721e284482Svaleriabarra   ierr = PetscFree(userO); CHKERRQ(ierr);
3731e284482Svaleriabarra   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
3741e284482Svaleriabarra 
3751e284482Svaleriabarra   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
3761e284482Svaleriabarra   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
3771e284482Svaleriabarra   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
3781e284482Svaleriabarra   CeedVectorDestroy(&target);
3791e284482Svaleriabarra   CeedQFunctionDestroy(&qferror);
3801e284482Svaleriabarra   CeedOperatorDestroy(&operror);
3811e284482Svaleriabarra   CeedDestroy(&ceed);
3821e284482Svaleriabarra   PetscFunctionReturn(0);
3831e284482Svaleriabarra }
3841e284482Svaleriabarra 
3855f284d84SJed Brown static PetscErrorCode Run(RunParams rp,
3865f284d84SJed Brown                           PetscInt num_resources, char *const *ceedresources,
3875f284d84SJed Brown                           PetscInt num_bpchoices, const bpType *bpchoices) {
3885f284d84SJed Brown   PetscInt ierr;
3895f284d84SJed Brown   DM  dm;
3905f284d84SJed Brown 
3915f284d84SJed Brown   PetscFunctionBeginUser;
3925f284d84SJed Brown   // Setup DM
3935f284d84SJed Brown   if (rp->read_mesh) {
3945f284d84SJed Brown     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, PETSC_TRUE, &dm);
3955f284d84SJed Brown     CHKERRQ(ierr);
3965f284d84SJed Brown   } else {
3975f284d84SJed Brown     if (rp->userlnodes) {
3985f284d84SJed Brown       // Find a nicely composite number of elements no less than global nodes
3995f284d84SJed Brown       PetscMPIInt size;
4005f284d84SJed Brown       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
4015f284d84SJed Brown       for (PetscInt gelem =
4025f284d84SJed Brown              PetscMax(1, size * rp->localnodes / PetscPowInt(rp->degree, rp->dim));
4035f284d84SJed Brown            ;
4045f284d84SJed Brown            gelem++) {
4055f284d84SJed Brown         Split3(gelem, rp->melem, true);
4065f284d84SJed Brown         if (Max3(rp->melem) / Min3(rp->melem) <= 2) break;
4075f284d84SJed Brown       }
4085f284d84SJed Brown     }
4095f284d84SJed Brown     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, PETSC_FALSE, rp->melem,
4105f284d84SJed Brown                                NULL, NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
4115f284d84SJed Brown   }
4125f284d84SJed Brown 
4135f284d84SJed Brown   {
4145f284d84SJed Brown     DM dmDist = NULL;
4155f284d84SJed Brown     PetscPartitioner part;
4165f284d84SJed Brown 
4175f284d84SJed Brown     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
4185f284d84SJed Brown     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
4195f284d84SJed Brown     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
4205f284d84SJed Brown     if (dmDist) {
4215f284d84SJed Brown       ierr = DMDestroy(&dm); CHKERRQ(ierr);
4225f284d84SJed Brown       dm  = dmDist;
4235f284d84SJed Brown     }
4245f284d84SJed Brown   }
4255f284d84SJed Brown 
4265f284d84SJed Brown   for (PetscInt b = 0; b < num_bpchoices; b++) {
4275f284d84SJed Brown     DM dm_deg;
4285f284d84SJed Brown     PetscInt qextra = rp->qextra;
4295f284d84SJed Brown     rp->bpchoice = bpchoices[b];
4305f284d84SJed Brown     rp->ncompu = bpOptions[rp->bpchoice].ncompu;
4315f284d84SJed Brown     rp->qextra = qextra < 0 ? bpOptions[rp->bpchoice].qextra : qextra;
4325f284d84SJed Brown     ierr = DMClone(dm, &dm_deg); CHKERRQ(ierr);
4335f284d84SJed Brown     // Create DM
4345f284d84SJed Brown     ierr = SetupDMByDegree(dm_deg, rp->degree, rp->ncompu, rp->bpchoice);
4355f284d84SJed Brown     CHKERRQ(ierr);
4365f284d84SJed Brown     for (PetscInt r = 0; r < num_resources; r++) {
4375f284d84SJed Brown       ierr = RunWithDM(rp, dm_deg, ceedresources[r]); CHKERRQ(ierr);
4385f284d84SJed Brown     }
4395f284d84SJed Brown     ierr = DMDestroy(&dm_deg); CHKERRQ(ierr);
4405f284d84SJed Brown     rp->qextra = qextra;
4415f284d84SJed Brown   }
4425f284d84SJed Brown 
4435f284d84SJed Brown   ierr = DMDestroy(&dm); CHKERRQ(ierr);
4445f284d84SJed Brown   PetscFunctionReturn(0);
4455f284d84SJed Brown }
4465f284d84SJed Brown 
4471e284482Svaleriabarra int main(int argc, char **argv) {
4481e284482Svaleriabarra   PetscInt ierr, commsize;
4491e284482Svaleriabarra   RunParams rp;
45053a0f73bSJed Brown   MPI_Comm comm;
45153a0f73bSJed Brown   char filename[PETSC_MAX_PATH_LEN];
452565a3730SJed Brown   char *ceedresources[30];
453565a3730SJed Brown   PetscInt num_ceedresources = 30;
45453a0f73bSJed Brown   char hostname[PETSC_MAX_PATH_LEN];
4551e284482Svaleriabarra 
456c36f77d8SJed Brown   PetscInt dim = 3, melem[3] = {3, 3, 3};
45753a0f73bSJed Brown   PetscInt num_degrees = 30, degree[30] = {}, num_localnodes = 2, localnodes[2] = {};
45853a0f73bSJed Brown   PetscMPIInt rankspernode;
459c36f77d8SJed Brown   PetscBool degree_set;
460565a3730SJed Brown   bpType bpchoices[10];
461565a3730SJed Brown   PetscInt num_bpchoices = 10;
4620c59ef15SJeremy L Thompson 
4639396343dSjeremylt   // Check PETSc CUDA support
464c36f77d8SJed Brown   PetscBool petschavecuda;
4659396343dSjeremylt   // *INDENT-OFF*
4669396343dSjeremylt   #ifdef PETSC_HAVE_CUDA
4679396343dSjeremylt   petschavecuda = PETSC_TRUE;
4689396343dSjeremylt   #else
4699396343dSjeremylt   petschavecuda = PETSC_FALSE;
4709396343dSjeremylt   #endif
4719396343dSjeremylt   // *INDENT-ON*
4729396343dSjeremylt 
4731e284482Svaleriabarra   // Initialize PETSc
4740c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
4750c59ef15SJeremy L Thompson   if (ierr) return ierr;
4760c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
4771e284482Svaleriabarra   ierr = MPI_Comm_size(comm, &commsize);
4781e284482Svaleriabarra   if (ierr != MPI_SUCCESS) return ierr;
47953a0f73bSJed Brown   #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
48053a0f73bSJed Brown   {
48153a0f73bSJed Brown     MPI_Comm splitcomm;
4821e284482Svaleriabarra     ierr = MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL,
4831e284482Svaleriabarra                                &splitcomm);
48453a0f73bSJed Brown     CHKERRQ(ierr);
48553a0f73bSJed Brown     ierr = MPI_Comm_size(splitcomm, &rankspernode); CHKERRQ(ierr);
48653a0f73bSJed Brown     ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr);
48753a0f73bSJed Brown   }
48853a0f73bSJed Brown   #else
48953a0f73bSJed Brown   rankspernode = -1; // Unknown
49053a0f73bSJed Brown   #endif
491cb32e2e7SValeria Barra 
492c36f77d8SJed Brown   // Setup all parameters needed in Run()
493c36f77d8SJed Brown   ierr = PetscMalloc1(1, &rp); CHKERRQ(ierr);
494c36f77d8SJed Brown   rp->comm = comm;
495c36f77d8SJed Brown 
49632d2ee49SValeria Barra   // Read command line options
4971e284482Svaleriabarra   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
4981e284482Svaleriabarra   CHKERRQ(ierr);
499565a3730SJed Brown   {
500565a3730SJed Brown     PetscBool set;
501565a3730SJed Brown     ierr = PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve",
502565a3730SJed Brown                                  NULL,
503565a3730SJed Brown                                  bpTypes, (PetscEnum *)bpchoices, &num_bpchoices, &set);
504565a3730SJed Brown     CHKERRQ(ierr);
505565a3730SJed Brown     if (!set) {
506565a3730SJed Brown       bpchoices[0] = CEED_BP1;
507565a3730SJed Brown       num_bpchoices = 1;
508565a3730SJed Brown     }
509565a3730SJed Brown   }
510c36f77d8SJed Brown   rp->test_mode = PETSC_FALSE;
5110c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
5120c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
513c36f77d8SJed Brown                           NULL, rp->test_mode, &rp->test_mode, NULL); CHKERRQ(ierr);
514c36f77d8SJed Brown   rp->write_solution = PETSC_FALSE;
51553a0f73bSJed Brown   ierr = PetscOptionsBool("-write_solution", "Write solution for visualization",
516c36f77d8SJed Brown                           NULL, rp->write_solution, &rp->write_solution, NULL);
5176bd9afcaSjeremylt   CHKERRQ(ierr);
518c36f77d8SJed Brown   degree[0] = rp->test_mode ? 3 : 2;
51953a0f73bSJed Brown   ierr = PetscOptionsIntArray("-degree",
52053a0f73bSJed Brown                               "Polynomial degree of tensor product basis", NULL,
52153a0f73bSJed Brown                               degree, &num_degrees, &degree_set); CHKERRQ(ierr);
52253a0f73bSJed Brown   if (!degree_set)
52353a0f73bSJed Brown     num_degrees = 1;
5245f284d84SJed Brown   rp->qextra = PETSC_DECIDE;
5255f284d84SJed Brown   ierr = PetscOptionsInt("-qextra",
5265f284d84SJed Brown                          "Number of extra quadrature points (-1 for auto)", NULL,
527c36f77d8SJed Brown                          rp->qextra, &rp->qextra, NULL); CHKERRQ(ierr);
528565a3730SJed Brown   {
529565a3730SJed Brown     PetscBool set;
530565a3730SJed Brown     ierr = PetscOptionsStringArray("-ceed",
531565a3730SJed Brown                                    "CEED resource specifier (comma-separated list)", NULL,
532565a3730SJed Brown                                    ceedresources, &num_ceedresources, &set); CHKERRQ(ierr);
533565a3730SJed Brown     if (!set) {
534565a3730SJed Brown       ierr = PetscStrallocpy( "/cpu/self", &ceedresources[0]); CHKERRQ(ierr);
535565a3730SJed Brown       num_ceedresources = 1;
536565a3730SJed Brown     }
537565a3730SJed Brown   }
53853a0f73bSJed Brown   ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
53953a0f73bSJed Brown   ierr = PetscOptionsString("-hostname", "Hostname for output", NULL, hostname,
54053a0f73bSJed Brown                             hostname, sizeof(hostname), NULL); CHKERRQ(ierr);
541c36f77d8SJed Brown   rp->read_mesh = PETSC_FALSE;
54253a0f73bSJed Brown   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, filename,
543c36f77d8SJed Brown                             filename, sizeof(filename), &rp->read_mesh);
5440c59ef15SJeremy L Thompson   CHKERRQ(ierr);
545c36f77d8SJed Brown   rp->filename = filename;
546c36f77d8SJed Brown   if (!rp->read_mesh) {
547cb32e2e7SValeria Barra     PetscInt tmp = dim;
548cb32e2e7SValeria Barra     ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL,
549cb32e2e7SValeria Barra                                 melem, &tmp, NULL); CHKERRQ(ierr);
550cb32e2e7SValeria Barra   }
551c36f77d8SJed Brown   rp->memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
55253a0f73bSJed Brown   ierr = PetscOptionsEnum("-memtype", "CEED MemType requested", NULL, memTypes,
553c36f77d8SJed Brown                           (PetscEnum)rp->memtyperequested,
554c36f77d8SJed Brown                           (PetscEnum *)&rp->memtyperequested, &rp->setmemtyperequest);
5559396343dSjeremylt   CHKERRQ(ierr);
55653a0f73bSJed Brown   localnodes[0] = 1000;
55753a0f73bSJed Brown   ierr = PetscOptionsIntArray("-local_nodes",
55853a0f73bSJed Brown                               "Target number of locally owned nodes per "
55953a0f73bSJed Brown                               "process (single value or min,max)",
560c36f77d8SJed Brown                               NULL, localnodes, &num_localnodes, &rp->userlnodes);
561153bcb04Svaleriabarra   CHKERRQ(ierr);
56253a0f73bSJed Brown   if (num_localnodes < 2)
56353a0f73bSJed Brown     localnodes[1] = 2 * localnodes[0];
564c36f77d8SJed Brown   {
565c36f77d8SJed Brown     PetscInt two = 2;
566c36f77d8SJed Brown     rp->ksp_max_it_clip[0] = 5;
567c36f77d8SJed Brown     rp->ksp_max_it_clip[1] = 20;
568c36f77d8SJed Brown     ierr = PetscOptionsIntArray("-ksp_max_it_clip",
569c36f77d8SJed Brown                                 "Min and max number of iterations to use during benchmarking",
570c36f77d8SJed Brown                                 NULL, rp->ksp_max_it_clip, &two, NULL); CHKERRQ(ierr);
571c36f77d8SJed Brown   }
572da9108adSvaleriabarra   if (!degree_set) {
57353a0f73bSJed Brown     PetscInt maxdegree = 8;
5745f284d84SJed Brown     ierr = PetscOptionsInt("-max_degree",
5755f284d84SJed Brown                            "Range of degrees [1, maxdegree] to run with",
57653a0f73bSJed Brown                            NULL, maxdegree, &maxdegree, NULL);
5771e284482Svaleriabarra     CHKERRQ(ierr);
57853a0f73bSJed Brown     for (PetscInt i = 0; i < maxdegree; i++)
57953a0f73bSJed Brown       degree[i] = i + 1;
58053a0f73bSJed Brown     num_degrees = maxdegree;
58153a0f73bSJed Brown   }
58253a0f73bSJed Brown   {
58353a0f73bSJed Brown     PetscBool flg;
58453a0f73bSJed Brown     PetscInt p = rankspernode;
58553a0f73bSJed Brown     ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL,
58653a0f73bSJed Brown                            p, &p, &flg);
58753a0f73bSJed Brown     CHKERRQ(ierr);
58853a0f73bSJed Brown     if (flg) rankspernode = p;
58953a0f73bSJed Brown   }
5909396343dSjeremylt 
59153a0f73bSJed Brown   ierr = PetscOptionsEnd();
59253a0f73bSJed Brown   CHKERRQ(ierr);
593cb32e2e7SValeria Barra 
5941e284482Svaleriabarra   // Register PETSc logging stage
595c36f77d8SJed Brown   ierr = PetscLogStageRegister("Solve Stage", &rp->solvestage);
59653a0f73bSJed Brown   CHKERRQ(ierr);
59709a940d7Sjeremylt 
5981e284482Svaleriabarra   rp->petschavecuda = petschavecuda;
59953a0f73bSJed Brown   rp->hostname = hostname;
6001e284482Svaleriabarra   rp->dim = dim;
6011e284482Svaleriabarra   rp->melem = melem;
60253a0f73bSJed Brown   rp->rankspernode = rankspernode;
603cb32e2e7SValeria Barra 
60453a0f73bSJed Brown   for (PetscInt d = 0; d < num_degrees; d++) {
60553a0f73bSJed Brown     PetscInt deg = degree[d];
60653a0f73bSJed Brown     for (PetscInt n = localnodes[0]; n < localnodes[1]; n *= 2) {
60753a0f73bSJed Brown       rp->degree = deg;
60853a0f73bSJed Brown       rp->localnodes = n;
6095f284d84SJed Brown       ierr = Run(rp, num_ceedresources, ceedresources,
6105f284d84SJed Brown                  num_bpchoices, bpchoices); CHKERRQ(ierr);
611565a3730SJed Brown     }
612565a3730SJed Brown   }
6131e284482Svaleriabarra   // Clear memory
6141e284482Svaleriabarra   ierr = PetscFree(rp); CHKERRQ(ierr);
615565a3730SJed Brown   for (PetscInt i=0; i<num_ceedresources; i++) {
616565a3730SJed Brown     ierr = PetscFree(ceedresources[i]); CHKERRQ(ierr);
617565a3730SJed Brown   }
6180c59ef15SJeremy L Thompson   return PetscFinalize();
6190c59ef15SJeremy L Thompson }
620