xref: /libCEED/examples/petsc/bps.c (revision 565a373045b9458f5f5657066be084ecad2ef16b)
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 //
37cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3
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;
791e284482Svaleriabarra   PetscBool test_mode, benchmark_mode, read_mesh, userlnodes, setmemtyperequest,
801e284482Svaleriabarra             petschavecuda, write_solution;
8153a0f73bSJed Brown   char *filename, *ceedresource, *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 // -----------------------------------------------------------------------------
931e284482Svaleriabarra 
941e284482Svaleriabarra static PetscErrorCode Run(RunParams rp) {
951e284482Svaleriabarra   PetscInt 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   DM  dm;
103cb32e2e7SValeria Barra   UserO userO;
1040c59ef15SJeremy L Thompson   Ceed ceed;
105cb32e2e7SValeria Barra   CeedData ceeddata;
106c70bd2a0Sjeremylt   CeedQFunction qferror;
107c70bd2a0Sjeremylt   CeedOperator operror;
108cb32e2e7SValeria Barra   CeedVector rhsceed, target;
1091e284482Svaleriabarra 
1101e284482Svaleriabarra   PetscFunctionBeginUser;
1111e284482Svaleriabarra   // Setup DM
1121e284482Svaleriabarra   if (rp->read_mesh) {
1131e284482Svaleriabarra     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, PETSC_TRUE, &dm);
1141e284482Svaleriabarra     CHKERRQ(ierr);
1151e284482Svaleriabarra   } else {
1161e284482Svaleriabarra     if (rp->userlnodes) {
1171e284482Svaleriabarra       // Find a nicely composite number of elements no less than global nodes
1181e284482Svaleriabarra       PetscMPIInt size;
1191e284482Svaleriabarra       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
1201e284482Svaleriabarra       for (PetscInt gelem =
1211e284482Svaleriabarra              PetscMax(1, size * rp->localnodes / PetscPowInt(rp->degree, rp->dim));
1221e284482Svaleriabarra            ;
1231e284482Svaleriabarra            gelem++) {
1241e284482Svaleriabarra         Split3(gelem, rp->melem, true);
1251e284482Svaleriabarra         if (Max3(rp->melem) / Min3(rp->melem) <= 2) break;
1261e284482Svaleriabarra       }
1271e284482Svaleriabarra     }
1281e284482Svaleriabarra     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, PETSC_FALSE, rp->melem,
1291e284482Svaleriabarra                                NULL, NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
1301e284482Svaleriabarra   }
1311e284482Svaleriabarra 
1321e284482Svaleriabarra   {
1331e284482Svaleriabarra     DM dmDist = NULL;
1341e284482Svaleriabarra     PetscPartitioner part;
1351e284482Svaleriabarra 
1361e284482Svaleriabarra     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1371e284482Svaleriabarra     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1381e284482Svaleriabarra     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1391e284482Svaleriabarra     if (dmDist) {
1401e284482Svaleriabarra       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1411e284482Svaleriabarra       dm  = dmDist;
1421e284482Svaleriabarra     }
1431e284482Svaleriabarra   }
1441e284482Svaleriabarra 
1451e284482Svaleriabarra   // Set up libCEED
1461e284482Svaleriabarra   CeedInit(rp->ceedresource, &ceed);
1471e284482Svaleriabarra   CeedMemType memtypebackend;
1481e284482Svaleriabarra   CeedGetPreferredMemType(ceed, &memtypebackend);
1491e284482Svaleriabarra 
1501e284482Svaleriabarra   // Check memtype compatibility
1511e284482Svaleriabarra   if (!rp->setmemtyperequest)
1521e284482Svaleriabarra     rp->memtyperequested = memtypebackend;
1531e284482Svaleriabarra   else if (!rp->petschavecuda && rp->memtyperequested == CEED_MEM_DEVICE)
1541e284482Svaleriabarra     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
1551e284482Svaleriabarra              "PETSc was not built with CUDA. "
1561e284482Svaleriabarra              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
1571e284482Svaleriabarra 
1581e284482Svaleriabarra   // Create DM
1591e284482Svaleriabarra   ierr = SetupDMByDegree(dm, rp->degree, rp->ncompu, rp->bpchoice);
1601e284482Svaleriabarra   CHKERRQ(ierr);
1611e284482Svaleriabarra 
1621e284482Svaleriabarra   // Create vectors
1631e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_DEVICE) {
1641e284482Svaleriabarra     ierr = DMSetVecType(dm, VECCUDA); CHKERRQ(ierr);
1651e284482Svaleriabarra   }
1661e284482Svaleriabarra   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
1671e284482Svaleriabarra   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
1681e284482Svaleriabarra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
1691e284482Svaleriabarra   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
1701e284482Svaleriabarra   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
1711e284482Svaleriabarra   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
1721e284482Svaleriabarra 
1731e284482Svaleriabarra   // Operator
1741e284482Svaleriabarra   ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr);
1751e284482Svaleriabarra   ierr = MatCreateShell(rp->comm, lsize, lsize, gsize, gsize,
1761e284482Svaleriabarra                         userO, &matO); CHKERRQ(ierr);
1771e284482Svaleriabarra   ierr = MatShellSetOperation(matO, MATOP_MULT,
1781e284482Svaleriabarra                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
1791e284482Svaleriabarra   ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL,
1801e284482Svaleriabarra                               (void(*)(void))MatGetDiag); CHKERRQ(ierr);
1811e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_DEVICE) {
1821e284482Svaleriabarra     ierr = MatShellSetVecType(matO, VECCUDA); CHKERRQ(ierr);
1831e284482Svaleriabarra   }
1841e284482Svaleriabarra 
1851e284482Svaleriabarra   // Print summary
1861e284482Svaleriabarra   if (!rp->test_mode) {
1871e284482Svaleriabarra     PetscInt P = rp->degree + 1, Q = P + rp->qextra;
1881e284482Svaleriabarra 
1891e284482Svaleriabarra     const char *usedresource;
1901e284482Svaleriabarra     CeedGetResource(ceed, &usedresource);
1911e284482Svaleriabarra 
1921e284482Svaleriabarra     VecType vectype;
1931e284482Svaleriabarra     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
1941e284482Svaleriabarra 
1951e284482Svaleriabarra     PetscInt cStart, cEnd;
1961e284482Svaleriabarra     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr);
19753a0f73bSJed Brown     PetscMPIInt comm_size;
19853a0f73bSJed Brown     ierr = MPI_Comm_size(rp->comm, &comm_size); CHKERRQ(ierr);
1991e284482Svaleriabarra     ierr = PetscPrintf(rp->comm,
2001e284482Svaleriabarra                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
20153a0f73bSJed Brown                        "  MPI:\n"
20253a0f73bSJed Brown                        "    Hostname                           : %s\n"
20353a0f73bSJed Brown                        "    Total ranks                        : %d\n"
20453a0f73bSJed Brown                        "    Ranks per compute node             : %d\n"
2051e284482Svaleriabarra                        "  PETSc:\n"
2061e284482Svaleriabarra                        "    PETSc Vec Type                     : %s\n"
2071e284482Svaleriabarra                        "  libCEED:\n"
2081e284482Svaleriabarra                        "    libCEED Backend                    : %s\n"
2091e284482Svaleriabarra                        "    libCEED Backend MemType            : %s\n"
2101e284482Svaleriabarra                        "    libCEED User Requested MemType     : %s\n"
2111e284482Svaleriabarra                        "  Mesh:\n"
2121e284482Svaleriabarra                        "    Number of 1D Basis Nodes (P)       : %d\n"
2131e284482Svaleriabarra                        "    Number of 1D Quadrature Points (Q) : %d\n"
2141e284482Svaleriabarra                        "    Global nodes                       : %D\n"
2151e284482Svaleriabarra                        "    Local Elements                     : %D\n"
2161e284482Svaleriabarra                        "    Owned nodes                        : %D\n"
2171e284482Svaleriabarra                        "    DoF per node                       : %D\n",
21853a0f73bSJed Brown                        rp->bpchoice+1,
21953a0f73bSJed Brown                        rp->hostname,
22053a0f73bSJed Brown                        comm_size, rp->rankspernode,
22153a0f73bSJed Brown                        vectype, usedresource,
2221e284482Svaleriabarra                        CeedMemTypes[memtypebackend],
2231e284482Svaleriabarra                        (rp->setmemtyperequest) ?
2241e284482Svaleriabarra                        CeedMemTypes[rp->memtyperequested] : "none",
2251e284482Svaleriabarra                        P, Q, gsize/rp->ncompu, cEnd - cStart, lsize/rp->ncompu,
2261e284482Svaleriabarra                        rp->ncompu);
2271e284482Svaleriabarra     CHKERRQ(ierr);
2281e284482Svaleriabarra   }
2291e284482Svaleriabarra 
2301e284482Svaleriabarra   // Create RHS vector
2311e284482Svaleriabarra   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
2321e284482Svaleriabarra   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
2331e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_HOST) {
2341e284482Svaleriabarra     ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
2351e284482Svaleriabarra   } else {
2361e284482Svaleriabarra     ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr);
2371e284482Svaleriabarra   }
2381e284482Svaleriabarra   CeedVectorCreate(ceed, xlsize, &rhsceed);
2391e284482Svaleriabarra   CeedVectorSetArray(rhsceed, rp->memtyperequested, CEED_USE_POINTER, r);
2401e284482Svaleriabarra 
2411e284482Svaleriabarra   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
2421e284482Svaleriabarra   ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->qextra,
2431e284482Svaleriabarra                               rp->ncompu, gsize, xlsize, rp->bpchoice, ceeddata,
2441e284482Svaleriabarra                               true, rhsceed, &target); CHKERRQ(ierr);
2451e284482Svaleriabarra 
2461e284482Svaleriabarra   // Gather RHS
2471e284482Svaleriabarra   CeedVectorSyncArray(rhsceed, rp->memtyperequested);
2481e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_HOST) {
2491e284482Svaleriabarra     ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
2501e284482Svaleriabarra   } else {
2511e284482Svaleriabarra     ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr);
2521e284482Svaleriabarra   }
2531e284482Svaleriabarra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
2541e284482Svaleriabarra   ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
2551e284482Svaleriabarra   ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
2561e284482Svaleriabarra   CeedVectorDestroy(&rhsceed);
2571e284482Svaleriabarra 
2581e284482Svaleriabarra   // Create the error Q-function
2591e284482Svaleriabarra   CeedQFunctionCreateInterior(ceed, 1, bpOptions[rp->bpchoice].error,
2601e284482Svaleriabarra                               bpOptions[rp->bpchoice].errorfname, &qferror);
2611e284482Svaleriabarra   CeedQFunctionAddInput(qferror, "u", rp->ncompu, CEED_EVAL_INTERP);
2621e284482Svaleriabarra   CeedQFunctionAddInput(qferror, "true_soln", rp->ncompu, CEED_EVAL_NONE);
2631e284482Svaleriabarra   CeedQFunctionAddOutput(qferror, "error", rp->ncompu, CEED_EVAL_NONE);
2641e284482Svaleriabarra 
2651e284482Svaleriabarra   // Create the error operator
2661e284482Svaleriabarra   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
2671e284482Svaleriabarra                      &operror);
2681e284482Svaleriabarra   CeedOperatorSetField(operror, "u", ceeddata->Erestrictu,
2691e284482Svaleriabarra                        ceeddata->basisu, CEED_VECTOR_ACTIVE);
2701e284482Svaleriabarra   CeedOperatorSetField(operror, "true_soln", ceeddata->Erestrictui,
2711e284482Svaleriabarra                        CEED_BASIS_COLLOCATED, target);
2721e284482Svaleriabarra   CeedOperatorSetField(operror, "error", ceeddata->Erestrictui,
2731e284482Svaleriabarra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
2741e284482Svaleriabarra 
2751e284482Svaleriabarra   // Set up Mat
2761e284482Svaleriabarra   userO->comm = rp->comm;
2771e284482Svaleriabarra   userO->dm = dm;
2781e284482Svaleriabarra   userO->Xloc = Xloc;
2791e284482Svaleriabarra   ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr);
2801e284482Svaleriabarra   userO->xceed = ceeddata->xceed;
2811e284482Svaleriabarra   userO->yceed = ceeddata->yceed;
2821e284482Svaleriabarra   userO->op = ceeddata->opapply;
2831e284482Svaleriabarra   userO->ceed = ceed;
2841e284482Svaleriabarra   userO->memtype = rp->memtyperequested;
2851e284482Svaleriabarra   if (rp->memtyperequested == CEED_MEM_HOST) {
2861e284482Svaleriabarra     userO->VecGetArray = VecGetArray;
2871e284482Svaleriabarra     userO->VecGetArrayRead = VecGetArrayRead;
2881e284482Svaleriabarra     userO->VecRestoreArray = VecRestoreArray;
2891e284482Svaleriabarra     userO->VecRestoreArrayRead = VecRestoreArrayRead;
2901e284482Svaleriabarra   } else {
2911e284482Svaleriabarra     userO->VecGetArray = VecCUDAGetArray;
2921e284482Svaleriabarra     userO->VecGetArrayRead = VecCUDAGetArrayRead;
2931e284482Svaleriabarra     userO->VecRestoreArray = VecCUDARestoreArray;
2941e284482Svaleriabarra     userO->VecRestoreArrayRead = VecCUDARestoreArrayRead;
2951e284482Svaleriabarra   }
2961e284482Svaleriabarra 
2971e284482Svaleriabarra   ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr);
2981e284482Svaleriabarra   {
2991e284482Svaleriabarra     PC pc;
3001e284482Svaleriabarra     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
3011e284482Svaleriabarra     if (rp->bpchoice == CEED_BP1 || rp->bpchoice == CEED_BP2) {
3021e284482Svaleriabarra       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
3031e284482Svaleriabarra       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
3041e284482Svaleriabarra     } else {
3051e284482Svaleriabarra       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
3061e284482Svaleriabarra     }
3071e284482Svaleriabarra     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
3081e284482Svaleriabarra     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
3091e284482Svaleriabarra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
3101e284482Svaleriabarra                             PETSC_DEFAULT); CHKERRQ(ierr);
3111e284482Svaleriabarra   }
3121e284482Svaleriabarra   ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr);
3131e284482Svaleriabarra 
3141e284482Svaleriabarra   // First run, if benchmarking
3151e284482Svaleriabarra   if (rp->benchmark_mode) {
3161e284482Svaleriabarra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
3171e284482Svaleriabarra     CHKERRQ(ierr);
3181e284482Svaleriabarra     my_rt_start = MPI_Wtime();
3191e284482Svaleriabarra     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
3201e284482Svaleriabarra     my_rt = MPI_Wtime() - my_rt_start;
3211e284482Svaleriabarra     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
3221e284482Svaleriabarra     CHKERRQ(ierr);
3231e284482Svaleriabarra     // Set maxits based on first iteration timing
3241e284482Svaleriabarra     if (my_rt > 0.02) {
325*565a3730SJed Brown       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
326*565a3730SJed Brown                               rp->ksp_max_it_clip[0]);
3271e284482Svaleriabarra       CHKERRQ(ierr);
3281e284482Svaleriabarra     } else {
329*565a3730SJed Brown       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
330*565a3730SJed Brown                               rp->ksp_max_it_clip[1]);
3311e284482Svaleriabarra       CHKERRQ(ierr);
3321e284482Svaleriabarra     }
3331e284482Svaleriabarra   }
3341e284482Svaleriabarra   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
3351e284482Svaleriabarra 
3361e284482Svaleriabarra   // Timed solve
3371e284482Svaleriabarra   ierr = VecZeroEntries(X); CHKERRQ(ierr);
3381e284482Svaleriabarra   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
3391e284482Svaleriabarra 
3401e284482Svaleriabarra   // -- Performance logging
3411e284482Svaleriabarra   ierr = PetscLogStagePush(rp->solvestage); CHKERRQ(ierr);
3421e284482Svaleriabarra 
3431e284482Svaleriabarra   // -- Solve
3441e284482Svaleriabarra   my_rt_start = MPI_Wtime();
3451e284482Svaleriabarra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
3461e284482Svaleriabarra   my_rt = MPI_Wtime() - my_rt_start;
3471e284482Svaleriabarra 
3481e284482Svaleriabarra   // -- Performance logging
3491e284482Svaleriabarra   ierr = PetscLogStagePop();
3501e284482Svaleriabarra 
3511e284482Svaleriabarra   // Output results
3521e284482Svaleriabarra   {
3531e284482Svaleriabarra     KSPType ksptype;
3541e284482Svaleriabarra     KSPConvergedReason reason;
3551e284482Svaleriabarra     PetscReal rnorm;
3561e284482Svaleriabarra     PetscInt its;
3571e284482Svaleriabarra     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
3581e284482Svaleriabarra     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
3591e284482Svaleriabarra     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
3601e284482Svaleriabarra     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
3611e284482Svaleriabarra     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
3621e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,
3631e284482Svaleriabarra                          "  KSP:\n"
3641e284482Svaleriabarra                          "    KSP Type                           : %s\n"
3651e284482Svaleriabarra                          "    KSP Convergence                    : %s\n"
3661e284482Svaleriabarra                          "    Total KSP Iterations               : %D\n"
3671e284482Svaleriabarra                          "    Final rnorm                        : %e\n",
3681e284482Svaleriabarra                          ksptype, KSPConvergedReasons[reason], its,
3691e284482Svaleriabarra                          (double)rnorm); CHKERRQ(ierr);
3701e284482Svaleriabarra     }
3711e284482Svaleriabarra     if (!rp->test_mode) {
3721e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,"  Performance:\n"); CHKERRQ(ierr);
3731e284482Svaleriabarra     }
3741e284482Svaleriabarra     {
3751e284482Svaleriabarra       PetscReal maxerror;
3761e284482Svaleriabarra       ierr = ComputeErrorMax(userO, operror, X, target, &maxerror);
3771e284482Svaleriabarra       CHKERRQ(ierr);
3781e284482Svaleriabarra       PetscReal tol = 5e-2;
3791e284482Svaleriabarra       if (!rp->test_mode || maxerror > tol) {
3801e284482Svaleriabarra         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm);
3811e284482Svaleriabarra         CHKERRQ(ierr);
3821e284482Svaleriabarra         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm);
3831e284482Svaleriabarra         CHKERRQ(ierr);
3841e284482Svaleriabarra         ierr = PetscPrintf(rp->comm,
3851e284482Svaleriabarra                            "    Pointwise Error (max)              : %e\n"
3861e284482Svaleriabarra                            "    CG Solve Time                      : %g (%g) sec\n",
3871e284482Svaleriabarra                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
3881e284482Svaleriabarra       }
3891e284482Svaleriabarra     }
3901e284482Svaleriabarra     if (rp->benchmark_mode && (!rp->test_mode)) {
3911e284482Svaleriabarra       ierr = PetscPrintf(rp->comm,
3921e284482Svaleriabarra                          "    DoFs/Sec in CG                     : %g (%g) million\n",
3931e284482Svaleriabarra                          1e-6*gsize*its/rt_max,
3941e284482Svaleriabarra                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
3951e284482Svaleriabarra     }
3961e284482Svaleriabarra   }
3971e284482Svaleriabarra 
3981e284482Svaleriabarra   if (rp->write_solution) {
3991e284482Svaleriabarra     PetscViewer vtkviewersoln;
4001e284482Svaleriabarra 
4011e284482Svaleriabarra     ierr = PetscViewerCreate(rp->comm, &vtkviewersoln); CHKERRQ(ierr);
4021e284482Svaleriabarra     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
4031e284482Svaleriabarra     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
4041e284482Svaleriabarra     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
4051e284482Svaleriabarra     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
4061e284482Svaleriabarra   }
4071e284482Svaleriabarra 
4081e284482Svaleriabarra   // Cleanup
4091e284482Svaleriabarra   ierr = VecDestroy(&X); CHKERRQ(ierr);
4101e284482Svaleriabarra   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
4111e284482Svaleriabarra   ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr);
4121e284482Svaleriabarra   ierr = MatDestroy(&matO); CHKERRQ(ierr);
4131e284482Svaleriabarra   ierr = PetscFree(userO); CHKERRQ(ierr);
4141e284482Svaleriabarra   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
4151e284482Svaleriabarra   ierr = DMDestroy(&dm); CHKERRQ(ierr);
4161e284482Svaleriabarra 
4171e284482Svaleriabarra   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
4181e284482Svaleriabarra   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
4191e284482Svaleriabarra   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
4201e284482Svaleriabarra   CeedVectorDestroy(&target);
4211e284482Svaleriabarra   CeedQFunctionDestroy(&qferror);
4221e284482Svaleriabarra   CeedOperatorDestroy(&operror);
4231e284482Svaleriabarra   CeedDestroy(&ceed);
4241e284482Svaleriabarra   PetscFunctionReturn(0);
4251e284482Svaleriabarra }
4261e284482Svaleriabarra 
4271e284482Svaleriabarra int main(int argc, char **argv) {
4281e284482Svaleriabarra   PetscInt ierr, commsize;
4291e284482Svaleriabarra   RunParams rp;
43053a0f73bSJed Brown   MPI_Comm comm;
43153a0f73bSJed Brown   char filename[PETSC_MAX_PATH_LEN];
432*565a3730SJed Brown   char *ceedresources[30];
433*565a3730SJed Brown   PetscInt num_ceedresources = 30;
43453a0f73bSJed Brown   char hostname[PETSC_MAX_PATH_LEN];
4351e284482Svaleriabarra 
436c36f77d8SJed Brown   PetscInt dim = 3, melem[3] = {3, 3, 3};
43753a0f73bSJed Brown   PetscInt num_degrees = 30, degree[30] = {}, num_localnodes = 2, localnodes[2] = {};
43853a0f73bSJed Brown   PetscMPIInt rankspernode;
439c36f77d8SJed Brown   PetscBool degree_set;
440*565a3730SJed Brown   bpType bpchoices[10];
441*565a3730SJed Brown   PetscInt num_bpchoices = 10;
4420c59ef15SJeremy L Thompson 
4439396343dSjeremylt   // Check PETSc CUDA support
444c36f77d8SJed Brown   PetscBool petschavecuda;
4459396343dSjeremylt   // *INDENT-OFF*
4469396343dSjeremylt   #ifdef PETSC_HAVE_CUDA
4479396343dSjeremylt   petschavecuda = PETSC_TRUE;
4489396343dSjeremylt   #else
4499396343dSjeremylt   petschavecuda = PETSC_FALSE;
4509396343dSjeremylt   #endif
4519396343dSjeremylt   // *INDENT-ON*
4529396343dSjeremylt 
4531e284482Svaleriabarra   // Initialize PETSc
4540c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
4550c59ef15SJeremy L Thompson   if (ierr) return ierr;
4560c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
4571e284482Svaleriabarra   ierr = MPI_Comm_size(comm, &commsize);
4581e284482Svaleriabarra   if (ierr != MPI_SUCCESS) return ierr;
45953a0f73bSJed Brown   #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
46053a0f73bSJed Brown   {
46153a0f73bSJed Brown     MPI_Comm splitcomm;
4621e284482Svaleriabarra     ierr = MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL,
4631e284482Svaleriabarra                                &splitcomm);
46453a0f73bSJed Brown     CHKERRQ(ierr);
46553a0f73bSJed Brown     ierr = MPI_Comm_size(splitcomm, &rankspernode); CHKERRQ(ierr);
46653a0f73bSJed Brown     ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr);
46753a0f73bSJed Brown   }
46853a0f73bSJed Brown   #else
46953a0f73bSJed Brown   rankspernode = -1; // Unknown
47053a0f73bSJed Brown   #endif
471cb32e2e7SValeria Barra 
472c36f77d8SJed Brown   // Setup all parameters needed in Run()
473c36f77d8SJed Brown   ierr = PetscMalloc1(1, &rp); CHKERRQ(ierr);
474c36f77d8SJed Brown   rp->comm = comm;
475c36f77d8SJed Brown 
47632d2ee49SValeria Barra   // Read command line options
4771e284482Svaleriabarra   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
4781e284482Svaleriabarra   CHKERRQ(ierr);
479*565a3730SJed Brown   {
480*565a3730SJed Brown     PetscBool set;
481*565a3730SJed Brown     ierr = PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve",
482*565a3730SJed Brown                                  NULL,
483*565a3730SJed Brown                                  bpTypes, (PetscEnum *)bpchoices, &num_bpchoices, &set);
484*565a3730SJed Brown     CHKERRQ(ierr);
485*565a3730SJed Brown     if (!set) {
486*565a3730SJed Brown       bpchoices[0] = CEED_BP1;
487*565a3730SJed Brown       num_bpchoices = 1;
488*565a3730SJed Brown     }
489*565a3730SJed Brown   }
490c36f77d8SJed Brown   rp->test_mode = PETSC_FALSE;
4910c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
4920c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
493c36f77d8SJed Brown                           NULL, rp->test_mode, &rp->test_mode, NULL); CHKERRQ(ierr);
494c36f77d8SJed Brown   rp->benchmark_mode = PETSC_FALSE;
4950c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-benchmark",
4960c59ef15SJeremy L Thompson                           "Benchmarking mode (prints benchmark statistics)",
497c36f77d8SJed Brown                           NULL, rp->benchmark_mode, &rp->benchmark_mode, NULL);
4980c59ef15SJeremy L Thompson   CHKERRQ(ierr);
499c36f77d8SJed Brown   rp->write_solution = PETSC_FALSE;
50053a0f73bSJed Brown   ierr = PetscOptionsBool("-write_solution", "Write solution for visualization",
501c36f77d8SJed Brown                           NULL, rp->write_solution, &rp->write_solution, NULL);
5026bd9afcaSjeremylt   CHKERRQ(ierr);
503c36f77d8SJed Brown   degree[0] = rp->test_mode ? 3 : 2;
50453a0f73bSJed Brown   ierr = PetscOptionsIntArray("-degree",
50553a0f73bSJed Brown                               "Polynomial degree of tensor product basis", NULL,
50653a0f73bSJed Brown                               degree, &num_degrees, &degree_set); CHKERRQ(ierr);
50753a0f73bSJed Brown   if (!degree_set)
50853a0f73bSJed Brown     num_degrees = 1;
50953a0f73bSJed Brown   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", NULL,
510c36f77d8SJed Brown                          rp->qextra, &rp->qextra, NULL); CHKERRQ(ierr);
511*565a3730SJed Brown   {
512*565a3730SJed Brown     PetscBool set;
513*565a3730SJed Brown     ierr = PetscOptionsStringArray("-ceed",
514*565a3730SJed Brown                                    "CEED resource specifier (comma-separated list)", NULL,
515*565a3730SJed Brown                                    ceedresources, &num_ceedresources, &set); CHKERRQ(ierr);
516*565a3730SJed Brown     if (!set) {
517*565a3730SJed Brown       ierr = PetscStrallocpy( "/cpu/self", &ceedresources[0]); CHKERRQ(ierr);
518*565a3730SJed Brown       num_ceedresources = 1;
519*565a3730SJed Brown     }
520*565a3730SJed Brown   }
52153a0f73bSJed Brown   ierr = PetscGetHostName(hostname, sizeof hostname); CHKERRQ(ierr);
52253a0f73bSJed Brown   ierr = PetscOptionsString("-hostname", "Hostname for output", NULL, hostname,
52353a0f73bSJed Brown                             hostname, sizeof(hostname), NULL); CHKERRQ(ierr);
524c36f77d8SJed Brown   rp->read_mesh = PETSC_FALSE;
52553a0f73bSJed Brown   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, filename,
526c36f77d8SJed Brown                             filename, sizeof(filename), &rp->read_mesh);
5270c59ef15SJeremy L Thompson   CHKERRQ(ierr);
528c36f77d8SJed Brown   rp->filename = filename;
529c36f77d8SJed Brown   if (!rp->read_mesh) {
530cb32e2e7SValeria Barra     PetscInt tmp = dim;
531cb32e2e7SValeria Barra     ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL,
532cb32e2e7SValeria Barra                                 melem, &tmp, NULL); CHKERRQ(ierr);
533cb32e2e7SValeria Barra   }
534c36f77d8SJed Brown   rp->memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
53553a0f73bSJed Brown   ierr = PetscOptionsEnum("-memtype", "CEED MemType requested", NULL, memTypes,
536c36f77d8SJed Brown                           (PetscEnum)rp->memtyperequested,
537c36f77d8SJed Brown                           (PetscEnum *)&rp->memtyperequested, &rp->setmemtyperequest);
5389396343dSjeremylt   CHKERRQ(ierr);
539c36f77d8SJed Brown 
54053a0f73bSJed Brown   localnodes[0] = 1000;
54153a0f73bSJed Brown   ierr = PetscOptionsIntArray("-local_nodes",
54253a0f73bSJed Brown                               "Target number of locally owned nodes per "
54353a0f73bSJed Brown                               "process (single value or min,max)",
544c36f77d8SJed Brown                               NULL, localnodes, &num_localnodes, &rp->userlnodes);
545153bcb04Svaleriabarra   CHKERRQ(ierr);
54653a0f73bSJed Brown   if (num_localnodes < 2)
54753a0f73bSJed Brown     localnodes[1] = 2 * localnodes[0];
548c36f77d8SJed Brown   {
549c36f77d8SJed Brown     PetscInt two = 2;
550c36f77d8SJed Brown     rp->ksp_max_it_clip[0] = 5;
551c36f77d8SJed Brown     rp->ksp_max_it_clip[1] = 20;
552c36f77d8SJed Brown     ierr = PetscOptionsIntArray("-ksp_max_it_clip",
553c36f77d8SJed Brown                                 "Min and max number of iterations to use during benchmarking",
554c36f77d8SJed Brown                                 NULL, rp->ksp_max_it_clip, &two, NULL); CHKERRQ(ierr);
555c36f77d8SJed Brown   }
556c36f77d8SJed Brown 
557c36f77d8SJed Brown   if (rp->benchmark_mode && !degree_set) {
55853a0f73bSJed Brown     PetscInt maxdegree = 8;
55953a0f73bSJed Brown     ierr = PetscOptionsInt("-max_degree", "Max degree to run in benchmark mode",
56053a0f73bSJed Brown                            NULL, maxdegree, &maxdegree, NULL);
5611e284482Svaleriabarra     CHKERRQ(ierr);
56253a0f73bSJed Brown     for (PetscInt i = 0; i < maxdegree; i++)
56353a0f73bSJed Brown       degree[i] = i + 1;
56453a0f73bSJed Brown     num_degrees = maxdegree;
56553a0f73bSJed Brown   }
56653a0f73bSJed Brown   {
56753a0f73bSJed Brown     PetscBool flg;
56853a0f73bSJed Brown     PetscInt p = rankspernode;
56953a0f73bSJed Brown     ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL,
57053a0f73bSJed Brown                            p, &p, &flg);
57153a0f73bSJed Brown     CHKERRQ(ierr);
57253a0f73bSJed Brown     if (flg) rankspernode = p;
57353a0f73bSJed Brown   }
5749396343dSjeremylt 
57553a0f73bSJed Brown   ierr = PetscOptionsEnd();
57653a0f73bSJed Brown   CHKERRQ(ierr);
577cb32e2e7SValeria Barra 
5781e284482Svaleriabarra   // Register PETSc logging stage
579c36f77d8SJed Brown   ierr = PetscLogStageRegister("Solve Stage", &rp->solvestage);
58053a0f73bSJed Brown   CHKERRQ(ierr);
58109a940d7Sjeremylt 
5821e284482Svaleriabarra   rp->petschavecuda = petschavecuda;
58353a0f73bSJed Brown   rp->hostname = hostname;
5841e284482Svaleriabarra   rp->dim = dim;
5851e284482Svaleriabarra   rp->melem = melem;
58653a0f73bSJed Brown   rp->rankspernode = rankspernode;
587cb32e2e7SValeria Barra 
588*565a3730SJed Brown   for (PetscInt b = 0; b < num_bpchoices; b++) {
589*565a3730SJed Brown     rp->bpchoice = bpchoices[b];
590*565a3730SJed Brown     rp->ncompu = bpOptions[rp->bpchoice].ncompu;
591*565a3730SJed Brown     rp->qextra = bpOptions[rp->bpchoice].qextra;
59253a0f73bSJed Brown     for (PetscInt d = 0; d < num_degrees; d++) {
59353a0f73bSJed Brown       PetscInt deg = degree[d];
59453a0f73bSJed Brown       for (PetscInt n = localnodes[0]; n < localnodes[1]; n *= 2) {
59553a0f73bSJed Brown         rp->degree = deg;
59653a0f73bSJed Brown         rp->localnodes = n;
597*565a3730SJed Brown         for (PetscInt c = 0; c < num_ceedresources; c++) {
598*565a3730SJed Brown           rp->ceedresource = ceedresources[c];
599*565a3730SJed Brown           ierr = Run(rp);
600*565a3730SJed Brown           CHKERRQ(ierr);
6010c59ef15SJeremy L Thompson         }
6020c59ef15SJeremy L Thompson       }
603*565a3730SJed Brown     }
604*565a3730SJed Brown   }
6051e284482Svaleriabarra   // Clear memory
6061e284482Svaleriabarra   ierr = PetscFree(rp); CHKERRQ(ierr);
607*565a3730SJed Brown   for (PetscInt i=0; i<num_ceedresources; i++) {
608*565a3730SJed Brown     ierr = PetscFree(ceedresources[i]); CHKERRQ(ierr);
609*565a3730SJed Brown   }
6100c59ef15SJeremy L Thompson   return PetscFinalize();
6110c59ef15SJeremy L Thompson }
612