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> 49*e83e87a5Sjeremylt #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) { 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; 79*e83e87a5Sjeremylt PetscBool test_mode, read_mesh, userlnodes, write_solution; 805f284d84SJed Brown char *filename, *hostname; 811e284482Svaleriabarra PetscInt localnodes, degree, qextra, dim, ncompu, *melem; 82c36f77d8SJed Brown PetscInt ksp_max_it_clip[2]; 8353a0f73bSJed Brown PetscMPIInt rankspernode; 841e284482Svaleriabarra bpType bpchoice; 8509a940d7Sjeremylt PetscLogStage solvestage; 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, 925f284d84SJed Brown const char *ceedresource) { 935f284d84SJed Brown PetscErrorCode ierr; 941e284482Svaleriabarra double my_rt_start, my_rt, rt_min, rt_max; 951e284482Svaleriabarra PetscInt xlsize, lsize, gsize; 961e284482Svaleriabarra PetscScalar *r; 97819eb1b3Sjeremylt Vec X, Xloc, rhs, rhsloc; 98cb32e2e7SValeria Barra Mat matO; 99819eb1b3Sjeremylt KSP ksp; 100cb32e2e7SValeria Barra UserO userO; 1010c59ef15SJeremy L Thompson Ceed ceed; 102cb32e2e7SValeria Barra CeedData ceeddata; 103*e83e87a5Sjeremylt CeedQFunction qfError; 104*e83e87a5Sjeremylt CeedOperator opError; 105cb32e2e7SValeria Barra CeedVector rhsceed, target; 106b68a8d79SJed Brown VecType vectype; 107b68a8d79SJed Brown PetscMemType memtype; 1081e284482Svaleriabarra 1091e284482Svaleriabarra PetscFunctionBeginUser; 1101e284482Svaleriabarra // Set up libCEED 1115f284d84SJed Brown CeedInit(ceedresource, &ceed); 1121e284482Svaleriabarra CeedMemType memtypebackend; 1131e284482Svaleriabarra CeedGetPreferredMemType(ceed, &memtypebackend); 1141e284482Svaleriabarra 115b68a8d79SJed Brown ierr = DMGetVecType(dm, &vectype); CHKERRQ(ierr); 116b68a8d79SJed Brown if (!vectype) { // Not yet set by user -dm_vec_type 117b68a8d79SJed Brown switch (memtypebackend) { 118b68a8d79SJed Brown case CEED_MEM_HOST: vectype = VECSTANDARD; break; 119b68a8d79SJed Brown case CEED_MEM_DEVICE: { 120b68a8d79SJed Brown const char *resolved; 121b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 122b68a8d79SJed Brown if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 123b68a8d79SJed Brown else if (strstr(resolved, "/gpu/hip/occa")) 124b68a8d79SJed Brown vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 125b68a8d79SJed Brown else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 126b68a8d79SJed Brown else vectype = VECSTANDARD; 1271e284482Svaleriabarra } 128b68a8d79SJed Brown } 129b68a8d79SJed Brown ierr = DMSetVecType(dm, vectype); CHKERRQ(ierr); 130b68a8d79SJed Brown } 131b68a8d79SJed Brown 132b68a8d79SJed Brown // Create global and local solution vectors 1331e284482Svaleriabarra ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr); 1341e284482Svaleriabarra ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr); 1351e284482Svaleriabarra ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 1361e284482Svaleriabarra ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr); 1371e284482Svaleriabarra ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr); 1381e284482Svaleriabarra ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr); 1391e284482Svaleriabarra 1401e284482Svaleriabarra // Operator 1411e284482Svaleriabarra ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr); 1421e284482Svaleriabarra ierr = MatCreateShell(rp->comm, lsize, lsize, gsize, gsize, 1431e284482Svaleriabarra userO, &matO); CHKERRQ(ierr); 1441e284482Svaleriabarra ierr = MatShellSetOperation(matO, MATOP_MULT, 1451e284482Svaleriabarra (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 1461e284482Svaleriabarra ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL, 1471e284482Svaleriabarra (void(*)(void))MatGetDiag); CHKERRQ(ierr); 148b68a8d79SJed Brown ierr = MatShellSetVecType(matO, vectype); CHKERRQ(ierr); 1491e284482Svaleriabarra 1501e284482Svaleriabarra // Print summary 1511e284482Svaleriabarra if (!rp->test_mode) { 1521e284482Svaleriabarra PetscInt P = rp->degree + 1, Q = P + rp->qextra; 1531e284482Svaleriabarra 1541e284482Svaleriabarra const char *usedresource; 1551e284482Svaleriabarra CeedGetResource(ceed, &usedresource); 1561e284482Svaleriabarra 1571e284482Svaleriabarra VecType vectype; 1581e284482Svaleriabarra ierr = VecGetType(X, &vectype); CHKERRQ(ierr); 1591e284482Svaleriabarra 1601e284482Svaleriabarra PetscInt cStart, cEnd; 1611e284482Svaleriabarra ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); 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", 18253a0f73bSJed Brown rp->bpchoice+1, 18353a0f73bSJed Brown rp->hostname, 18453a0f73bSJed Brown comm_size, rp->rankspernode, 18553a0f73bSJed Brown vectype, usedresource, 1861e284482Svaleriabarra CeedMemTypes[memtypebackend], 1871e284482Svaleriabarra P, Q, gsize/rp->ncompu, cEnd - cStart, lsize/rp->ncompu, 1881e284482Svaleriabarra rp->ncompu); 1891e284482Svaleriabarra CHKERRQ(ierr); 1901e284482Svaleriabarra } 1911e284482Svaleriabarra 1921e284482Svaleriabarra // Create RHS vector 1931e284482Svaleriabarra ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 1941e284482Svaleriabarra ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 195b68a8d79SJed Brown ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr); 1961e284482Svaleriabarra CeedVectorCreate(ceed, xlsize, &rhsceed); 197b68a8d79SJed Brown CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r); 1981e284482Svaleriabarra 1991e284482Svaleriabarra ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr); 2001e284482Svaleriabarra ierr = SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->qextra, 201*e83e87a5Sjeremylt rp->dim, rp->ncompu, gsize, xlsize, bpOptions[rp->bpchoice], 202*e83e87a5Sjeremylt ceeddata, true, rhsceed, &target); CHKERRQ(ierr); 2031e284482Svaleriabarra 2041e284482Svaleriabarra // Gather RHS 205b68a8d79SJed Brown CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL); 206b68a8d79SJed Brown ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr); 2071e284482Svaleriabarra ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 208b68a8d79SJed Brown ierr = DMLocalToGlobal(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 2091e284482Svaleriabarra CeedVectorDestroy(&rhsceed); 2101e284482Svaleriabarra 2116a6c615bSJeremy L Thompson // Create the error QFunction 2121e284482Svaleriabarra CeedQFunctionCreateInterior(ceed, 1, bpOptions[rp->bpchoice].error, 213*e83e87a5Sjeremylt bpOptions[rp->bpchoice].errorfname, &qfError); 214*e83e87a5Sjeremylt CeedQFunctionAddInput(qfError, "u", rp->ncompu, CEED_EVAL_INTERP); 215*e83e87a5Sjeremylt CeedQFunctionAddInput(qfError, "true_soln", rp->ncompu, CEED_EVAL_NONE); 216*e83e87a5Sjeremylt CeedQFunctionAddOutput(qfError, "error", rp->ncompu, CEED_EVAL_NONE); 2171e284482Svaleriabarra 2181e284482Svaleriabarra // Create the error operator 219*e83e87a5Sjeremylt CeedOperatorCreate(ceed, qfError, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 220*e83e87a5Sjeremylt &opError); 221*e83e87a5Sjeremylt CeedOperatorSetField(opError, "u", ceeddata->Erestrictu, 2221e284482Svaleriabarra ceeddata->basisu, CEED_VECTOR_ACTIVE); 223*e83e87a5Sjeremylt CeedOperatorSetField(opError, "true_soln", ceeddata->Erestrictui, 2241e284482Svaleriabarra CEED_BASIS_COLLOCATED, target); 225*e83e87a5Sjeremylt CeedOperatorSetField(opError, "error", ceeddata->Erestrictui, 2261e284482Svaleriabarra CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 2271e284482Svaleriabarra 2281e284482Svaleriabarra // Set up Mat 2291e284482Svaleriabarra userO->comm = rp->comm; 2301e284482Svaleriabarra userO->dm = dm; 2311e284482Svaleriabarra userO->Xloc = Xloc; 2321e284482Svaleriabarra ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr); 233*e83e87a5Sjeremylt userO->Xceed = ceeddata->Xceed; 234*e83e87a5Sjeremylt userO->Yceed = ceeddata->Yceed; 235*e83e87a5Sjeremylt userO->op = ceeddata->opApply; 2361e284482Svaleriabarra userO->ceed = ceed; 2371e284482Svaleriabarra 2381e284482Svaleriabarra ierr = KSPCreate(rp->comm, &ksp); CHKERRQ(ierr); 2391e284482Svaleriabarra { 2401e284482Svaleriabarra PC pc; 2411e284482Svaleriabarra ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 2421e284482Svaleriabarra if (rp->bpchoice == CEED_BP1 || rp->bpchoice == CEED_BP2) { 2431e284482Svaleriabarra ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 2441e284482Svaleriabarra ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 2451e284482Svaleriabarra } else { 2461e284482Svaleriabarra ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 2471e284482Svaleriabarra } 2481e284482Svaleriabarra ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 2491e284482Svaleriabarra ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 2501e284482Svaleriabarra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 2511e284482Svaleriabarra PETSC_DEFAULT); CHKERRQ(ierr); 2521e284482Svaleriabarra } 2531e284482Svaleriabarra ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr); 2541e284482Svaleriabarra 255da9108adSvaleriabarra // First run's performance log is not considered for benchmarking purposes 2561e284482Svaleriabarra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 2571e284482Svaleriabarra CHKERRQ(ierr); 2581e284482Svaleriabarra my_rt_start = MPI_Wtime(); 2591e284482Svaleriabarra ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 2601e284482Svaleriabarra my_rt = MPI_Wtime() - my_rt_start; 2611e284482Svaleriabarra ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm); 2621e284482Svaleriabarra CHKERRQ(ierr); 2631e284482Svaleriabarra // Set maxits based on first iteration timing 2641e284482Svaleriabarra if (my_rt > 0.02) { 265565a3730SJed Brown ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 266565a3730SJed Brown rp->ksp_max_it_clip[0]); 2671e284482Svaleriabarra CHKERRQ(ierr); 2681e284482Svaleriabarra } else { 269565a3730SJed Brown ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 270565a3730SJed Brown rp->ksp_max_it_clip[1]); 2711e284482Svaleriabarra CHKERRQ(ierr); 2721e284482Svaleriabarra } 2731e284482Svaleriabarra ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 2741e284482Svaleriabarra 2751e284482Svaleriabarra // Timed solve 2761e284482Svaleriabarra ierr = VecZeroEntries(X); CHKERRQ(ierr); 2771e284482Svaleriabarra ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 2781e284482Svaleriabarra 2791e284482Svaleriabarra // -- Performance logging 2801e284482Svaleriabarra ierr = PetscLogStagePush(rp->solvestage); CHKERRQ(ierr); 2811e284482Svaleriabarra 2821e284482Svaleriabarra // -- Solve 2831e284482Svaleriabarra my_rt_start = MPI_Wtime(); 2841e284482Svaleriabarra ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 2851e284482Svaleriabarra my_rt = MPI_Wtime() - my_rt_start; 2861e284482Svaleriabarra 2871e284482Svaleriabarra // -- Performance logging 2881e284482Svaleriabarra ierr = PetscLogStagePop(); 2891e284482Svaleriabarra 2901e284482Svaleriabarra // Output results 2911e284482Svaleriabarra { 2921e284482Svaleriabarra KSPType ksptype; 2931e284482Svaleriabarra KSPConvergedReason reason; 2941e284482Svaleriabarra PetscReal rnorm; 2951e284482Svaleriabarra PetscInt its; 2961e284482Svaleriabarra ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 2971e284482Svaleriabarra ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 2981e284482Svaleriabarra ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 2991e284482Svaleriabarra ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 3001e284482Svaleriabarra if (!rp->test_mode || reason < 0 || rnorm > 1e-8) { 3011e284482Svaleriabarra ierr = PetscPrintf(rp->comm, 3021e284482Svaleriabarra " KSP:\n" 3031e284482Svaleriabarra " KSP Type : %s\n" 3041e284482Svaleriabarra " KSP Convergence : %s\n" 3051e284482Svaleriabarra " Total KSP Iterations : %D\n" 3061e284482Svaleriabarra " Final rnorm : %e\n", 3071e284482Svaleriabarra ksptype, KSPConvergedReasons[reason], its, 3081e284482Svaleriabarra (double)rnorm); CHKERRQ(ierr); 3091e284482Svaleriabarra } 3101e284482Svaleriabarra if (!rp->test_mode) { 3111e284482Svaleriabarra ierr = PetscPrintf(rp->comm," Performance:\n"); CHKERRQ(ierr); 3121e284482Svaleriabarra } 3131e284482Svaleriabarra { 3141e284482Svaleriabarra PetscReal maxerror; 315*e83e87a5Sjeremylt ierr = ComputeErrorMax(userO, opError, X, target, &maxerror); 3161e284482Svaleriabarra CHKERRQ(ierr); 3171e284482Svaleriabarra PetscReal tol = 5e-2; 3181e284482Svaleriabarra if (!rp->test_mode || maxerror > tol) { 3191e284482Svaleriabarra ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm); 3201e284482Svaleriabarra CHKERRQ(ierr); 3211e284482Svaleriabarra ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm); 3221e284482Svaleriabarra CHKERRQ(ierr); 3231e284482Svaleriabarra ierr = PetscPrintf(rp->comm, 3241e284482Svaleriabarra " Pointwise Error (max) : %e\n" 3251e284482Svaleriabarra " CG Solve Time : %g (%g) sec\n", 3261e284482Svaleriabarra (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 3271e284482Svaleriabarra } 3281e284482Svaleriabarra } 3294c583f1fSvaleriabarra if (!rp->test_mode) { 3301e284482Svaleriabarra ierr = PetscPrintf(rp->comm, 3311e284482Svaleriabarra " DoFs/Sec in CG : %g (%g) million\n", 3321e284482Svaleriabarra 1e-6*gsize*its/rt_max, 3331e284482Svaleriabarra 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 3341e284482Svaleriabarra } 3351e284482Svaleriabarra } 3361e284482Svaleriabarra 3371e284482Svaleriabarra if (rp->write_solution) { 3381e284482Svaleriabarra PetscViewer vtkviewersoln; 3391e284482Svaleriabarra 3401e284482Svaleriabarra ierr = PetscViewerCreate(rp->comm, &vtkviewersoln); CHKERRQ(ierr); 3411e284482Svaleriabarra ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 342a854c6c0Sjeremylt ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtu"); CHKERRQ(ierr); 3431e284482Svaleriabarra ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 3441e284482Svaleriabarra ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 3451e284482Svaleriabarra } 3461e284482Svaleriabarra 3471e284482Svaleriabarra // Cleanup 3481e284482Svaleriabarra ierr = VecDestroy(&X); CHKERRQ(ierr); 3491e284482Svaleriabarra ierr = VecDestroy(&Xloc); CHKERRQ(ierr); 3501e284482Svaleriabarra ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr); 3511e284482Svaleriabarra ierr = MatDestroy(&matO); CHKERRQ(ierr); 3521e284482Svaleriabarra ierr = PetscFree(userO); CHKERRQ(ierr); 3531e284482Svaleriabarra ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr); 3541e284482Svaleriabarra 3551e284482Svaleriabarra ierr = VecDestroy(&rhs); CHKERRQ(ierr); 3561e284482Svaleriabarra ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 3571e284482Svaleriabarra ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 3581e284482Svaleriabarra CeedVectorDestroy(&target); 359*e83e87a5Sjeremylt CeedQFunctionDestroy(&qfError); 360*e83e87a5Sjeremylt CeedOperatorDestroy(&opError); 3611e284482Svaleriabarra CeedDestroy(&ceed); 3621e284482Svaleriabarra PetscFunctionReturn(0); 3631e284482Svaleriabarra } 3641e284482Svaleriabarra 365*e83e87a5Sjeremylt static PetscErrorCode Run(RunParams rp, PetscInt num_resources, 366*e83e87a5Sjeremylt char *const *ceedresources, PetscInt num_bpchoices, 367*e83e87a5Sjeremylt const bpType *bpchoices) { 3685f284d84SJed Brown PetscInt ierr; 3695f284d84SJed Brown DM dm; 3705f284d84SJed Brown 3715f284d84SJed Brown PetscFunctionBeginUser; 3725f284d84SJed Brown // Setup DM 3735f284d84SJed Brown if (rp->read_mesh) { 3745f284d84SJed Brown ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, PETSC_TRUE, &dm); 3755f284d84SJed Brown CHKERRQ(ierr); 3765f284d84SJed Brown } else { 3775f284d84SJed Brown if (rp->userlnodes) { 3785f284d84SJed Brown // Find a nicely composite number of elements no less than global nodes 3795f284d84SJed Brown PetscMPIInt size; 3805f284d84SJed Brown ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr); 3815f284d84SJed Brown for (PetscInt gelem = 3825f284d84SJed Brown PetscMax(1, size * rp->localnodes / PetscPowInt(rp->degree, rp->dim)); 3835f284d84SJed Brown ; 3845f284d84SJed Brown gelem++) { 3855f284d84SJed Brown Split3(gelem, rp->melem, true); 3865f284d84SJed Brown if (Max3(rp->melem) / Min3(rp->melem) <= 2) break; 3875f284d84SJed Brown } 3885f284d84SJed Brown } 3895f284d84SJed Brown ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, PETSC_FALSE, rp->melem, 3905f284d84SJed Brown NULL, NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr); 3915f284d84SJed Brown } 3925f284d84SJed Brown 3935f284d84SJed Brown { 3945f284d84SJed Brown DM dmDist = NULL; 3955f284d84SJed Brown PetscPartitioner part; 3965f284d84SJed Brown 3975f284d84SJed Brown ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 3985f284d84SJed Brown ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 3995f284d84SJed Brown ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 4005f284d84SJed Brown if (dmDist) { 4015f284d84SJed Brown ierr = DMDestroy(&dm); CHKERRQ(ierr); 4025f284d84SJed Brown dm = dmDist; 4035f284d84SJed Brown } 4045f284d84SJed Brown } 405b68a8d79SJed Brown // Disable default VECSTANDARD *after* distribution (which creates a Vec) 406b68a8d79SJed Brown ierr = DMSetVecType(dm, NULL); CHKERRQ(ierr); 4075f284d84SJed Brown 4085f284d84SJed Brown for (PetscInt b = 0; b < num_bpchoices; b++) { 4095f284d84SJed Brown DM dm_deg; 410b68a8d79SJed Brown VecType vectype; 4115f284d84SJed Brown PetscInt qextra = rp->qextra; 4125f284d84SJed Brown rp->bpchoice = bpchoices[b]; 4135f284d84SJed Brown rp->ncompu = bpOptions[rp->bpchoice].ncompu; 4145f284d84SJed Brown rp->qextra = qextra < 0 ? bpOptions[rp->bpchoice].qextra : qextra; 4155f284d84SJed Brown ierr = DMClone(dm, &dm_deg); CHKERRQ(ierr); 416b68a8d79SJed Brown ierr = DMGetVecType(dm, &vectype); CHKERRQ(ierr); 417b68a8d79SJed Brown ierr = DMSetVecType(dm_deg, vectype); CHKERRQ(ierr); 4185f284d84SJed Brown // Create DM 419*e83e87a5Sjeremylt PetscInt dim; 420*e83e87a5Sjeremylt ierr = DMGetDimension(dm_deg, &dim); CHKERRQ(ierr); 421*e83e87a5Sjeremylt ierr = SetupDMByDegree(dm_deg, rp->degree, rp->ncompu, dim, 422*e83e87a5Sjeremylt bpOptions[rp->bpchoice].enforcebc, 423*e83e87a5Sjeremylt bpOptions[rp->bpchoice].bcsfunc); CHKERRQ(ierr); 4245f284d84SJed Brown for (PetscInt r = 0; r < num_resources; r++) { 4255f284d84SJed Brown ierr = RunWithDM(rp, dm_deg, ceedresources[r]); CHKERRQ(ierr); 4265f284d84SJed Brown } 4275f284d84SJed Brown ierr = DMDestroy(&dm_deg); CHKERRQ(ierr); 4285f284d84SJed Brown rp->qextra = qextra; 4295f284d84SJed Brown } 4305f284d84SJed Brown 4315f284d84SJed Brown ierr = DMDestroy(&dm); CHKERRQ(ierr); 4325f284d84SJed Brown PetscFunctionReturn(0); 4335f284d84SJed Brown } 4345f284d84SJed Brown 4351e284482Svaleriabarra int main(int argc, char **argv) { 4361e284482Svaleriabarra PetscInt ierr, commsize; 4371e284482Svaleriabarra RunParams rp; 43853a0f73bSJed Brown MPI_Comm comm; 43953a0f73bSJed Brown char filename[PETSC_MAX_PATH_LEN]; 440565a3730SJed Brown char *ceedresources[30]; 441565a3730SJed Brown PetscInt num_ceedresources = 30; 44253a0f73bSJed Brown char hostname[PETSC_MAX_PATH_LEN]; 4431e284482Svaleriabarra 444c36f77d8SJed Brown PetscInt dim = 3, melem[3] = {3, 3, 3}; 44553a0f73bSJed Brown PetscInt num_degrees = 30, degree[30] = {}, num_localnodes = 2, localnodes[2] = {}; 44653a0f73bSJed Brown PetscMPIInt rankspernode; 447c36f77d8SJed Brown PetscBool degree_set; 448565a3730SJed Brown bpType bpchoices[10]; 449565a3730SJed Brown PetscInt num_bpchoices = 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; 4551e284482Svaleriabarra ierr = MPI_Comm_size(comm, &commsize); 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); 46353a0f73bSJed Brown ierr = MPI_Comm_size(splitcomm, &rankspernode); CHKERRQ(ierr); 46453a0f73bSJed Brown ierr = MPI_Comm_free(&splitcomm); CHKERRQ(ierr); 46553a0f73bSJed Brown } 46653a0f73bSJed Brown #else 46753a0f73bSJed Brown rankspernode = -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, 481565a3730SJed Brown bpTypes, (PetscEnum *)bpchoices, &num_bpchoices, &set); 482565a3730SJed Brown CHKERRQ(ierr); 483565a3730SJed Brown if (!set) { 484565a3730SJed Brown bpchoices[0] = CEED_BP1; 485565a3730SJed Brown num_bpchoices = 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, °ree_set); CHKERRQ(ierr); 50053a0f73bSJed Brown if (!degree_set) 50153a0f73bSJed Brown num_degrees = 1; 5025f284d84SJed Brown rp->qextra = PETSC_DECIDE; 5035f284d84SJed Brown ierr = PetscOptionsInt("-qextra", 5045f284d84SJed Brown "Number of extra quadrature points (-1 for auto)", NULL, 505c36f77d8SJed Brown rp->qextra, &rp->qextra, NULL); CHKERRQ(ierr); 506565a3730SJed Brown { 507565a3730SJed Brown PetscBool set; 508565a3730SJed Brown ierr = PetscOptionsStringArray("-ceed", 509565a3730SJed Brown "CEED resource specifier (comma-separated list)", NULL, 510565a3730SJed Brown ceedresources, &num_ceedresources, &set); CHKERRQ(ierr); 511565a3730SJed Brown if (!set) { 512565a3730SJed Brown ierr = PetscStrallocpy( "/cpu/self", &ceedresources[0]); CHKERRQ(ierr); 513565a3730SJed Brown num_ceedresources = 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, 527cb32e2e7SValeria Barra melem, &tmp, NULL); CHKERRQ(ierr); 528cb32e2e7SValeria Barra } 52953a0f73bSJed Brown localnodes[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)", 533c36f77d8SJed Brown NULL, localnodes, &num_localnodes, &rp->userlnodes); 534153bcb04Svaleriabarra CHKERRQ(ierr); 53553a0f73bSJed Brown if (num_localnodes < 2) 53653a0f73bSJed Brown localnodes[1] = 2 * localnodes[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) { 54653a0f73bSJed Brown PetscInt maxdegree = 8; 5475f284d84SJed Brown ierr = PetscOptionsInt("-max_degree", 5485f284d84SJed Brown "Range of degrees [1, maxdegree] to run with", 54953a0f73bSJed Brown NULL, maxdegree, &maxdegree, NULL); 5501e284482Svaleriabarra CHKERRQ(ierr); 55153a0f73bSJed Brown for (PetscInt i = 0; i < maxdegree; i++) 55253a0f73bSJed Brown degree[i] = i + 1; 55353a0f73bSJed Brown num_degrees = maxdegree; 55453a0f73bSJed Brown } 55553a0f73bSJed Brown { 55653a0f73bSJed Brown PetscBool flg; 55753a0f73bSJed Brown PetscInt p = rankspernode; 55853a0f73bSJed Brown ierr = PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, 55953a0f73bSJed Brown p, &p, &flg); 56053a0f73bSJed Brown CHKERRQ(ierr); 56153a0f73bSJed Brown if (flg) rankspernode = p; 56253a0f73bSJed Brown } 5639396343dSjeremylt 56453a0f73bSJed Brown ierr = PetscOptionsEnd(); 56553a0f73bSJed Brown CHKERRQ(ierr); 566cb32e2e7SValeria Barra 5671e284482Svaleriabarra // Register PETSc logging stage 568c36f77d8SJed Brown ierr = PetscLogStageRegister("Solve Stage", &rp->solvestage); 56953a0f73bSJed Brown CHKERRQ(ierr); 57009a940d7Sjeremylt 57153a0f73bSJed Brown rp->hostname = hostname; 5721e284482Svaleriabarra rp->dim = dim; 5731e284482Svaleriabarra rp->melem = melem; 57453a0f73bSJed Brown rp->rankspernode = rankspernode; 575cb32e2e7SValeria Barra 57653a0f73bSJed Brown for (PetscInt d = 0; d < num_degrees; d++) { 57753a0f73bSJed Brown PetscInt deg = degree[d]; 57853a0f73bSJed Brown for (PetscInt n = localnodes[0]; n < localnodes[1]; n *= 2) { 57953a0f73bSJed Brown rp->degree = deg; 58053a0f73bSJed Brown rp->localnodes = n; 5815f284d84SJed Brown ierr = Run(rp, num_ceedresources, ceedresources, 5825f284d84SJed Brown num_bpchoices, bpchoices); CHKERRQ(ierr); 583565a3730SJed Brown } 584565a3730SJed Brown } 5851e284482Svaleriabarra // Clear memory 5861e284482Svaleriabarra ierr = PetscFree(rp); CHKERRQ(ierr); 587565a3730SJed Brown for (PetscInt i=0; i<num_ceedresources; i++) { 588565a3730SJed Brown ierr = PetscFree(ceedresources[i]); CHKERRQ(ierr); 589565a3730SJed Brown } 5900c59ef15SJeremy L Thompson return PetscFinalize(); 5910c59ef15SJeremy L Thompson } 592