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 3128688798Sjeremylt // multigrid -problem bp4 3228688798Sjeremylt // multigrid -problem bp5 -ceed /cpu/self 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; 5409a940d7Sjeremylt PetscLogStage solvestage; 55199551f5Sjeremylt DM *dm, dmorig; 56c70bd2a0Sjeremylt SNES snesdummy; 576c5df90dSjeremylt KSP ksp; 586c5df90dSjeremylt PC pc; 59199551f5Sjeremylt Mat *matO, *matPR, matcoarse; 60cdb3667fSjeremylt Vec *X, *Xloc, *mult, rhs, rhsloc; 61b68a8d79SJed Brown PetscMemType memtype; 626c5df90dSjeremylt UserO *userO; 63a97643b0Sjeremylt UserProlongRestr *userPR; 646c5df90dSjeremylt Ceed ceed; 656c5df90dSjeremylt CeedData *ceeddata; 665e81177dSjeremylt CeedVector rhsceed, target; 67c70bd2a0Sjeremylt CeedQFunction qferror, qfrestrict, qfprolong; 68c70bd2a0Sjeremylt CeedOperator operror; 69199551f5Sjeremylt bpType bpchoice; 706c5df90dSjeremylt coarsenType coarsen; 716c5df90dSjeremylt 726c5df90dSjeremylt ierr = PetscInitialize(&argc, &argv, NULL, help); 736c5df90dSjeremylt if (ierr) return ierr; 746c5df90dSjeremylt comm = PETSC_COMM_WORLD; 756c5df90dSjeremylt 766c5df90dSjeremylt // Parse command line options 776c5df90dSjeremylt ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 78199551f5Sjeremylt bpchoice = CEED_BP3; 796c5df90dSjeremylt ierr = PetscOptionsEnum("-problem", 806c5df90dSjeremylt "CEED benchmark problem to solve", NULL, 81199551f5Sjeremylt bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice, 826c5df90dSjeremylt NULL); CHKERRQ(ierr); 83199551f5Sjeremylt ncompu = bpOptions[bpchoice].ncompu; 846c5df90dSjeremylt test_mode = PETSC_FALSE; 856c5df90dSjeremylt ierr = PetscOptionsBool("-test", 866c5df90dSjeremylt "Testing mode (do not print unless error is large)", 876c5df90dSjeremylt NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 886c5df90dSjeremylt benchmark_mode = PETSC_FALSE; 896c5df90dSjeremylt ierr = PetscOptionsBool("-benchmark", 906c5df90dSjeremylt "Benchmarking mode (prints benchmark statistics)", 916c5df90dSjeremylt NULL, benchmark_mode, &benchmark_mode, NULL); 926c5df90dSjeremylt CHKERRQ(ierr); 936c5df90dSjeremylt write_solution = PETSC_FALSE; 946c5df90dSjeremylt ierr = PetscOptionsBool("-write_solution", 956c5df90dSjeremylt "Write solution for visualization", 966c5df90dSjeremylt NULL, write_solution, &write_solution, NULL); 976c5df90dSjeremylt CHKERRQ(ierr); 986c5df90dSjeremylt degree = test_mode ? 3 : 2; 996c5df90dSjeremylt ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 1006c5df90dSjeremylt NULL, degree, °ree, NULL); CHKERRQ(ierr); 1016c5df90dSjeremylt if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 1026c5df90dSjeremylt "-degree %D must be at least 1", degree); 103199551f5Sjeremylt qextra = bpOptions[bpchoice].qextra; 1046c5df90dSjeremylt ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 1056c5df90dSjeremylt NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 1066c5df90dSjeremylt ierr = PetscOptionsString("-ceed", "CEED resource specifier", 1076c5df90dSjeremylt NULL, ceedresource, ceedresource, 1086c5df90dSjeremylt sizeof(ceedresource), NULL); CHKERRQ(ierr); 1096c5df90dSjeremylt coarsen = COARSEN_UNIFORM; 1106c5df90dSjeremylt ierr = PetscOptionsEnum("-coarsen", 1116c5df90dSjeremylt "Coarsening strategy to use", NULL, 1126c5df90dSjeremylt coarsenTypes, (PetscEnum)coarsen, 1136c5df90dSjeremylt (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr); 114cb32e2e7SValeria Barra read_mesh = PETSC_FALSE; 1156c5df90dSjeremylt ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 1166c5df90dSjeremylt filename, filename, sizeof(filename), &read_mesh); 1176c5df90dSjeremylt CHKERRQ(ierr); 1186c5df90dSjeremylt if (!read_mesh) { 1196c5df90dSjeremylt PetscInt tmp = dim; 1206c5df90dSjeremylt ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL, 1216c5df90dSjeremylt melem, &tmp, NULL); CHKERRQ(ierr); 1226c5df90dSjeremylt } 1236c5df90dSjeremylt ierr = PetscOptionsEnd(); CHKERRQ(ierr); 1246c5df90dSjeremylt 1259396343dSjeremylt // Set up libCEED 1269396343dSjeremylt CeedInit(ceedresource, &ceed); 1279396343dSjeremylt CeedMemType memtypebackend; 1289396343dSjeremylt CeedGetPreferredMemType(ceed, &memtypebackend); 1299396343dSjeremylt 1306c5df90dSjeremylt // Setup DM 1316c5df90dSjeremylt if (read_mesh) { 132199551f5Sjeremylt ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmorig); 1336c5df90dSjeremylt CHKERRQ(ierr); 1346c5df90dSjeremylt } else { 1356c5df90dSjeremylt ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL, 136199551f5Sjeremylt NULL, NULL, PETSC_TRUE,&dmorig); CHKERRQ(ierr); 1376c5df90dSjeremylt } 1386c5df90dSjeremylt 1396c5df90dSjeremylt { 1406c5df90dSjeremylt DM dmDist = NULL; 1416c5df90dSjeremylt PetscPartitioner part; 1426c5df90dSjeremylt 143199551f5Sjeremylt ierr = DMPlexGetPartitioner(dmorig, &part); CHKERRQ(ierr); 1446c5df90dSjeremylt ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 145199551f5Sjeremylt ierr = DMPlexDistribute(dmorig, 0, NULL, &dmDist); CHKERRQ(ierr); 1466c5df90dSjeremylt if (dmDist) { 147199551f5Sjeremylt ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 148199551f5Sjeremylt dmorig = dmDist; 1496c5df90dSjeremylt } 1506c5df90dSjeremylt } 1516c5df90dSjeremylt 152b68a8d79SJed Brown VecType vectype; 153b68a8d79SJed Brown switch (memtypebackend) { 154b68a8d79SJed Brown case CEED_MEM_HOST: vectype = VECSTANDARD; break; 155b68a8d79SJed Brown case CEED_MEM_DEVICE: { 156b68a8d79SJed Brown const char *resolved; 157b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 158b68a8d79SJed Brown if (strstr(resolved, "/gpu/cuda")) vectype = VECCUDA; 159*ca2d516cSJed Brown else if (strstr(resolved, "/gpu/hip/occa")) 160*ca2d516cSJed Brown vectype = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 161b68a8d79SJed Brown else if (strstr(resolved, "/gpu/hip")) vectype = VECHIP; 162b68a8d79SJed Brown else vectype = VECSTANDARD; 163b68a8d79SJed Brown } 164b68a8d79SJed Brown } 165b68a8d79SJed Brown ierr = DMSetVecType(dmorig, vectype); CHKERRQ(ierr); 166b68a8d79SJed Brown ierr = DMSetFromOptions(dmorig); CHKERRQ(ierr); 167b68a8d79SJed Brown 1686c5df90dSjeremylt // Allocate arrays for PETSc objects for each level 1696c5df90dSjeremylt switch (coarsen) { 1706c5df90dSjeremylt case COARSEN_UNIFORM: 1716c5df90dSjeremylt numlevels = degree; 1726c5df90dSjeremylt break; 173dc7d240cSValeria Barra case COARSEN_LOGARITHMIC: 1746c5df90dSjeremylt numlevels = ceil(log(degree)/log(2)) + 1; 1756c5df90dSjeremylt break; 1766c5df90dSjeremylt } 1776c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr); 17861608365Sjeremylt fineLevel = numlevels - 1; 17961608365Sjeremylt 1806c5df90dSjeremylt switch (coarsen) { 1816c5df90dSjeremylt case COARSEN_UNIFORM: 1826c5df90dSjeremylt for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1; 1836c5df90dSjeremylt break; 184dc7d240cSValeria Barra case COARSEN_LOGARITHMIC: 1856c5df90dSjeremylt for (int i=0; i<numlevels - 1; i++) leveldegrees[i] = pow(2,i); 18661608365Sjeremylt leveldegrees[fineLevel] = degree; 1876c5df90dSjeremylt break; 1886c5df90dSjeremylt } 1896c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr); 1906c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr); 1916c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr); 1926c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr); 1936c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr); 194a97643b0Sjeremylt ierr = PetscMalloc1(numlevels, &userPR); CHKERRQ(ierr); 1956c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr); 196a97643b0Sjeremylt ierr = PetscMalloc1(numlevels, &matPR); CHKERRQ(ierr); 1976c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr); 1986c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr); 1996c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr); 2006c5df90dSjeremylt 2016c5df90dSjeremylt // Setup DM and Operator Mat Shells for each level 2026c5df90dSjeremylt for (CeedInt i=0; i<numlevels; i++) { 2036c5df90dSjeremylt // Create DM 204199551f5Sjeremylt ierr = DMClone(dmorig, &dm[i]); CHKERRQ(ierr); 205b68a8d79SJed Brown ierr = DMGetVecType(dmorig, &vectype); CHKERRQ(ierr); 206b68a8d79SJed Brown ierr = DMSetVecType(dm[i], vectype); CHKERRQ(ierr); 207199551f5Sjeremylt ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpchoice); 2086c5df90dSjeremylt CHKERRQ(ierr); 2096c5df90dSjeremylt 2106c5df90dSjeremylt // Create vectors 2116c5df90dSjeremylt ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr); 2126c5df90dSjeremylt ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr); 2136c5df90dSjeremylt ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr); 2146c5df90dSjeremylt ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr); 2156c5df90dSjeremylt ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr); 2166c5df90dSjeremylt 2176c5df90dSjeremylt // Operator 2186c5df90dSjeremylt ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr); 2196c5df90dSjeremylt ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i], 2206c5df90dSjeremylt userO[i], &matO[i]); CHKERRQ(ierr); 2216c5df90dSjeremylt ierr = MatShellSetOperation(matO[i], MATOP_MULT, 222ce74dcefSjeremylt (void(*)(void))MatMult_Ceed); CHKERRQ(ierr); 2236c5df90dSjeremylt ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL, 224ce74dcefSjeremylt (void(*)(void))MatGetDiag); CHKERRQ(ierr); 225b68a8d79SJed Brown ierr = MatShellSetVecType(matO[i], vectype); CHKERRQ(ierr); 2266c5df90dSjeremylt 2276c5df90dSjeremylt // Level transfers 2286c5df90dSjeremylt if (i > 0) { 2296c5df90dSjeremylt // Interp 230a97643b0Sjeremylt ierr = PetscMalloc1(1, &userPR[i]); CHKERRQ(ierr); 2316c5df90dSjeremylt ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1], 232a97643b0Sjeremylt userPR[i], &matPR[i]); CHKERRQ(ierr); 233a97643b0Sjeremylt ierr = MatShellSetOperation(matPR[i], MATOP_MULT, 234a97643b0Sjeremylt (void(*)(void))MatMult_Prolong); 2356c5df90dSjeremylt CHKERRQ(ierr); 236a97643b0Sjeremylt ierr = MatShellSetOperation(matPR[i], MATOP_MULT_TRANSPOSE, 2376c5df90dSjeremylt (void(*)(void))MatMult_Restrict); 2386c5df90dSjeremylt CHKERRQ(ierr); 239b68a8d79SJed Brown ierr = MatShellSetVecType(matPR[i], vectype); CHKERRQ(ierr); 2406c5df90dSjeremylt } 2416c5df90dSjeremylt } 24261608365Sjeremylt ierr = VecDuplicate(X[fineLevel], &rhs); CHKERRQ(ierr); 2436c5df90dSjeremylt 2446c5df90dSjeremylt // Print global grid information 2456c5df90dSjeremylt if (!test_mode) { 2466c5df90dSjeremylt PetscInt P = degree + 1, Q = P + qextra; 2479396343dSjeremylt 2482d03409cSjeremylt const char *usedresource; 2492d03409cSjeremylt CeedGetResource(ceed, &usedresource); 2509396343dSjeremylt 2519396343dSjeremylt ierr = VecGetType(X[0], &vectype); CHKERRQ(ierr); 2529396343dSjeremylt 2536c5df90dSjeremylt ierr = PetscPrintf(comm, 2546c5df90dSjeremylt "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n" 2559396343dSjeremylt " PETSc:\n" 2569396343dSjeremylt " PETSc Vec Type : %s\n" 2576c5df90dSjeremylt " libCEED:\n" 2586c5df90dSjeremylt " libCEED Backend : %s\n" 2599396343dSjeremylt " libCEED Backend MemType : %s\n" 2606c5df90dSjeremylt " Mesh:\n" 2616c5df90dSjeremylt " Number of 1D Basis Nodes (p) : %d\n" 2626c5df90dSjeremylt " Number of 1D Quadrature Points (q) : %d\n" 2636c5df90dSjeremylt " Global Nodes : %D\n" 2646c5df90dSjeremylt " Owned Nodes : %D\n" 265db419314Sjeremylt " DoF per node : %D\n" 2666c5df90dSjeremylt " Multigrid:\n" 2676c5df90dSjeremylt " Number of Levels : %d\n", 2689396343dSjeremylt bpchoice+1, vectype, usedresource, 2699396343dSjeremylt CeedMemTypes[memtypebackend], 2709396343dSjeremylt P, Q, gsize[fineLevel]/ncompu, lsize[fineLevel]/ncompu, 271db419314Sjeremylt ncompu, numlevels); CHKERRQ(ierr); 2726c5df90dSjeremylt } 2736c5df90dSjeremylt 2746c5df90dSjeremylt // Create RHS vector 27561608365Sjeremylt ierr = VecDuplicate(Xloc[fineLevel], &rhsloc); CHKERRQ(ierr); 2766c5df90dSjeremylt ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 277b68a8d79SJed Brown ierr = VecGetArrayAndMemType(rhsloc, &r, &memtype); CHKERRQ(ierr); 27861608365Sjeremylt CeedVectorCreate(ceed, xlsize[fineLevel], &rhsceed); 279b68a8d79SJed Brown CeedVectorSetArray(rhsceed, MemTypeP2C(memtype), CEED_USE_POINTER, r); 2806c5df90dSjeremylt 2816c5df90dSjeremylt // Set up libCEED operators on each level 2826c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr); 2836c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 2846c5df90dSjeremylt // Print level information 28561608365Sjeremylt if (!test_mode && (i == 0 || i == fineLevel)) { 2866c5df90dSjeremylt ierr = PetscPrintf(comm," Level %D (%s):\n" 2876c5df90dSjeremylt " Number of 1D Basis Nodes (p) : %d\n" 2886c5df90dSjeremylt " Global Nodes : %D\n" 2896c5df90dSjeremylt " Owned Nodes : %D\n", 2906c5df90dSjeremylt i, (i? "fine" : "coarse"), leveldegrees[i] + 1, 2916c5df90dSjeremylt gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr); 2926c5df90dSjeremylt } 2936c5df90dSjeremylt ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr); 2946c5df90dSjeremylt ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra, 295199551f5Sjeremylt ncompu, gsize[i], xlsize[i], bpchoice, 296f9342789SJeremy L Thompson ceeddata[i], i==(fineLevel), rhsceed, &target); 297f9342789SJeremy L Thompson CHKERRQ(ierr); 2986c5df90dSjeremylt } 2996c5df90dSjeremylt 3006c5df90dSjeremylt // Gather RHS 301b68a8d79SJed Brown CeedVectorTakeArray(rhsceed, MemTypeP2C(memtype), NULL); 302b68a8d79SJed Brown ierr = VecRestoreArrayAndMemType(rhsloc, &r); CHKERRQ(ierr); 3036c5df90dSjeremylt ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 3049396343dSjeremylt ierr = DMLocalToGlobal(dm[fineLevel], rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr); 3056c5df90dSjeremylt CeedVectorDestroy(&rhsceed); 3066c5df90dSjeremylt 30743eb8658SJeremy L Thompson // Create the restriction/interpolation QFunction 30860f77c51Sjeremylt CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_NONE, CEED_EVAL_INTERP, 309c70bd2a0Sjeremylt &qfrestrict); 31060f77c51Sjeremylt CeedQFunctionCreateIdentity(ceed, ncompu, CEED_EVAL_INTERP, CEED_EVAL_NONE, 311c70bd2a0Sjeremylt &qfprolong); 3126c5df90dSjeremylt 3136c5df90dSjeremylt // Set up libCEED level transfer operators 314199551f5Sjeremylt ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpchoice, ceeddata, 315c70bd2a0Sjeremylt leveldegrees, qfrestrict, qfprolong); 3166c5df90dSjeremylt CHKERRQ(ierr); 3176c5df90dSjeremylt 31843eb8658SJeremy L Thompson // Create the error QFunction 319199551f5Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error, 320199551f5Sjeremylt bpOptions[bpchoice].errorfname, &qferror); 321c70bd2a0Sjeremylt CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP); 322c70bd2a0Sjeremylt CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE); 323c70bd2a0Sjeremylt CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE); 3246c5df90dSjeremylt 3256c5df90dSjeremylt // Create the error operator 326c70bd2a0Sjeremylt CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 327c70bd2a0Sjeremylt &operror); 328c70bd2a0Sjeremylt CeedOperatorSetField(operror, "u", ceeddata[fineLevel]->Erestrictu, 32961608365Sjeremylt ceeddata[fineLevel]->basisu, CEED_VECTOR_ACTIVE); 330c70bd2a0Sjeremylt CeedOperatorSetField(operror, "true_soln", ceeddata[fineLevel]->Erestrictui, 331a8d32208Sjeremylt CEED_BASIS_COLLOCATED, target); 332c70bd2a0Sjeremylt CeedOperatorSetField(operror, "error", ceeddata[fineLevel]->Erestrictui, 333a8d32208Sjeremylt CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 3346c5df90dSjeremylt 3356c5df90dSjeremylt // Calculate multiplicity 3366c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 3376c5df90dSjeremylt PetscScalar *x; 3386c5df90dSjeremylt 3396c5df90dSjeremylt // CEED vector 340f9342789SJeremy L Thompson ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 3416c5df90dSjeremylt ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr); 3426c5df90dSjeremylt CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 3436c5df90dSjeremylt 3446c5df90dSjeremylt // Multiplicity 345a8d32208Sjeremylt CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu, 3466c5df90dSjeremylt ceeddata[i]->xceed); 347f9342789SJeremy L Thompson CeedVectorSyncArray(ceeddata[i]->xceed, CEED_MEM_HOST); 3486c5df90dSjeremylt 3496c5df90dSjeremylt // Restore vector 3506c5df90dSjeremylt ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr); 3516c5df90dSjeremylt 3526c5df90dSjeremylt // Creat mult vector 3536c5df90dSjeremylt ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr); 3546c5df90dSjeremylt 3556c5df90dSjeremylt // Local-to-global 3566c5df90dSjeremylt ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 357483f8b0dSjeremylt ierr = DMLocalToGlobal(dm[i], Xloc[i], ADD_VALUES, X[i]); 3586c5df90dSjeremylt CHKERRQ(ierr); 3596c5df90dSjeremylt ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 3606c5df90dSjeremylt 3616c5df90dSjeremylt // Global-to-local 362483f8b0dSjeremylt ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]); 3636c5df90dSjeremylt CHKERRQ(ierr); 3646c5df90dSjeremylt ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 3656c5df90dSjeremylt 3666c5df90dSjeremylt // Multiplicity scaling 3676c5df90dSjeremylt ierr = VecReciprocal(mult[i]); 3686c5df90dSjeremylt } 3696c5df90dSjeremylt 3706c5df90dSjeremylt // Set up Mat 3716c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 372226c3a8fSjeremylt // User Operator 3736c5df90dSjeremylt userO[i]->comm = comm; 3746c5df90dSjeremylt userO[i]->dm = dm[i]; 3756c5df90dSjeremylt userO[i]->Xloc = Xloc[i]; 3766c5df90dSjeremylt ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr); 3776c5df90dSjeremylt userO[i]->xceed = ceeddata[i]->xceed; 3786c5df90dSjeremylt userO[i]->yceed = ceeddata[i]->yceed; 379c70bd2a0Sjeremylt userO[i]->op = ceeddata[i]->opapply; 3806c5df90dSjeremylt userO[i]->ceed = ceed; 3816c5df90dSjeremylt 3826c5df90dSjeremylt if (i > 0) { 383a97643b0Sjeremylt // Prolongation/Restriction Operator 384a97643b0Sjeremylt userPR[i]->comm = comm; 385199551f5Sjeremylt userPR[i]->dmf = dm[i]; 386199551f5Sjeremylt userPR[i]->dmc = dm[i-1]; 387199551f5Sjeremylt userPR[i]->locvecc = Xloc[i-1]; 388199551f5Sjeremylt userPR[i]->locvecf = userO[i]->Yloc; 389199551f5Sjeremylt userPR[i]->multvec = mult[i]; 390199551f5Sjeremylt userPR[i]->ceedvecc = userO[i-1]->xceed; 391199551f5Sjeremylt userPR[i]->ceedvecf = userO[i]->yceed; 392c70bd2a0Sjeremylt userPR[i]->opprolong = ceeddata[i]->opprolong; 393c70bd2a0Sjeremylt userPR[i]->oprestrict = ceeddata[i]->oprestrict; 394a97643b0Sjeremylt userPR[i]->ceed = ceed; 3956c5df90dSjeremylt } 3966c5df90dSjeremylt } 3976c5df90dSjeremylt 39815ce0ef0Sjeremylt // Setup dummy SNES for AMG coarse solve 399c70bd2a0Sjeremylt ierr = SNESCreate(comm, &snesdummy); CHKERRQ(ierr); 400c70bd2a0Sjeremylt ierr = SNESSetDM(snesdummy, dm[0]); CHKERRQ(ierr); 401c70bd2a0Sjeremylt ierr = SNESSetSolution(snesdummy, X[0]); CHKERRQ(ierr); 40215ce0ef0Sjeremylt 40315ce0ef0Sjeremylt // -- Jacobian matrix 40415ce0ef0Sjeremylt ierr = DMSetMatType(dm[0], MATAIJ); CHKERRQ(ierr); 405199551f5Sjeremylt ierr = DMCreateMatrix(dm[0], &matcoarse); CHKERRQ(ierr); 406199551f5Sjeremylt ierr = SNESSetJacobian(snesdummy, matcoarse, matcoarse, NULL, 40715ce0ef0Sjeremylt NULL); CHKERRQ(ierr); 40815ce0ef0Sjeremylt 40915ce0ef0Sjeremylt // -- Residual evaluation function 410c70bd2a0Sjeremylt ierr = SNESSetFunction(snesdummy, X[0], FormResidual_Ceed, 41115ce0ef0Sjeremylt userO[0]); CHKERRQ(ierr); 41215ce0ef0Sjeremylt 41315ce0ef0Sjeremylt // -- Form Jacobian 414c70bd2a0Sjeremylt ierr = SNESComputeJacobianDefaultColor(snesdummy, X[0], matO[0], 415199551f5Sjeremylt matcoarse, NULL); CHKERRQ(ierr); 41615ce0ef0Sjeremylt 4176c5df90dSjeremylt // Set up KSP 4186c5df90dSjeremylt ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 4196c5df90dSjeremylt { 4206c5df90dSjeremylt ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 4216c5df90dSjeremylt ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 4226c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 4236c5df90dSjeremylt PETSC_DEFAULT); CHKERRQ(ierr); 4246c5df90dSjeremylt } 4256c5df90dSjeremylt ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 42661608365Sjeremylt ierr = KSPSetOperators(ksp, matO[fineLevel], matO[fineLevel]); 4276c5df90dSjeremylt CHKERRQ(ierr); 4286c5df90dSjeremylt 4296c5df90dSjeremylt // Set up PCMG 4306c5df90dSjeremylt ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 4316c5df90dSjeremylt PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V; 4326c5df90dSjeremylt { 4336c5df90dSjeremylt ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 4346c5df90dSjeremylt 4356c5df90dSjeremylt // PCMG levels 4366c5df90dSjeremylt ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr); 4376c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 4386c5df90dSjeremylt // Smoother 4396c5df90dSjeremylt KSP smoother; 4406c5df90dSjeremylt PC smoother_pc; 4416c5df90dSjeremylt ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 4426c5df90dSjeremylt ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 4436c5df90dSjeremylt ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 4446c5df90dSjeremylt ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 4456c5df90dSjeremylt ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr); 4466c5df90dSjeremylt ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 4476c5df90dSjeremylt ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 4486c5df90dSjeremylt ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 4496c5df90dSjeremylt 4506c5df90dSjeremylt // Work vector 4516c5df90dSjeremylt if (i < numlevels - 1) { 4526c5df90dSjeremylt ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 4536c5df90dSjeremylt } 4546c5df90dSjeremylt 4556c5df90dSjeremylt // Level transfers 4566c5df90dSjeremylt if (i > 0) { 4576c5df90dSjeremylt // Interpolation 458a97643b0Sjeremylt ierr = PCMGSetInterpolation(pc, i, matPR[i]); CHKERRQ(ierr); 4596c5df90dSjeremylt } 4606c5df90dSjeremylt 4616c5df90dSjeremylt // Coarse solve 4626c5df90dSjeremylt KSP coarse; 4636c5df90dSjeremylt PC coarse_pc; 4646c5df90dSjeremylt ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 46515ce0ef0Sjeremylt ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr); 466199551f5Sjeremylt ierr = KSPSetOperators(coarse, matcoarse, matcoarse); CHKERRQ(ierr); 46715ce0ef0Sjeremylt 4686c5df90dSjeremylt ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 46915ce0ef0Sjeremylt ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr); 47015ce0ef0Sjeremylt 47115ce0ef0Sjeremylt ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr); 47215ce0ef0Sjeremylt ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr); 47315ce0ef0Sjeremylt ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr); 47415ce0ef0Sjeremylt ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr); 4756c5df90dSjeremylt } 4766c5df90dSjeremylt 4776c5df90dSjeremylt // PCMG options 4786c5df90dSjeremylt ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 4796c5df90dSjeremylt ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 4806c5df90dSjeremylt ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr); 4816c5df90dSjeremylt } 4826c5df90dSjeremylt 4836c5df90dSjeremylt // First run, if benchmarking 4846c5df90dSjeremylt if (benchmark_mode) { 4856c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 4866c5df90dSjeremylt CHKERRQ(ierr); 48761608365Sjeremylt ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 4886c5df90dSjeremylt my_rt_start = MPI_Wtime(); 48961608365Sjeremylt ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 4906c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 4916c5df90dSjeremylt ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 4926c5df90dSjeremylt CHKERRQ(ierr); 4936c5df90dSjeremylt // Set maxits based on first iteration timing 4946c5df90dSjeremylt if (my_rt > 0.02) { 4956c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 4966c5df90dSjeremylt CHKERRQ(ierr); 4976c5df90dSjeremylt } else { 4986c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 4996c5df90dSjeremylt CHKERRQ(ierr); 5006c5df90dSjeremylt } 5016c5df90dSjeremylt } 5026c5df90dSjeremylt 5036c5df90dSjeremylt // Timed solve 50461608365Sjeremylt ierr = VecZeroEntries(X[fineLevel]); CHKERRQ(ierr); 5056c5df90dSjeremylt ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 50609a940d7Sjeremylt 50709a940d7Sjeremylt // -- Performance logging 50809a940d7Sjeremylt ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr); 50909a940d7Sjeremylt ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr); 51009a940d7Sjeremylt 51109a940d7Sjeremylt // -- Solve 5126c5df90dSjeremylt my_rt_start = MPI_Wtime(); 51361608365Sjeremylt ierr = KSPSolve(ksp, rhs, X[fineLevel]); CHKERRQ(ierr); 5146c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 5156c5df90dSjeremylt 51609a940d7Sjeremylt 51709a940d7Sjeremylt // -- Performance logging 51809a940d7Sjeremylt ierr = PetscLogStagePop(); 51909a940d7Sjeremylt 5206c5df90dSjeremylt // Output results 5216c5df90dSjeremylt { 5226c5df90dSjeremylt KSPType ksptype; 5236c5df90dSjeremylt PCMGType pcmgtype; 5246c5df90dSjeremylt KSPConvergedReason reason; 5256c5df90dSjeremylt PetscReal rnorm; 5266c5df90dSjeremylt PetscInt its; 5276c5df90dSjeremylt ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 5286c5df90dSjeremylt ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 5296c5df90dSjeremylt ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 5306c5df90dSjeremylt ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 5316c5df90dSjeremylt ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr); 5326c5df90dSjeremylt if (!test_mode || reason < 0 || rnorm > 1e-8) { 5336c5df90dSjeremylt ierr = PetscPrintf(comm, 5346c5df90dSjeremylt " KSP:\n" 5356c5df90dSjeremylt " KSP Type : %s\n" 5366c5df90dSjeremylt " KSP Convergence : %s\n" 5376c5df90dSjeremylt " Total KSP Iterations : %D\n" 5386c5df90dSjeremylt " Final rnorm : %e\n", 5396c5df90dSjeremylt ksptype, KSPConvergedReasons[reason], its, 5406c5df90dSjeremylt (double)rnorm); CHKERRQ(ierr); 5416c5df90dSjeremylt ierr = PetscPrintf(comm, 5426c5df90dSjeremylt " PCMG:\n" 5436c5df90dSjeremylt " PCMG Type : %s\n" 5446c5df90dSjeremylt " PCMG Cycle Type : %s\n", 5456c5df90dSjeremylt PCMGTypes[pcmgtype], 5466c5df90dSjeremylt PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr); 5476c5df90dSjeremylt } 5486c5df90dSjeremylt if (!test_mode) { 5496c5df90dSjeremylt ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 5506c5df90dSjeremylt } 5516c5df90dSjeremylt { 5526c5df90dSjeremylt PetscReal maxerror; 553c70bd2a0Sjeremylt ierr = ComputeErrorMax(userO[fineLevel], operror, X[fineLevel], target, 5546c5df90dSjeremylt &maxerror); CHKERRQ(ierr); 5556c5df90dSjeremylt PetscReal tol = 5e-2; 5566c5df90dSjeremylt if (!test_mode || maxerror > tol) { 5576c5df90dSjeremylt ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 5586c5df90dSjeremylt CHKERRQ(ierr); 5596c5df90dSjeremylt ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 5606c5df90dSjeremylt CHKERRQ(ierr); 5616c5df90dSjeremylt ierr = PetscPrintf(comm, 5626c5df90dSjeremylt " Pointwise Error (max) : %e\n" 5636c5df90dSjeremylt " CG Solve Time : %g (%g) sec\n", 5646c5df90dSjeremylt (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 5656c5df90dSjeremylt } 5666c5df90dSjeremylt } 5676c5df90dSjeremylt if (benchmark_mode && (!test_mode)) { 5686c5df90dSjeremylt ierr = PetscPrintf(comm, 5696c5df90dSjeremylt " DoFs/Sec in CG : %g (%g) million\n", 57061608365Sjeremylt 1e-6*gsize[fineLevel]*its/rt_max, 57161608365Sjeremylt 1e-6*gsize[fineLevel]*its/rt_min); 5726c5df90dSjeremylt CHKERRQ(ierr); 5736c5df90dSjeremylt } 5746c5df90dSjeremylt } 5756c5df90dSjeremylt 5766c5df90dSjeremylt if (write_solution) { 5776c5df90dSjeremylt PetscViewer vtkviewersoln; 5786c5df90dSjeremylt 5796c5df90dSjeremylt ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 5806c5df90dSjeremylt ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 5816c5df90dSjeremylt ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 58261608365Sjeremylt ierr = VecView(X[fineLevel], vtkviewersoln); CHKERRQ(ierr); 5836c5df90dSjeremylt ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 5846c5df90dSjeremylt } 5856c5df90dSjeremylt 5866c5df90dSjeremylt // Cleanup 5876c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 5886c5df90dSjeremylt ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 5896c5df90dSjeremylt ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr); 5906c5df90dSjeremylt ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 5916c5df90dSjeremylt ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr); 5926c5df90dSjeremylt ierr = MatDestroy(&matO[i]); CHKERRQ(ierr); 5936c5df90dSjeremylt ierr = PetscFree(userO[i]); CHKERRQ(ierr); 5946c5df90dSjeremylt if (i > 0) { 595a97643b0Sjeremylt ierr = MatDestroy(&matPR[i]); CHKERRQ(ierr); 596a97643b0Sjeremylt ierr = PetscFree(userPR[i]); CHKERRQ(ierr); 5976c5df90dSjeremylt } 5986c5df90dSjeremylt ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr); 5996c5df90dSjeremylt ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 6006c5df90dSjeremylt } 6016c5df90dSjeremylt ierr = PetscFree(leveldegrees); CHKERRQ(ierr); 6026c5df90dSjeremylt ierr = PetscFree(dm); CHKERRQ(ierr); 6036c5df90dSjeremylt ierr = PetscFree(X); CHKERRQ(ierr); 6046c5df90dSjeremylt ierr = PetscFree(Xloc); CHKERRQ(ierr); 6056c5df90dSjeremylt ierr = PetscFree(mult); CHKERRQ(ierr); 6066c5df90dSjeremylt ierr = PetscFree(matO); CHKERRQ(ierr); 607a97643b0Sjeremylt ierr = PetscFree(matPR); CHKERRQ(ierr); 6086c5df90dSjeremylt ierr = PetscFree(ceeddata); CHKERRQ(ierr); 6096c5df90dSjeremylt ierr = PetscFree(userO); CHKERRQ(ierr); 610a97643b0Sjeremylt ierr = PetscFree(userPR); CHKERRQ(ierr); 6116c5df90dSjeremylt ierr = PetscFree(lsize); CHKERRQ(ierr); 6126c5df90dSjeremylt ierr = PetscFree(xlsize); CHKERRQ(ierr); 6136c5df90dSjeremylt ierr = PetscFree(gsize); CHKERRQ(ierr); 6146c5df90dSjeremylt ierr = VecDestroy(&rhs); CHKERRQ(ierr); 6156c5df90dSjeremylt ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 616199551f5Sjeremylt ierr = MatDestroy(&matcoarse); CHKERRQ(ierr); 6176c5df90dSjeremylt ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 618c70bd2a0Sjeremylt ierr = SNESDestroy(&snesdummy); CHKERRQ(ierr); 619199551f5Sjeremylt ierr = DMDestroy(&dmorig); CHKERRQ(ierr); 6206c5df90dSjeremylt CeedVectorDestroy(&target); 621c70bd2a0Sjeremylt CeedQFunctionDestroy(&qferror); 622c70bd2a0Sjeremylt CeedQFunctionDestroy(&qfrestrict); 623c70bd2a0Sjeremylt CeedQFunctionDestroy(&qfprolong); 624c70bd2a0Sjeremylt CeedOperatorDestroy(&operror); 6256c5df90dSjeremylt CeedDestroy(&ceed); 6266c5df90dSjeremylt return PetscFinalize(); 6276c5df90dSjeremylt } 628