16c5df90dSjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 26c5df90dSjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 36c5df90dSjeremylt // reserved. See files LICENSE and NOTICE for details. 46c5df90dSjeremylt // 56c5df90dSjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software 66c5df90dSjeremylt // libraries and APIs for efficient high-order finite element and spectral 76c5df90dSjeremylt // element discretizations for exascale applications. For more information and 86c5df90dSjeremylt // source code availability see http://github.com/ceed. 96c5df90dSjeremylt // 106c5df90dSjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 116c5df90dSjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office 126c5df90dSjeremylt // of Science and the National Nuclear Security Administration) responsible for 136c5df90dSjeremylt // the planning and preparation of a capable exascale ecosystem, including 146c5df90dSjeremylt // software, applications, hardware, advanced system engineering and early 156c5df90dSjeremylt // testbed platforms, in support of the nation's exascale computing imperative. 166c5df90dSjeremylt 176c5df90dSjeremylt // libCEED + PETSc Example: CEED BPs 3-6 with Multigrid 186c5df90dSjeremylt // 196c5df90dSjeremylt // This example demonstrates a simple usage of libCEED with PETSc to solve the 206c5df90dSjeremylt // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps. 216c5df90dSjeremylt // 226c5df90dSjeremylt // The code uses higher level communication protocols in DMPlex. 236c5df90dSjeremylt // 246c5df90dSjeremylt // Build with: 256c5df90dSjeremylt // 266c5df90dSjeremylt // make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 276c5df90dSjeremylt // 286c5df90dSjeremylt // Sample runs: 296c5df90dSjeremylt // 306c5df90dSjeremylt // multigrid -problem bp3 316c5df90dSjeremylt // multigrid -problem bp4 -ceed /cpu/self 326c5df90dSjeremylt // multigrid -problem bp5 -ceed /cpu/occa 336c5df90dSjeremylt // multigrid -problem bp6 -ceed /gpu/cuda 346c5df90dSjeremylt // 356c5df90dSjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 366c5df90dSjeremylt 376c5df90dSjeremylt /// @file 386c5df90dSjeremylt /// CEED BPs 1-6 multigrid example using PETSc 396c5df90dSjeremylt const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n"; 406c5df90dSjeremylt 416c5df90dSjeremylt #define multigrid 426c5df90dSjeremylt #include "setup.h" 436c5df90dSjeremylt 446c5df90dSjeremylt int main(int argc, char **argv) { 456c5df90dSjeremylt PetscInt ierr; 466c5df90dSjeremylt MPI_Comm comm; 47cb0b5415Sjeremylt char filename[PETSC_MAX_PATH_LEN], 48cb0b5415Sjeremylt ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 496c5df90dSjeremylt double my_rt_start, my_rt, rt_min, rt_max; 5061608365Sjeremylt PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3, fineLevel, 516c5df90dSjeremylt melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees; 526c5df90dSjeremylt PetscScalar *r; 536c5df90dSjeremylt PetscBool test_mode, benchmark_mode, read_mesh, write_solution; 54199551f5Sjeremylt DM *dm, dmorig; 55c70bd2a0Sjeremylt SNES snesdummy; 566c5df90dSjeremylt KSP ksp; 576c5df90dSjeremylt PC pc; 58199551f5Sjeremylt Mat *matO, *matPR, matcoarse; 59cdb3667fSjeremylt Vec *X, *Xloc, *mult, rhs, rhsloc; 606c5df90dSjeremylt UserO *userO; 61a97643b0Sjeremylt UserProlongRestr *userPR; 626c5df90dSjeremylt Ceed ceed; 636c5df90dSjeremylt CeedData *ceeddata; 64*9396343dSjeremylt CeedMemType memtyperequested; 655e81177dSjeremylt CeedVector rhsceed, target; 66c70bd2a0Sjeremylt CeedQFunction qferror, qfrestrict, qfprolong; 67c70bd2a0Sjeremylt CeedOperator operror; 68199551f5Sjeremylt bpType bpchoice; 696c5df90dSjeremylt coarsenType coarsen; 706c5df90dSjeremylt 716c5df90dSjeremylt ierr = PetscInitialize(&argc, &argv, NULL, help); 726c5df90dSjeremylt if (ierr) return ierr; 736c5df90dSjeremylt comm = PETSC_COMM_WORLD; 746c5df90dSjeremylt 75*9396343dSjeremylt // Check PETSc CUDA avaliability 76*9396343dSjeremylt PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE; 77*9396343dSjeremylt // *INDENT-OFF* 78*9396343dSjeremylt #ifdef PETSC_HAVE_CUDA 79*9396343dSjeremylt petschavecuda = PETSC_TRUE; 80*9396343dSjeremylt #else 81*9396343dSjeremylt petschavecuda = PETSC_FALSE; 82*9396343dSjeremylt #endif 83*9396343dSjeremylt // *INDENT-ON* 84*9396343dSjeremylt 856c5df90dSjeremylt // Parse command line options 866c5df90dSjeremylt ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 87199551f5Sjeremylt bpchoice = CEED_BP3; 886c5df90dSjeremylt ierr = PetscOptionsEnum("-problem", 896c5df90dSjeremylt "CEED benchmark problem to solve", NULL, 90199551f5Sjeremylt bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice, 916c5df90dSjeremylt NULL); CHKERRQ(ierr); 92199551f5Sjeremylt ncompu = bpOptions[bpchoice].ncompu; 936c5df90dSjeremylt test_mode = PETSC_FALSE; 946c5df90dSjeremylt ierr = PetscOptionsBool("-test", 956c5df90dSjeremylt "Testing mode (do not print unless error is large)", 966c5df90dSjeremylt NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 976c5df90dSjeremylt benchmark_mode = PETSC_FALSE; 986c5df90dSjeremylt ierr = PetscOptionsBool("-benchmark", 996c5df90dSjeremylt "Benchmarking mode (prints benchmark statistics)", 1006c5df90dSjeremylt NULL, benchmark_mode, &benchmark_mode, NULL); 1016c5df90dSjeremylt CHKERRQ(ierr); 1026c5df90dSjeremylt write_solution = PETSC_FALSE; 1036c5df90dSjeremylt ierr = PetscOptionsBool("-write_solution", 1046c5df90dSjeremylt "Write solution for visualization", 1056c5df90dSjeremylt NULL, write_solution, &write_solution, NULL); 1066c5df90dSjeremylt CHKERRQ(ierr); 1076c5df90dSjeremylt degree = test_mode ? 3 : 2; 1086c5df90dSjeremylt ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 1096c5df90dSjeremylt NULL, degree, °ree, NULL); CHKERRQ(ierr); 1106c5df90dSjeremylt if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 1116c5df90dSjeremylt "-degree %D must be at least 1", degree); 112199551f5Sjeremylt qextra = bpOptions[bpchoice].qextra; 1136c5df90dSjeremylt ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 1146c5df90dSjeremylt NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 1156c5df90dSjeremylt ierr = PetscOptionsString("-ceed", "CEED resource specifier", 1166c5df90dSjeremylt NULL, ceedresource, ceedresource, 1176c5df90dSjeremylt sizeof(ceedresource), NULL); CHKERRQ(ierr); 1186c5df90dSjeremylt coarsen = COARSEN_UNIFORM; 1196c5df90dSjeremylt ierr = PetscOptionsEnum("-coarsen", 1206c5df90dSjeremylt "Coarsening strategy to use", NULL, 1216c5df90dSjeremylt coarsenTypes, (PetscEnum)coarsen, 1226c5df90dSjeremylt (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr); 123cb32e2e7SValeria Barra read_mesh = PETSC_FALSE; 1246c5df90dSjeremylt ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 1256c5df90dSjeremylt filename, filename, sizeof(filename), &read_mesh); 1266c5df90dSjeremylt CHKERRQ(ierr); 1276c5df90dSjeremylt if (!read_mesh) { 1286c5df90dSjeremylt PetscInt tmp = dim; 1296c5df90dSjeremylt ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL, 1306c5df90dSjeremylt melem, &tmp, NULL); CHKERRQ(ierr); 1316c5df90dSjeremylt } 132*9396343dSjeremylt memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 133*9396343dSjeremylt ierr = PetscOptionsEnum("-memtype", 134*9396343dSjeremylt "CEED MemType requested", NULL, 135*9396343dSjeremylt memTypes, (PetscEnum)memtyperequested, 136*9396343dSjeremylt (PetscEnum *)&memtyperequested, &setmemtyperequest); 137*9396343dSjeremylt CHKERRQ(ierr); 1386c5df90dSjeremylt ierr = PetscOptionsEnd(); CHKERRQ(ierr); 1396c5df90dSjeremylt 140*9396343dSjeremylt // Set up libCEED 141*9396343dSjeremylt CeedInit(ceedresource, &ceed); 142*9396343dSjeremylt CeedMemType memtypebackend; 143*9396343dSjeremylt CeedGetPreferredMemType(ceed, &memtypebackend); 144*9396343dSjeremylt 145*9396343dSjeremylt // Check memtype compatibility 146*9396343dSjeremylt if (!setmemtyperequest) 147*9396343dSjeremylt memtyperequested = memtypebackend; 148*9396343dSjeremylt else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE) 149*9396343dSjeremylt SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 150*9396343dSjeremylt "PETSc was not built with CUDA. " 151*9396343dSjeremylt "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 152*9396343dSjeremylt 1536c5df90dSjeremylt // Setup DM 1546c5df90dSjeremylt if (read_mesh) { 155199551f5Sjeremylt ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig); 1566c5df90dSjeremylt CHKERRQ(ierr); 1576c5df90dSjeremylt } else { 1586c5df90dSjeremylt ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL, 159199551f5Sjeremylt NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr); 1606c5df90dSjeremylt } 1616c5df90dSjeremylt 1626c5df90dSjeremylt { 1636c5df90dSjeremylt DM dmDist = NULL; 1646c5df90dSjeremylt PetscPartitioner part; 1656c5df90dSjeremylt 166199551f5Sjeremylt ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr); 1676c5df90dSjeremylt ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 168199551f5Sjeremylt ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr); 1696c5df90dSjeremylt if (dmDist) { 170199551f5Sjeremylt ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 171199551f5Sjeremylt dmorig = dmDist; 1726c5df90dSjeremylt } 1736c5df90dSjeremylt } 1746c5df90dSjeremylt 1756c5df90dSjeremylt // Allocate arrays for PETSc objects for each level 1766c5df90dSjeremylt switch (coarsen) { 1776c5df90dSjeremylt case COARSEN_UNIFORM: 1786c5df90dSjeremylt numlevels = degree; 1796c5df90dSjeremylt break; 180dc7d240cSValeria Barra case COARSEN_LOGARITHMIC: 1816c5df90dSjeremylt numlevels = ceil(log(degree)/log(2)) + 1; 1826c5df90dSjeremylt break; 1836c5df90dSjeremylt } 1846c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr); 18561608365Sjeremylt fineLevel = numlevels - 1; 18661608365Sjeremylt 1876c5df90dSjeremylt switch (coarsen) { 1886c5df90dSjeremylt case COARSEN_UNIFORM: 1896c5df90dSjeremylt for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1; 1906c5df90dSjeremylt break; 191dc7d240cSValeria Barra case COARSEN_LOGARITHMIC: 1926c5df90dSjeremylt for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i); 19361608365Sjeremylt leveldegrees[fineLevel] = degree; 1946c5df90dSjeremylt break; 1956c5df90dSjeremylt } 1966c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr); 1976c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr); 1986c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr); 1996c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr); 2006c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr); 201a97643b0Sjeremylt ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr); 2026c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr); 203a97643b0Sjeremylt ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr); 2046c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr); 2056c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr); 2066c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr); 2076c5df90dSjeremylt 2086c5df90dSjeremylt // Setup DM and Operator Mat Shells for each level 2096c5df90dSjeremylt for (CeedInt i=0; i<numlevels; i++) { 2106c5df90dSjeremylt // Create DM 211199551f5Sjeremylt ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr); 212199551f5Sjeremylt ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpchoice); 2136c5df90dSjeremylt CHKERRQ(ierr); 2146c5df90dSjeremylt 2156c5df90dSjeremylt // Create vectors 216*9396343dSjeremylt if (memtyperequested == CEED_MEM_DEVICE) { 217*9396343dSjeremylt ierr = DMSetVecType(dm[i], VECCUDA); CHKERRQ(ierr); 218*9396343dSjeremylt } 2196c5df90dSjeremylt ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr); 2206c5df90dSjeremylt ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr); 2216c5df90dSjeremylt ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr); 2226c5df90dSjeremylt ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr); 2236c5df90dSjeremylt ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr); 2246c5df90dSjeremylt 2256c5df90dSjeremylt // Operator 2266c5df90dSjeremylt ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr); 2276c5df90dSjeremylt ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i], 2286c5df90dSjeremylt userO[i], &matO[i]); CHKERRQ(ierr); 2296c5df90dSjeremylt ierr = MatShellSetOperation(matO[i], MATOP_MULT, 230ce74dcefSjeremylt (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 2316c5df90dSjeremylt ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL, 232ce74dcefSjeremylt (void(*)(void))MatGetDiag); CHKERRQ(ierr); 233*9396343dSjeremylt if (memtyperequested == CEED_MEM_DEVICE) { 234*9396343dSjeremylt ierr = MatShellSetVecType(matO[i], VECCUDA); CHKERRQ(ierr); 235*9396343dSjeremylt } 2366c5df90dSjeremylt 2376c5df90dSjeremylt // Level transfers 2386c5df90dSjeremylt if (i > 0) { 2396c5df90dSjeremylt // Interp 240a97643b0Sjeremylt ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr); 2416c5df90dSjeremylt ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1], 242a97643b0Sjeremylt userPR[i], &matPR[i]); CHKERRQ(ierr); 243a97643b0Sjeremylt ierr = MatShellSetOperation(matPR[i], MATOP_MULT, 244a97643b0Sjeremylt (void(*)(void))MatMult_Prolong); 2456c5df90dSjeremylt CHKERRQ(ierr); 246a97643b0Sjeremylt ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE, 2476c5df90dSjeremylt (void(*)(void))MatMult_Restrict); 2486c5df90dSjeremylt CHKERRQ(ierr); 249*9396343dSjeremylt if (memtyperequested == CEED_MEM_DEVICE) { 250*9396343dSjeremylt ierr = MatShellSetVecType(matPR[i], VECCUDA); CHKERRQ(ierr); 251*9396343dSjeremylt } 2526c5df90dSjeremylt } 2536c5df90dSjeremylt } 25461608365Sjeremylt ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr); 2556c5df90dSjeremylt 2562d03409cSjeremylt // Set up libCEED 2572d03409cSjeremylt CeedInit(ceedresource, &ceed); 2582d03409cSjeremylt 2596c5df90dSjeremylt // Print global grid information 2606c5df90dSjeremylt if (!test_mode) { 2616c5df90dSjeremylt PetscInt P = degree + 1, Q = P + qextra; 262*9396343dSjeremylt 2632d03409cSjeremylt const char *usedresource; 2642d03409cSjeremylt CeedGetResource(ceed, &usedresource); 265*9396343dSjeremylt 266*9396343dSjeremylt VecType vectype; 267*9396343dSjeremylt ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr); 268*9396343dSjeremylt 2696c5df90dSjeremylt ierr = PetscPrintf(comm, 2706c5df90dSjeremylt "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n" 271*9396343dSjeremylt " PETSc:\n" 272*9396343dSjeremylt " PETSc Vec Type : %s\n" 2736c5df90dSjeremylt " libCEED:\n" 2746c5df90dSjeremylt " libCEED Backend : %s\n" 275*9396343dSjeremylt " libCEED Backend MemType : %s\n" 276*9396343dSjeremylt " libCEED User Requested MemType : %s\n" 2776c5df90dSjeremylt " Mesh:\n" 2786c5df90dSjeremylt " Number of 1D Basis Nodes (p) : %d\n" 2796c5df90dSjeremylt " Number of 1D Quadrature Points (q) : %d\n" 2806c5df90dSjeremylt " Global Nodes : %D\n" 2816c5df90dSjeremylt " Owned Nodes : %D\n" 282db419314Sjeremylt " DoF per node : %D\n" 2836c5df90dSjeremylt " Multigrid:\n" 2846c5df90dSjeremylt " Number of Levels : %d\n", 285*9396343dSjeremylt bpchoice+1, vectype, usedresource, 286*9396343dSjeremylt CeedMemTypes[memtypebackend], 287*9396343dSjeremylt (setmemtyperequest) ? 288*9396343dSjeremylt CeedMemTypes[memtyperequested] : "none", 289*9396343dSjeremylt P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu, 290db419314Sjeremylt ncompu, numlevels); CHKERRQ(ierr); 2916c5df90dSjeremylt } 2926c5df90dSjeremylt 2936c5df90dSjeremylt // Create RHS vector 29461608365Sjeremylt ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr); 2956c5df90dSjeremylt ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 296*9396343dSjeremylt if (memtyperequested == CEED_MEM_HOST) { 2976c5df90dSjeremylt ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 298*9396343dSjeremylt } else { 299*9396343dSjeremylt ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr); 300*9396343dSjeremylt } 30161608365Sjeremylt CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed); 3026c5df90dSjeremylt CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r); 3036c5df90dSjeremylt 3046c5df90dSjeremylt // Set up libCEED operators on each level 3056c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr); 3066c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 3076c5df90dSjeremylt // Print level information 30861608365Sjeremylt if (!test_mode && (i == 0 || i == fineLevel)) { 3096c5df90dSjeremylt ierr = PetscPrintf(comm," Level %D (%s):\n" 3106c5df90dSjeremylt " Number of 1D Basis Nodes (p) : %d\n" 3116c5df90dSjeremylt " Global Nodes : %D\n" 3126c5df90dSjeremylt " Owned Nodes : %D\n", 3136c5df90dSjeremylt i, (i? "fine" : "coarse"), leveldegrees[i] + 1, 3146c5df90dSjeremylt gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr); 3156c5df90dSjeremylt } 3166c5df90dSjeremylt ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr); 3176c5df90dSjeremylt ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra, 318199551f5Sjeremylt ncompu, gsize[i], xlsize[i], bpchoice, 31961608365Sjeremylt ceeddata[i], i==(fineLevel), rhsceed, 3207f823360Sjeremylt &target); CHKERRQ(ierr); 3216c5df90dSjeremylt } 3226c5df90dSjeremylt 3236c5df90dSjeremylt // Gather RHS 324*9396343dSjeremylt CeedVectorSyncArray(rhsceed, memtyperequested); 325*9396343dSjeremylt if (memtyperequested == CEED_MEM_HOST) { 3266c5df90dSjeremylt ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 327*9396343dSjeremylt } else { 328*9396343dSjeremylt ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr); 329*9396343dSjeremylt } 3306c5df90dSjeremylt ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 331*9396343dSjeremylt ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 3326c5df90dSjeremylt CeedVectorDestroy(&rhsceed); 3336c5df90dSjeremylt 3346c5df90dSjeremylt // Create the restriction/interpolation Q-function 33560f77c51Sjeremylt CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 336c70bd2a0Sjeremylt &qfrestrict); 33760f77c51Sjeremylt CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 338c70bd2a0Sjeremylt &qfprolong); 3396c5df90dSjeremylt 3406c5df90dSjeremylt // Set up libCEED level transfer operators 341199551f5Sjeremylt ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpchoice, ceeddata, 342c70bd2a0Sjeremylt leveldegrees, qfrestrict, qfprolong); 3436c5df90dSjeremylt CHKERRQ(ierr); 3446c5df90dSjeremylt 3456c5df90dSjeremylt // Create the error Q-function 346199551f5Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error, 347199551f5Sjeremylt bpOptions[bpchoice].errorfname, &qferror); 348c70bd2a0Sjeremylt CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP); 349c70bd2a0Sjeremylt CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE); 350c70bd2a0Sjeremylt CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE); 3516c5df90dSjeremylt 3526c5df90dSjeremylt // Create the error operator 353c70bd2a0Sjeremylt CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 354c70bd2a0Sjeremylt &operror); 355c70bd2a0Sjeremylt CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu, 35661608365Sjeremylt ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE); 357c70bd2a0Sjeremylt CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui, 358a8d32208Sjeremylt CEED_BASIS_COLLOCATED, target); 359c70bd2a0Sjeremylt CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui, 360a8d32208Sjeremylt CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 3616c5df90dSjeremylt 3626c5df90dSjeremylt // Calculate multiplicity 3636c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 3646c5df90dSjeremylt PetscScalar *x; 3656c5df90dSjeremylt 3666c5df90dSjeremylt // CEED vector 3676c5df90dSjeremylt ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr); 3686c5df90dSjeremylt CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 3696c5df90dSjeremylt 3706c5df90dSjeremylt // Multiplicity 371a8d32208Sjeremylt CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu, 3726c5df90dSjeremylt ceeddata[i]->xceed); 3736c5df90dSjeremylt 3746c5df90dSjeremylt // Restore vector 3756c5df90dSjeremylt ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr); 3766c5df90dSjeremylt 3776c5df90dSjeremylt // Creat mult vector 3786c5df90dSjeremylt ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr); 3796c5df90dSjeremylt 3806c5df90dSjeremylt // Local-to-global 3816c5df90dSjeremylt ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 382483f8b0dSjeremylt ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]); 3836c5df90dSjeremylt CHKERRQ(ierr); 3846c5df90dSjeremylt ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 3856c5df90dSjeremylt 3866c5df90dSjeremylt // Global-to-local 387483f8b0dSjeremylt ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]); 3886c5df90dSjeremylt CHKERRQ(ierr); 3896c5df90dSjeremylt ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 3906c5df90dSjeremylt 3916c5df90dSjeremylt // Multiplicity scaling 3926c5df90dSjeremylt ierr = VecReciprocal(mult[i]); 3936c5df90dSjeremylt } 3946c5df90dSjeremylt 3956c5df90dSjeremylt // Set up Mat 3966c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 397226c3a8fSjeremylt // User Operator 3986c5df90dSjeremylt userO[i]->comm = comm; 3996c5df90dSjeremylt userO[i]->dm = dm[i]; 4006c5df90dSjeremylt userO[i]->Xloc = Xloc[i]; 4016c5df90dSjeremylt ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr); 4026c5df90dSjeremylt userO[i]->xceed = ceeddata[i]->xceed; 4036c5df90dSjeremylt userO[i]->yceed = ceeddata[i]->yceed; 404c70bd2a0Sjeremylt userO[i]->op = ceeddata[i]->opapply; 4056c5df90dSjeremylt userO[i]->ceed = ceed; 406*9396343dSjeremylt userO[i]->memtype = memtyperequested; 407*9396343dSjeremylt if (memtyperequested == CEED_MEM_HOST) { 408*9396343dSjeremylt userO[i]->VecGetArray = VecGetArray; 409*9396343dSjeremylt userO[i]->VecGetArrayRead = VecGetArrayRead; 410*9396343dSjeremylt userO[i]->VecRestoreArray = VecRestoreArray; 411*9396343dSjeremylt userO[i]->VecRestoreArrayRead = VecRestoreArrayRead; 412*9396343dSjeremylt } else { 413*9396343dSjeremylt userO[i]->VecGetArray = VecCUDAGetArray; 414*9396343dSjeremylt userO[i]->VecGetArrayRead = VecCUDAGetArrayRead; 415*9396343dSjeremylt userO[i]->VecRestoreArray = VecCUDARestoreArray; 416*9396343dSjeremylt userO[i]->VecRestoreArrayRead = VecCUDARestoreArrayRead; 417*9396343dSjeremylt } 4186c5df90dSjeremylt 4196c5df90dSjeremylt if (i > 0) { 420a97643b0Sjeremylt // Prolongation/Restriction Operator 421a97643b0Sjeremylt userPR[i]->comm = comm; 422199551f5Sjeremylt userPR[i]->dmf = dm[i]; 423199551f5Sjeremylt userPR[i]->dmc = dm[i-1]; 424199551f5Sjeremylt userPR[i]->locvecc = Xloc[i-1]; 425199551f5Sjeremylt userPR[i]->locvecf = userO[i]->Yloc; 426199551f5Sjeremylt userPR[i]->multvec = mult[i]; 427199551f5Sjeremylt userPR[i]->ceedvecc = userO[i-1]->xceed; 428199551f5Sjeremylt userPR[i]->ceedvecf = userO[i]->yceed; 429c70bd2a0Sjeremylt userPR[i]->opprolong = ceeddata[i]->opprolong; 430c70bd2a0Sjeremylt userPR[i]->oprestrict = ceeddata[i]->oprestrict; 431a97643b0Sjeremylt userPR[i]->ceed = ceed; 432*9396343dSjeremylt userPR[i]->memtype = userO[i]->memtype; 433*9396343dSjeremylt userPR[i]->VecGetArray = userO[i]->VecGetArray; 434*9396343dSjeremylt userPR[i]->VecGetArrayRead = userO[i]->VecGetArrayRead; 435*9396343dSjeremylt userPR[i]->VecRestoreArray = userO[i]->VecRestoreArray; 436*9396343dSjeremylt userPR[i]->VecRestoreArrayRead = userO[i]->VecRestoreArrayRead; 4376c5df90dSjeremylt } 4386c5df90dSjeremylt } 4396c5df90dSjeremylt 44015ce0ef0Sjeremylt // Setup dummy SNES for AMG coarse solve 441c70bd2a0Sjeremylt ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr); 442c70bd2a0Sjeremylt ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr); 443c70bd2a0Sjeremylt ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr); 44415ce0ef0Sjeremylt 44515ce0ef0Sjeremylt // -- Jacobian matrix 44615ce0ef0Sjeremylt ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr); 447199551f5Sjeremylt ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr); 448199551f5Sjeremylt ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL, 44915ce0ef0Sjeremylt NULL); CHKERRQ(ierr); 45015ce0ef0Sjeremylt 45115ce0ef0Sjeremylt // -- Residual evaluation function 452c70bd2a0Sjeremylt ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed, 45315ce0ef0Sjeremylt userO[0]); CHKERRQ(ierr); 45415ce0ef0Sjeremylt 45515ce0ef0Sjeremylt // -- Form Jacobian 456c70bd2a0Sjeremylt ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0], 457199551f5Sjeremylt matcoarse, NULL); CHKERRQ(ierr); 45815ce0ef0Sjeremylt 4596c5df90dSjeremylt // Set up KSP 4606c5df90dSjeremylt ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 4616c5df90dSjeremylt { 4626c5df90dSjeremylt ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 4636c5df90dSjeremylt ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 4646c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 4656c5df90dSjeremylt PETSC_DEFAULT); CHKERRQ(ierr); 4666c5df90dSjeremylt } 4676c5df90dSjeremylt ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 46861608365Sjeremylt ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]); 4696c5df90dSjeremylt CHKERRQ(ierr); 4706c5df90dSjeremylt 4716c5df90dSjeremylt // Set up PCMG 4726c5df90dSjeremylt ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 4736c5df90dSjeremylt PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V; 4746c5df90dSjeremylt { 4756c5df90dSjeremylt ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 4766c5df90dSjeremylt 4776c5df90dSjeremylt // PCMG levels 4786c5df90dSjeremylt ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr); 4796c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 4806c5df90dSjeremylt // Smoother 4816c5df90dSjeremylt KSP smoother; 4826c5df90dSjeremylt PC smoother_pc; 4836c5df90dSjeremylt ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 4846c5df90dSjeremylt ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 4856c5df90dSjeremylt ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 4866c5df90dSjeremylt ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 4876c5df90dSjeremylt ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr); 4886c5df90dSjeremylt ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 4896c5df90dSjeremylt ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 4906c5df90dSjeremylt ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 4916c5df90dSjeremylt 4926c5df90dSjeremylt // Work vector 4936c5df90dSjeremylt if (i < numlevels - 1) { 4946c5df90dSjeremylt ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 4956c5df90dSjeremylt } 4966c5df90dSjeremylt 4976c5df90dSjeremylt // Level transfers 4986c5df90dSjeremylt if (i > 0) { 4996c5df90dSjeremylt // Interpolation 500a97643b0Sjeremylt ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr); 5016c5df90dSjeremylt } 5026c5df90dSjeremylt 5036c5df90dSjeremylt // Coarse solve 5046c5df90dSjeremylt KSP coarse; 5056c5df90dSjeremylt PC coarse_pc; 5066c5df90dSjeremylt ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 50715ce0ef0Sjeremylt ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr); 508199551f5Sjeremylt ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr); 50915ce0ef0Sjeremylt 5106c5df90dSjeremylt ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 51115ce0ef0Sjeremylt ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr); 51215ce0ef0Sjeremylt 51315ce0ef0Sjeremylt ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr); 51415ce0ef0Sjeremylt ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr); 51515ce0ef0Sjeremylt ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr); 51615ce0ef0Sjeremylt ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr); 5176c5df90dSjeremylt } 5186c5df90dSjeremylt 5196c5df90dSjeremylt // PCMG options 5206c5df90dSjeremylt ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 5216c5df90dSjeremylt ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 5226c5df90dSjeremylt ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr); 5236c5df90dSjeremylt } 5246c5df90dSjeremylt 5256c5df90dSjeremylt // First run, if benchmarking 5266c5df90dSjeremylt if (benchmark_mode) { 5276c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 5286c5df90dSjeremylt CHKERRQ(ierr); 52961608365Sjeremylt ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 5306c5df90dSjeremylt my_rt_start = MPI_Wtime(); 53161608365Sjeremylt ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 5326c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 5336c5df90dSjeremylt ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 5346c5df90dSjeremylt CHKERRQ(ierr); 5356c5df90dSjeremylt // Set maxits based on first iteration timing 5366c5df90dSjeremylt if (my_rt > 0.02) { 5376c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 5386c5df90dSjeremylt CHKERRQ(ierr); 5396c5df90dSjeremylt } else { 5406c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 5416c5df90dSjeremylt CHKERRQ(ierr); 5426c5df90dSjeremylt } 5436c5df90dSjeremylt } 5446c5df90dSjeremylt 5456c5df90dSjeremylt // Timed solve 54661608365Sjeremylt ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 5476c5df90dSjeremylt ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 5486c5df90dSjeremylt my_rt_start = MPI_Wtime(); 54961608365Sjeremylt ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 5506c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 5516c5df90dSjeremylt 5526c5df90dSjeremylt // Output results 5536c5df90dSjeremylt { 5546c5df90dSjeremylt KSPType ksptype; 5556c5df90dSjeremylt PCMGType pcmgtype; 5566c5df90dSjeremylt KSPConvergedReason reason; 5576c5df90dSjeremylt PetscReal rnorm; 5586c5df90dSjeremylt PetscInt its; 5596c5df90dSjeremylt ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 5606c5df90dSjeremylt ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 5616c5df90dSjeremylt ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 5626c5df90dSjeremylt ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 5636c5df90dSjeremylt ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr); 5646c5df90dSjeremylt if (!test_mode || reason < 0 || rnorm > 1e-8) { 5656c5df90dSjeremylt ierr = PetscPrintf(comm, 5666c5df90dSjeremylt " KSP:\n" 5676c5df90dSjeremylt " KSP Type : %s\n" 5686c5df90dSjeremylt " KSP Convergence : %s\n" 5696c5df90dSjeremylt " Total KSP Iterations : %D\n" 5706c5df90dSjeremylt " Final rnorm : %e\n", 5716c5df90dSjeremylt ksptype, KSPConvergedReasons[reason], its, 5726c5df90dSjeremylt (double)rnorm); CHKERRQ(ierr); 5736c5df90dSjeremylt ierr = PetscPrintf(comm, 5746c5df90dSjeremylt " PCMG:\n" 5756c5df90dSjeremylt " PCMG Type : %s\n" 5766c5df90dSjeremylt " PCMG Cycle Type : %s\n", 5776c5df90dSjeremylt PCMGTypes[pcmgtype], 5786c5df90dSjeremylt PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr); 5796c5df90dSjeremylt } 5806c5df90dSjeremylt if (!test_mode) { 5816c5df90dSjeremylt ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 5826c5df90dSjeremylt } 5836c5df90dSjeremylt { 5846c5df90dSjeremylt PetscReal maxerror; 585c70bd2a0Sjeremylt ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target, 5866c5df90dSjeremylt &maxerror); CHKERRQ(ierr); 5876c5df90dSjeremylt PetscReal tol = 5e-2; 5886c5df90dSjeremylt if (!test_mode || maxerror > tol) { 5896c5df90dSjeremylt ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 5906c5df90dSjeremylt CHKERRQ(ierr); 5916c5df90dSjeremylt ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 5926c5df90dSjeremylt CHKERRQ(ierr); 5936c5df90dSjeremylt ierr = PetscPrintf(comm, 5946c5df90dSjeremylt " Pointwise Error (max) : %e\n" 5956c5df90dSjeremylt " CG Solve Time : %g (%g) sec\n", 5966c5df90dSjeremylt (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 5976c5df90dSjeremylt } 5986c5df90dSjeremylt } 5996c5df90dSjeremylt if (benchmark_mode && (!test_mode)) { 6006c5df90dSjeremylt ierr = PetscPrintf(comm, 6016c5df90dSjeremylt " DoFs/Sec in CG : %g (%g) million\n", 60261608365Sjeremylt 1e-6*gsize[fineLevel]*its/rt_max, 60361608365Sjeremylt 1e-6*gsize[fineLevel]*its/rt_min); 6046c5df90dSjeremylt CHKERRQ(ierr); 6056c5df90dSjeremylt } 6066c5df90dSjeremylt } 6076c5df90dSjeremylt 6086c5df90dSjeremylt if (write_solution) { 6096c5df90dSjeremylt PetscViewer vtkviewersoln; 6106c5df90dSjeremylt 6116c5df90dSjeremylt ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 6126c5df90dSjeremylt ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 6136c5df90dSjeremylt ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 61461608365Sjeremylt ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr); 6156c5df90dSjeremylt ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 6166c5df90dSjeremylt } 6176c5df90dSjeremylt 6186c5df90dSjeremylt // Cleanup 6196c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 6206c5df90dSjeremylt ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 6216c5df90dSjeremylt ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr); 6226c5df90dSjeremylt ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 6236c5df90dSjeremylt ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr); 6246c5df90dSjeremylt ierr = MatDestroy(&matO[i]); CHKERRQ(ierr); 6256c5df90dSjeremylt ierr = PetscFree(userO[i]); CHKERRQ(ierr); 6266c5df90dSjeremylt if (i > 0) { 627a97643b0Sjeremylt ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr); 628a97643b0Sjeremylt ierr = PetscFree(userPR[i]); CHKERRQ(ierr); 6296c5df90dSjeremylt } 6306c5df90dSjeremylt ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr); 6316c5df90dSjeremylt ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 6326c5df90dSjeremylt } 6336c5df90dSjeremylt ierr = PetscFree(leveldegrees); CHKERRQ(ierr); 6346c5df90dSjeremylt ierr = PetscFree(dm); CHKERRQ(ierr); 6356c5df90dSjeremylt ierr = PetscFree(X); CHKERRQ(ierr); 6366c5df90dSjeremylt ierr = PetscFree(Xloc); CHKERRQ(ierr); 6376c5df90dSjeremylt ierr = PetscFree(mult); CHKERRQ(ierr); 6386c5df90dSjeremylt ierr = PetscFree(matO); CHKERRQ(ierr); 639a97643b0Sjeremylt ierr = PetscFree(matPR); CHKERRQ(ierr); 6406c5df90dSjeremylt ierr = PetscFree(ceeddata); CHKERRQ(ierr); 6416c5df90dSjeremylt ierr = PetscFree(userO); CHKERRQ(ierr); 642a97643b0Sjeremylt ierr = PetscFree(userPR); CHKERRQ(ierr); 6436c5df90dSjeremylt ierr = PetscFree(lsize); CHKERRQ(ierr); 6446c5df90dSjeremylt ierr = PetscFree(xlsize); CHKERRQ(ierr); 6456c5df90dSjeremylt ierr = PetscFree(gsize); CHKERRQ(ierr); 6466c5df90dSjeremylt ierr = VecDestroy(&rhs); CHKERRQ(ierr); 6476c5df90dSjeremylt ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 648199551f5Sjeremylt ierr = MatDestroy(&matcoarse); CHKERRQ(ierr); 6496c5df90dSjeremylt ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 650c70bd2a0Sjeremylt ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr); 651199551f5Sjeremylt ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 6526c5df90dSjeremylt CeedVectorDestroy(&target); 653c70bd2a0Sjeremylt CeedQFunctionDestroy(&qferror); 654c70bd2a0Sjeremylt CeedQFunctionDestroy(&qfrestrict); 655c70bd2a0Sjeremylt CeedQFunctionDestroy(&qfprolong); 656c70bd2a0Sjeremylt CeedOperatorDestroy(&operror); 6576c5df90dSjeremylt CeedDestroy(&ceed); 6586c5df90dSjeremylt return PetscFinalize(); 6596c5df90dSjeremylt } 660