1*6c5df90dSjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2*6c5df90dSjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3*6c5df90dSjeremylt // reserved. See files LICENSE and NOTICE for details. 4*6c5df90dSjeremylt // 5*6c5df90dSjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software 6*6c5df90dSjeremylt // libraries and APIs for efficient high-order finite element and spectral 7*6c5df90dSjeremylt // element discretizations for exascale applications. For more information and 8*6c5df90dSjeremylt // source code availability see http://github.com/ceed. 9*6c5df90dSjeremylt // 10*6c5df90dSjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11*6c5df90dSjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office 12*6c5df90dSjeremylt // of Science and the National Nuclear Security Administration) responsible for 13*6c5df90dSjeremylt // the planning and preparation of a capable exascale ecosystem, including 14*6c5df90dSjeremylt // software, applications, hardware, advanced system engineering and early 15*6c5df90dSjeremylt // testbed platforms, in support of the nation's exascale computing imperative. 16*6c5df90dSjeremylt 17*6c5df90dSjeremylt // libCEED + PETSc Example: CEED BPs 3-6 with Multigrid 18*6c5df90dSjeremylt // 19*6c5df90dSjeremylt // This example demonstrates a simple usage of libCEED with PETSc to solve the 20*6c5df90dSjeremylt // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps. 21*6c5df90dSjeremylt // 22*6c5df90dSjeremylt // The code uses higher level communication protocols in DMPlex. 23*6c5df90dSjeremylt // 24*6c5df90dSjeremylt // Build with: 25*6c5df90dSjeremylt // 26*6c5df90dSjeremylt // make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 27*6c5df90dSjeremylt // 28*6c5df90dSjeremylt // Sample runs: 29*6c5df90dSjeremylt // 30*6c5df90dSjeremylt // multigrid -problem bp3 31*6c5df90dSjeremylt // multigrid -problem bp4 -ceed /cpu/self 32*6c5df90dSjeremylt // multigrid -problem bp5 -ceed /cpu/occa 33*6c5df90dSjeremylt // multigrid -problem bp6 -ceed /gpu/cuda 34*6c5df90dSjeremylt // 35*6c5df90dSjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 36*6c5df90dSjeremylt 37*6c5df90dSjeremylt /// @file 38*6c5df90dSjeremylt /// CEED BPs 1-6 multigrid example using PETSc 39*6c5df90dSjeremylt const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n"; 40*6c5df90dSjeremylt 41*6c5df90dSjeremylt #define multigrid 42*6c5df90dSjeremylt #include "setup.h" 43*6c5df90dSjeremylt 44*6c5df90dSjeremylt int main(int argc, char **argv) { 45*6c5df90dSjeremylt PetscInt ierr; 46*6c5df90dSjeremylt MPI_Comm comm; 47*6c5df90dSjeremylt char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self", 48*6c5df90dSjeremylt filename[PETSC_MAX_PATH_LEN]; 49*6c5df90dSjeremylt double my_rt_start, my_rt, rt_min, rt_max; 50*6c5df90dSjeremylt PetscInt degree = 3, qextra, *lsize, *xlsize, *gsize, dim = 3, 51*6c5df90dSjeremylt melem[3] = {3, 3, 3}, ncompu = 1, numlevels = degree, *leveldegrees; 52*6c5df90dSjeremylt PetscScalar *r; 53*6c5df90dSjeremylt PetscBool test_mode, benchmark_mode, read_mesh, write_solution; 54*6c5df90dSjeremylt DM *dm, dmOrig; 55*6c5df90dSjeremylt KSP ksp; 56*6c5df90dSjeremylt PC pc; 57*6c5df90dSjeremylt Mat *matO, *matI, *matR; 58*6c5df90dSjeremylt Vec *X, *Xloc, *mult, rhs, rhsloc, diagloc; 59*6c5df90dSjeremylt UserO *userO; 60*6c5df90dSjeremylt UserIR *userI, *userR; 61*6c5df90dSjeremylt Ceed ceed; 62*6c5df90dSjeremylt CeedData *ceeddata; 63*6c5df90dSjeremylt CeedVector rhsceed, diagceed, target; 64*6c5df90dSjeremylt CeedQFunction qf_error, qf_restrict, qf_prolong; 65*6c5df90dSjeremylt CeedOperator op_error; 66*6c5df90dSjeremylt bpType bpChoice; 67*6c5df90dSjeremylt coarsenType coarsen; 68*6c5df90dSjeremylt 69*6c5df90dSjeremylt ierr = PetscInitialize(&argc, &argv, NULL, help); 70*6c5df90dSjeremylt if (ierr) return ierr; 71*6c5df90dSjeremylt comm = PETSC_COMM_WORLD; 72*6c5df90dSjeremylt 73*6c5df90dSjeremylt // Parse command line options 74*6c5df90dSjeremylt ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 75*6c5df90dSjeremylt bpChoice = CEED_BP3; 76*6c5df90dSjeremylt ierr = PetscOptionsEnum("-problem", 77*6c5df90dSjeremylt "CEED benchmark problem to solve", NULL, 78*6c5df90dSjeremylt bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice, 79*6c5df90dSjeremylt NULL); CHKERRQ(ierr); 80*6c5df90dSjeremylt ncompu = bpOptions[bpChoice].ncompu; 81*6c5df90dSjeremylt test_mode = PETSC_FALSE; 82*6c5df90dSjeremylt ierr = PetscOptionsBool("-test", 83*6c5df90dSjeremylt "Testing mode (do not print unless error is large)", 84*6c5df90dSjeremylt NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 85*6c5df90dSjeremylt benchmark_mode = PETSC_FALSE; 86*6c5df90dSjeremylt ierr = PetscOptionsBool("-benchmark", 87*6c5df90dSjeremylt "Benchmarking mode (prints benchmark statistics)", 88*6c5df90dSjeremylt NULL, benchmark_mode, &benchmark_mode, NULL); 89*6c5df90dSjeremylt CHKERRQ(ierr); 90*6c5df90dSjeremylt write_solution = PETSC_FALSE; 91*6c5df90dSjeremylt ierr = PetscOptionsBool("-write_solution", 92*6c5df90dSjeremylt "Write solution for visualization", 93*6c5df90dSjeremylt NULL, write_solution, &write_solution, NULL); 94*6c5df90dSjeremylt CHKERRQ(ierr); 95*6c5df90dSjeremylt degree = test_mode ? 3 : 2; 96*6c5df90dSjeremylt ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 97*6c5df90dSjeremylt NULL, degree, °ree, NULL); CHKERRQ(ierr); 98*6c5df90dSjeremylt if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 99*6c5df90dSjeremylt "-degree %D must be at least 1", degree); 100*6c5df90dSjeremylt qextra = bpOptions[bpChoice].qextra; 101*6c5df90dSjeremylt ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 102*6c5df90dSjeremylt NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 103*6c5df90dSjeremylt ierr = PetscOptionsString("-ceed", "CEED resource specifier", 104*6c5df90dSjeremylt NULL, ceedresource, ceedresource, 105*6c5df90dSjeremylt sizeof(ceedresource), NULL); CHKERRQ(ierr); 106*6c5df90dSjeremylt coarsen = COARSEN_UNIFORM; 107*6c5df90dSjeremylt ierr = PetscOptionsEnum("-coarsen", 108*6c5df90dSjeremylt "Coarsening strategy to use", NULL, 109*6c5df90dSjeremylt coarsenTypes, (PetscEnum)coarsen, 110*6c5df90dSjeremylt (PetscEnum*)&coarsen, NULL); CHKERRQ(ierr); 111*6c5df90dSjeremylt ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 112*6c5df90dSjeremylt filename, filename, sizeof(filename), &read_mesh); 113*6c5df90dSjeremylt CHKERRQ(ierr); 114*6c5df90dSjeremylt if (!read_mesh) { 115*6c5df90dSjeremylt PetscInt tmp = dim; 116*6c5df90dSjeremylt ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL, 117*6c5df90dSjeremylt melem, &tmp, NULL); CHKERRQ(ierr); 118*6c5df90dSjeremylt } 119*6c5df90dSjeremylt ierr = PetscOptionsEnd(); CHKERRQ(ierr); 120*6c5df90dSjeremylt 121*6c5df90dSjeremylt // Setup DM 122*6c5df90dSjeremylt if (read_mesh) { 123*6c5df90dSjeremylt ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dmOrig); 124*6c5df90dSjeremylt CHKERRQ(ierr); 125*6c5df90dSjeremylt } else { 126*6c5df90dSjeremylt ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL, 127*6c5df90dSjeremylt NULL, NULL, PETSC_TRUE,&dmOrig); CHKERRQ(ierr); 128*6c5df90dSjeremylt } 129*6c5df90dSjeremylt 130*6c5df90dSjeremylt { 131*6c5df90dSjeremylt DM dmDist = NULL; 132*6c5df90dSjeremylt PetscPartitioner part; 133*6c5df90dSjeremylt 134*6c5df90dSjeremylt ierr = DMPlexGetPartitioner(dmOrig, &part); CHKERRQ(ierr); 135*6c5df90dSjeremylt ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 136*6c5df90dSjeremylt ierr = DMPlexDistribute(dmOrig, 0, NULL, &dmDist); CHKERRQ(ierr); 137*6c5df90dSjeremylt if (dmDist) { 138*6c5df90dSjeremylt ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 139*6c5df90dSjeremylt dmOrig = dmDist; 140*6c5df90dSjeremylt } 141*6c5df90dSjeremylt } 142*6c5df90dSjeremylt 143*6c5df90dSjeremylt // Allocate arrays for PETSc objects for each level 144*6c5df90dSjeremylt switch (coarsen) { 145*6c5df90dSjeremylt case COARSEN_UNIFORM: 146*6c5df90dSjeremylt numlevels = degree; 147*6c5df90dSjeremylt break; 148*6c5df90dSjeremylt case COARSEN_LOGRITHMIC: 149*6c5df90dSjeremylt numlevels = ceil(log(degree)/log(2)) + 1; 150*6c5df90dSjeremylt break; 151*6c5df90dSjeremylt } 152*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &leveldegrees); CHKERRQ(ierr); 153*6c5df90dSjeremylt switch (coarsen) { 154*6c5df90dSjeremylt case COARSEN_UNIFORM: 155*6c5df90dSjeremylt for (int i=0; i<numlevels; i++) leveldegrees[i] = i + 1; 156*6c5df90dSjeremylt break; 157*6c5df90dSjeremylt case COARSEN_LOGRITHMIC: 158*6c5df90dSjeremylt for (int i=0; i<numlevels-1; i++) leveldegrees[i] = pow(2,i); 159*6c5df90dSjeremylt leveldegrees[numlevels-1] = degree; 160*6c5df90dSjeremylt break; 161*6c5df90dSjeremylt } 162*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &dm); CHKERRQ(ierr); 163*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &X); CHKERRQ(ierr); 164*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &Xloc); CHKERRQ(ierr); 165*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &mult); CHKERRQ(ierr); 166*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &userO); CHKERRQ(ierr); 167*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &userI); CHKERRQ(ierr); 168*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &userR); CHKERRQ(ierr); 169*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &matO); CHKERRQ(ierr); 170*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &matI); CHKERRQ(ierr); 171*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &matR); CHKERRQ(ierr); 172*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &lsize); CHKERRQ(ierr); 173*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &xlsize); CHKERRQ(ierr); 174*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &gsize); CHKERRQ(ierr); 175*6c5df90dSjeremylt 176*6c5df90dSjeremylt // Setup DM and Operator Mat Shells for each level 177*6c5df90dSjeremylt for (CeedInt i=0; i<numlevels; i++) { 178*6c5df90dSjeremylt // Create DM 179*6c5df90dSjeremylt ierr = DMClone(dmOrig, &dm[i]); CHKERRQ(ierr); 180*6c5df90dSjeremylt ierr = SetupDMByDegree(dm[i], leveldegrees[i], ncompu, bpChoice); 181*6c5df90dSjeremylt CHKERRQ(ierr); 182*6c5df90dSjeremylt 183*6c5df90dSjeremylt // Create vectors 184*6c5df90dSjeremylt ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr); 185*6c5df90dSjeremylt ierr = VecGetLocalSize(X[i], &lsize[i]); CHKERRQ(ierr); 186*6c5df90dSjeremylt ierr = VecGetSize(X[i], &gsize[i]); CHKERRQ(ierr); 187*6c5df90dSjeremylt ierr = DMCreateLocalVector(dm[i], &Xloc[i]); CHKERRQ(ierr); 188*6c5df90dSjeremylt ierr = VecGetSize(Xloc[i], &xlsize[i]); CHKERRQ(ierr); 189*6c5df90dSjeremylt 190*6c5df90dSjeremylt // Operator 191*6c5df90dSjeremylt ierr = PetscMalloc1(1, &userO[i]); CHKERRQ(ierr); 192*6c5df90dSjeremylt ierr = MatCreateShell(comm, lsize[i], lsize[i], gsize[i], gsize[i], 193*6c5df90dSjeremylt userO[i], &matO[i]); CHKERRQ(ierr); 194*6c5df90dSjeremylt ierr = MatShellSetOperation(matO[i], MATOP_MULT, 195*6c5df90dSjeremylt (void(*)(void))MatMult_Ceed); 196*6c5df90dSjeremylt ierr = MatShellSetOperation(matO[i], MATOP_GET_DIAGONAL, 197*6c5df90dSjeremylt (void(*)(void))MatGetDiag); 198*6c5df90dSjeremylt CHKERRQ(ierr); 199*6c5df90dSjeremylt 200*6c5df90dSjeremylt // Level transfers 201*6c5df90dSjeremylt if (i > 0) { 202*6c5df90dSjeremylt // Interp 203*6c5df90dSjeremylt ierr = PetscMalloc1(1, &userI[i]); CHKERRQ(ierr); 204*6c5df90dSjeremylt ierr = MatCreateShell(comm, lsize[i], lsize[i-1], gsize[i], gsize[i-1], 205*6c5df90dSjeremylt userI[i], &matI[i]); CHKERRQ(ierr); 206*6c5df90dSjeremylt ierr = MatShellSetOperation(matI[i], MATOP_MULT, 207*6c5df90dSjeremylt (void(*)(void))MatMult_Interp); 208*6c5df90dSjeremylt CHKERRQ(ierr); 209*6c5df90dSjeremylt 210*6c5df90dSjeremylt // Restrict 211*6c5df90dSjeremylt ierr = PetscMalloc1(1, &userR[i]); CHKERRQ(ierr); 212*6c5df90dSjeremylt ierr = MatCreateShell(comm, lsize[i-1], lsize[i], gsize[i-1], gsize[i], 213*6c5df90dSjeremylt userR[i], &matR[i]); CHKERRQ(ierr); 214*6c5df90dSjeremylt ierr = MatShellSetOperation(matR[i], MATOP_MULT, 215*6c5df90dSjeremylt (void(*)(void))MatMult_Restrict); 216*6c5df90dSjeremylt CHKERRQ(ierr); 217*6c5df90dSjeremylt } 218*6c5df90dSjeremylt } 219*6c5df90dSjeremylt ierr = VecDuplicate(X[numlevels-1], &rhs); CHKERRQ(ierr); 220*6c5df90dSjeremylt 221*6c5df90dSjeremylt // Print global grid information 222*6c5df90dSjeremylt if (!test_mode) { 223*6c5df90dSjeremylt PetscInt P = degree + 1, Q = P + qextra; 224*6c5df90dSjeremylt ierr = PetscPrintf(comm, 225*6c5df90dSjeremylt "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n" 226*6c5df90dSjeremylt " libCEED:\n" 227*6c5df90dSjeremylt " libCEED Backend : %s\n" 228*6c5df90dSjeremylt " Mesh:\n" 229*6c5df90dSjeremylt " Number of 1D Basis Nodes (p) : %d\n" 230*6c5df90dSjeremylt " Number of 1D Quadrature Points (q) : %d\n" 231*6c5df90dSjeremylt " Global Nodes : %D\n" 232*6c5df90dSjeremylt " Owned Nodes : %D\n" 233*6c5df90dSjeremylt " Multigrid:\n" 234*6c5df90dSjeremylt " Number of Levels : %d\n", 235*6c5df90dSjeremylt bpChoice+1, ceedresource, P, Q, 236*6c5df90dSjeremylt gsize[numlevels-1]/ncompu, lsize[numlevels-1]/ncompu, 237*6c5df90dSjeremylt numlevels); CHKERRQ(ierr); 238*6c5df90dSjeremylt } 239*6c5df90dSjeremylt 240*6c5df90dSjeremylt // Set up libCEED 241*6c5df90dSjeremylt CeedInit(ceedresource, &ceed); 242*6c5df90dSjeremylt 243*6c5df90dSjeremylt // Create RHS vector 244*6c5df90dSjeremylt ierr = VecDuplicate(Xloc[numlevels-1], &rhsloc); CHKERRQ(ierr); 245*6c5df90dSjeremylt ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 246*6c5df90dSjeremylt ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 247*6c5df90dSjeremylt CeedVectorCreate(ceed, xlsize[numlevels-1], &rhsceed); 248*6c5df90dSjeremylt CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r); 249*6c5df90dSjeremylt 250*6c5df90dSjeremylt // Set up libCEED operators on each level 251*6c5df90dSjeremylt ierr = PetscMalloc1(numlevels, &ceeddata); CHKERRQ(ierr); 252*6c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 253*6c5df90dSjeremylt // Print level information 254*6c5df90dSjeremylt if (!test_mode && (i == 0 || i == numlevels-1)) { 255*6c5df90dSjeremylt ierr = PetscPrintf(comm," Level %D (%s):\n" 256*6c5df90dSjeremylt " Number of 1D Basis Nodes (p) : %d\n" 257*6c5df90dSjeremylt " Global Nodes : %D\n" 258*6c5df90dSjeremylt " Owned Nodes : %D\n", 259*6c5df90dSjeremylt i, (i? "fine" : "coarse"), leveldegrees[i] + 1, 260*6c5df90dSjeremylt gsize[i]/ncompu, lsize[i]/ncompu); CHKERRQ(ierr); 261*6c5df90dSjeremylt } 262*6c5df90dSjeremylt ierr = PetscMalloc1(1, &ceeddata[i]); CHKERRQ(ierr); 263*6c5df90dSjeremylt ierr = SetupLibceedByDegree(dm[i], ceed, leveldegrees[i], dim, qextra, 264*6c5df90dSjeremylt ncompu, gsize[i], xlsize[i], bpChoice, ceeddata[i], 265*6c5df90dSjeremylt i==(numlevels-1), rhsceed, &target); CHKERRQ(ierr); 266*6c5df90dSjeremylt } 267*6c5df90dSjeremylt 268*6c5df90dSjeremylt // Gather RHS 269*6c5df90dSjeremylt ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 270*6c5df90dSjeremylt ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 271*6c5df90dSjeremylt ierr = DMLocalToGlobalBegin(dm[numlevels-1], rhsloc, ADD_VALUES, rhs); 272*6c5df90dSjeremylt CHKERRQ(ierr); 273*6c5df90dSjeremylt ierr = DMLocalToGlobalEnd(dm[numlevels-1], rhsloc, ADD_VALUES, rhs); 274*6c5df90dSjeremylt CHKERRQ(ierr); 275*6c5df90dSjeremylt CeedVectorDestroy(&rhsceed); 276*6c5df90dSjeremylt 277*6c5df90dSjeremylt // Create the restriction/interpolation Q-function 278*6c5df90dSjeremylt CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident, 279*6c5df90dSjeremylt bpOptions[bpChoice].identfname, &qf_restrict); 280*6c5df90dSjeremylt CeedQFunctionAddInput(qf_restrict, "uin", ncompu, CEED_EVAL_NONE); 281*6c5df90dSjeremylt CeedQFunctionAddOutput(qf_restrict, "uout", ncompu, CEED_EVAL_INTERP); 282*6c5df90dSjeremylt CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].ident, 283*6c5df90dSjeremylt bpOptions[bpChoice].identfname, &qf_prolong); 284*6c5df90dSjeremylt CeedQFunctionAddInput(qf_prolong, "uin", ncompu, CEED_EVAL_INTERP); 285*6c5df90dSjeremylt CeedQFunctionAddOutput(qf_prolong, "uout", ncompu, CEED_EVAL_NONE); 286*6c5df90dSjeremylt 287*6c5df90dSjeremylt // Set up libCEED level transfer operators 288*6c5df90dSjeremylt ierr = CeedLevelTransferSetup(ceed, numlevels, ncompu, bpChoice, ceeddata, 289*6c5df90dSjeremylt leveldegrees, qf_restrict, qf_prolong); 290*6c5df90dSjeremylt CHKERRQ(ierr); 291*6c5df90dSjeremylt 292*6c5df90dSjeremylt // Create the error Q-function 293*6c5df90dSjeremylt CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error, 294*6c5df90dSjeremylt bpOptions[bpChoice].errorfname, &qf_error); 295*6c5df90dSjeremylt CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP); 296*6c5df90dSjeremylt CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE); 297*6c5df90dSjeremylt CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE); 298*6c5df90dSjeremylt 299*6c5df90dSjeremylt // Create the error operator 300*6c5df90dSjeremylt CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 301*6c5df90dSjeremylt CeedOperatorSetField(op_error, "u", ceeddata[numlevels-1]->Erestrictu, 302*6c5df90dSjeremylt CEED_TRANSPOSE, ceeddata[numlevels-1]->basisu, 303*6c5df90dSjeremylt CEED_VECTOR_ACTIVE); 304*6c5df90dSjeremylt CeedOperatorSetField(op_error, "true_soln", ceeddata[numlevels-1]->Erestrictui, 305*6c5df90dSjeremylt CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, target); 306*6c5df90dSjeremylt CeedOperatorSetField(op_error, "error", ceeddata[numlevels-1]->Erestrictui, 307*6c5df90dSjeremylt CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, 308*6c5df90dSjeremylt CEED_VECTOR_ACTIVE); 309*6c5df90dSjeremylt 310*6c5df90dSjeremylt // Calculate multiplicity 311*6c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 312*6c5df90dSjeremylt PetscScalar *x; 313*6c5df90dSjeremylt 314*6c5df90dSjeremylt // CEED vector 315*6c5df90dSjeremylt ierr = VecGetArray(Xloc[i], &x); CHKERRQ(ierr); 316*6c5df90dSjeremylt CeedVectorSetArray(ceeddata[i]->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 317*6c5df90dSjeremylt 318*6c5df90dSjeremylt // Multiplicity 319*6c5df90dSjeremylt CeedElemRestrictionGetMultiplicity(ceeddata[i]->Erestrictu, 320*6c5df90dSjeremylt ceeddata[i]->xceed); 321*6c5df90dSjeremylt 322*6c5df90dSjeremylt // Restore vector 323*6c5df90dSjeremylt ierr = VecRestoreArray(Xloc[i], &x); CHKERRQ(ierr); 324*6c5df90dSjeremylt 325*6c5df90dSjeremylt // Creat mult vector 326*6c5df90dSjeremylt ierr = VecDuplicate(Xloc[i], &mult[i]); CHKERRQ(ierr); 327*6c5df90dSjeremylt 328*6c5df90dSjeremylt // Local-to-global 329*6c5df90dSjeremylt ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 330*6c5df90dSjeremylt ierr = DMLocalToGlobalBegin(dm[i], Xloc[i], ADD_VALUES, X[i]); 331*6c5df90dSjeremylt CHKERRQ(ierr); 332*6c5df90dSjeremylt ierr = DMLocalToGlobalEnd(dm[i], Xloc[i], ADD_VALUES, X[i]); 333*6c5df90dSjeremylt CHKERRQ(ierr); 334*6c5df90dSjeremylt ierr = VecZeroEntries(Xloc[i]); CHKERRQ(ierr); 335*6c5df90dSjeremylt 336*6c5df90dSjeremylt // Global-to-local 337*6c5df90dSjeremylt ierr = DMGlobalToLocalBegin(dm[i], X[i], INSERT_VALUES, mult[i]); 338*6c5df90dSjeremylt CHKERRQ(ierr); 339*6c5df90dSjeremylt ierr = DMGlobalToLocalEnd(dm[i], X[i], INSERT_VALUES, mult[i]); 340*6c5df90dSjeremylt CHKERRQ(ierr); 341*6c5df90dSjeremylt ierr = VecZeroEntries(X[i]); CHKERRQ(ierr); 342*6c5df90dSjeremylt 343*6c5df90dSjeremylt // Multiplicity scaling 344*6c5df90dSjeremylt ierr = VecReciprocal(mult[i]); 345*6c5df90dSjeremylt } 346*6c5df90dSjeremylt 347*6c5df90dSjeremylt // Set up Mat 348*6c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 349*6c5df90dSjeremylt // Diff Operator 350*6c5df90dSjeremylt userO[i]->comm = comm; 351*6c5df90dSjeremylt userO[i]->dm = dm[i]; 352*6c5df90dSjeremylt userO[i]->Xloc = Xloc[i]; 353*6c5df90dSjeremylt ierr = VecDuplicate(Xloc[i], &userO[i]->Yloc); CHKERRQ(ierr); 354*6c5df90dSjeremylt userO[i]->xceed = ceeddata[i]->xceed; 355*6c5df90dSjeremylt userO[i]->yceed = ceeddata[i]->yceed; 356*6c5df90dSjeremylt userO[i]->op = ceeddata[i]->op_apply; 357*6c5df90dSjeremylt userO[i]->ceed = ceed; 358*6c5df90dSjeremylt 359*6c5df90dSjeremylt // Set up diagonal 360*6c5df90dSjeremylt const CeedScalar *ceedarray; 361*6c5df90dSjeremylt PetscScalar *petscarray; 362*6c5df90dSjeremylt CeedInt length; 363*6c5df90dSjeremylt 364*6c5df90dSjeremylt ierr = VecDuplicate(X[i], &userO[i]->diag); CHKERRQ(ierr); 365*6c5df90dSjeremylt ierr = VecDuplicate(Xloc[i], &diagloc); CHKERRQ(ierr); 366*6c5df90dSjeremylt 367*6c5df90dSjeremylt // -- Local diagonal 368*6c5df90dSjeremylt CeedOperatorAssembleLinearDiagonal(userO[i]->op, &diagceed, 369*6c5df90dSjeremylt CEED_REQUEST_IMMEDIATE); 370*6c5df90dSjeremylt 371*6c5df90dSjeremylt // -- Copy values 372*6c5df90dSjeremylt CeedVectorGetArrayRead(diagceed, CEED_MEM_HOST, &ceedarray); 373*6c5df90dSjeremylt ierr = VecGetArray(diagloc, &petscarray); CHKERRQ(ierr); 374*6c5df90dSjeremylt CeedVectorGetLength(diagceed, &length); 375*6c5df90dSjeremylt for (CeedInt i=0; i<length; i++) 376*6c5df90dSjeremylt petscarray[i] = ceedarray[i]; 377*6c5df90dSjeremylt CeedVectorRestoreArrayRead(diagceed, &ceedarray); 378*6c5df90dSjeremylt ierr = VecRestoreArray(diagloc, &petscarray); CHKERRQ(ierr); 379*6c5df90dSjeremylt 380*6c5df90dSjeremylt // -- Global diagonal 381*6c5df90dSjeremylt ierr = VecZeroEntries(userO[i]->diag); CHKERRQ(ierr); 382*6c5df90dSjeremylt ierr = DMLocalToGlobalBegin(userO[i]->dm, diagloc, ADD_VALUES, 383*6c5df90dSjeremylt userO[i]->diag); CHKERRQ(ierr); 384*6c5df90dSjeremylt ierr = DMLocalToGlobalEnd(userO[i]->dm, diagloc, ADD_VALUES, 385*6c5df90dSjeremylt userO[i]->diag); CHKERRQ(ierr); 386*6c5df90dSjeremylt 387*6c5df90dSjeremylt // -- Cleanup 388*6c5df90dSjeremylt ierr = VecDestroy(&diagloc); CHKERRQ(ierr); 389*6c5df90dSjeremylt CeedVectorDestroy(&diagceed); 390*6c5df90dSjeremylt 391*6c5df90dSjeremylt if (i > 0) { 392*6c5df90dSjeremylt // Interp Operator 393*6c5df90dSjeremylt userI[i]->comm = comm; 394*6c5df90dSjeremylt userI[i]->dmc = dm[i-1]; 395*6c5df90dSjeremylt userI[i]->dmf = dm[i]; 396*6c5df90dSjeremylt userI[i]->Xloc = Xloc[i-1]; 397*6c5df90dSjeremylt userI[i]->Yloc = userO[i]->Yloc; 398*6c5df90dSjeremylt userI[i]->mult = mult[i]; 399*6c5df90dSjeremylt userI[i]->ceedvecc = userO[i-1]->xceed; 400*6c5df90dSjeremylt userI[i]->ceedvecf = userO[i]->yceed; 401*6c5df90dSjeremylt userI[i]->op = ceeddata[i]->op_interp; 402*6c5df90dSjeremylt userI[i]->ceed = ceed; 403*6c5df90dSjeremylt 404*6c5df90dSjeremylt // Restrict Operator 405*6c5df90dSjeremylt userR[i]->comm = comm; 406*6c5df90dSjeremylt userR[i]->dmc = dm[i-1]; 407*6c5df90dSjeremylt userR[i]->dmf = dm[i]; 408*6c5df90dSjeremylt userR[i]->Xloc = Xloc[i]; 409*6c5df90dSjeremylt userR[i]->Yloc = userO[i-1]->Yloc; 410*6c5df90dSjeremylt userR[i]->mult = mult[i]; 411*6c5df90dSjeremylt userR[i]->ceedvecf = userO[i]->xceed; 412*6c5df90dSjeremylt userR[i]->ceedvecc = userO[i-1]->yceed; 413*6c5df90dSjeremylt userR[i]->op = ceeddata[i]->op_restrict; 414*6c5df90dSjeremylt userR[i]->ceed = ceed; 415*6c5df90dSjeremylt } 416*6c5df90dSjeremylt } 417*6c5df90dSjeremylt 418*6c5df90dSjeremylt // Set up KSP 419*6c5df90dSjeremylt ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 420*6c5df90dSjeremylt { 421*6c5df90dSjeremylt ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 422*6c5df90dSjeremylt ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 423*6c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 424*6c5df90dSjeremylt PETSC_DEFAULT); CHKERRQ(ierr); 425*6c5df90dSjeremylt } 426*6c5df90dSjeremylt ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 427*6c5df90dSjeremylt ierr = KSPSetOperators(ksp, matO[numlevels-1], matO[numlevels-1]); 428*6c5df90dSjeremylt CHKERRQ(ierr); 429*6c5df90dSjeremylt 430*6c5df90dSjeremylt // Set up PCMG 431*6c5df90dSjeremylt ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 432*6c5df90dSjeremylt PCMGCycleType pcgmcycletype = PC_MG_CYCLE_V; 433*6c5df90dSjeremylt { 434*6c5df90dSjeremylt ierr = PCSetType(pc, PCMG); CHKERRQ(ierr); 435*6c5df90dSjeremylt 436*6c5df90dSjeremylt // PCMG levels 437*6c5df90dSjeremylt ierr = PCMGSetLevels(pc, numlevels, NULL); CHKERRQ(ierr); 438*6c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 439*6c5df90dSjeremylt // Smoother 440*6c5df90dSjeremylt KSP smoother; 441*6c5df90dSjeremylt PC smoother_pc; 442*6c5df90dSjeremylt ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr); 443*6c5df90dSjeremylt ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr); 444*6c5df90dSjeremylt ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr); 445*6c5df90dSjeremylt ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr); 446*6c5df90dSjeremylt ierr = KSPSetOperators(smoother, matO[i], matO[i]); CHKERRQ(ierr); 447*6c5df90dSjeremylt ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr); 448*6c5df90dSjeremylt ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr); 449*6c5df90dSjeremylt ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 450*6c5df90dSjeremylt 451*6c5df90dSjeremylt // Work vector 452*6c5df90dSjeremylt if (i < numlevels-1) { 453*6c5df90dSjeremylt ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr); 454*6c5df90dSjeremylt } 455*6c5df90dSjeremylt 456*6c5df90dSjeremylt // Level transfers 457*6c5df90dSjeremylt if (i > 0) { 458*6c5df90dSjeremylt // Interpolation 459*6c5df90dSjeremylt ierr = PCMGSetInterpolation(pc, i, matI[i]); CHKERRQ(ierr); 460*6c5df90dSjeremylt 461*6c5df90dSjeremylt // Restriction 462*6c5df90dSjeremylt ierr = PCMGSetRestriction(pc, i, matR[i]); CHKERRQ(ierr); 463*6c5df90dSjeremylt } 464*6c5df90dSjeremylt 465*6c5df90dSjeremylt // Coarse solve 466*6c5df90dSjeremylt KSP coarse; 467*6c5df90dSjeremylt PC coarse_pc; 468*6c5df90dSjeremylt ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr); 469*6c5df90dSjeremylt ierr = KSPSetType(coarse, KSPCG); CHKERRQ(ierr); 470*6c5df90dSjeremylt ierr = KSPSetOperators(coarse, matO[0], matO[0]); CHKERRQ(ierr); 471*6c5df90dSjeremylt ierr = KSPSetTolerances(coarse, 1e-10, 1e-10, PETSC_DEFAULT, 472*6c5df90dSjeremylt PETSC_DEFAULT); CHKERRQ(ierr); 473*6c5df90dSjeremylt ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr); 474*6c5df90dSjeremylt ierr = PCSetType(coarse_pc, PCJACOBI); CHKERRQ(ierr); 475*6c5df90dSjeremylt ierr = PCJacobiSetType(coarse_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr); 476*6c5df90dSjeremylt } 477*6c5df90dSjeremylt 478*6c5df90dSjeremylt // PCMG options 479*6c5df90dSjeremylt ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr); 480*6c5df90dSjeremylt ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr); 481*6c5df90dSjeremylt ierr = PCMGSetCycleType(pc, pcgmcycletype); CHKERRQ(ierr); 482*6c5df90dSjeremylt } 483*6c5df90dSjeremylt 484*6c5df90dSjeremylt // First run, if benchmarking 485*6c5df90dSjeremylt if (benchmark_mode) { 486*6c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 487*6c5df90dSjeremylt CHKERRQ(ierr); 488*6c5df90dSjeremylt ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr); 489*6c5df90dSjeremylt my_rt_start = MPI_Wtime(); 490*6c5df90dSjeremylt ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr); 491*6c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 492*6c5df90dSjeremylt ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 493*6c5df90dSjeremylt CHKERRQ(ierr); 494*6c5df90dSjeremylt // Set maxits based on first iteration timing 495*6c5df90dSjeremylt if (my_rt > 0.02) { 496*6c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 497*6c5df90dSjeremylt CHKERRQ(ierr); 498*6c5df90dSjeremylt } else { 499*6c5df90dSjeremylt ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 500*6c5df90dSjeremylt CHKERRQ(ierr); 501*6c5df90dSjeremylt } 502*6c5df90dSjeremylt } 503*6c5df90dSjeremylt 504*6c5df90dSjeremylt // Timed solve 505*6c5df90dSjeremylt ierr = VecZeroEntries(X[numlevels-1]); CHKERRQ(ierr); 506*6c5df90dSjeremylt ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 507*6c5df90dSjeremylt my_rt_start = MPI_Wtime(); 508*6c5df90dSjeremylt ierr = KSPSolve(ksp, rhs, X[numlevels-1]); CHKERRQ(ierr); 509*6c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 510*6c5df90dSjeremylt 511*6c5df90dSjeremylt // Output results 512*6c5df90dSjeremylt { 513*6c5df90dSjeremylt KSPType ksptype; 514*6c5df90dSjeremylt PCMGType pcmgtype; 515*6c5df90dSjeremylt KSPConvergedReason reason; 516*6c5df90dSjeremylt PetscReal rnorm; 517*6c5df90dSjeremylt PetscInt its; 518*6c5df90dSjeremylt ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 519*6c5df90dSjeremylt ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 520*6c5df90dSjeremylt ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 521*6c5df90dSjeremylt ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 522*6c5df90dSjeremylt ierr = PCMGGetType(pc, &pcmgtype); CHKERRQ(ierr); 523*6c5df90dSjeremylt if (!test_mode || reason < 0 || rnorm > 1e-8) { 524*6c5df90dSjeremylt ierr = PetscPrintf(comm, 525*6c5df90dSjeremylt " KSP:\n" 526*6c5df90dSjeremylt " KSP Type : %s\n" 527*6c5df90dSjeremylt " KSP Convergence : %s\n" 528*6c5df90dSjeremylt " Total KSP Iterations : %D\n" 529*6c5df90dSjeremylt " Final rnorm : %e\n", 530*6c5df90dSjeremylt ksptype, KSPConvergedReasons[reason], its, 531*6c5df90dSjeremylt (double)rnorm); CHKERRQ(ierr); 532*6c5df90dSjeremylt ierr = PetscPrintf(comm, 533*6c5df90dSjeremylt " PCMG:\n" 534*6c5df90dSjeremylt " PCMG Type : %s\n" 535*6c5df90dSjeremylt " PCMG Cycle Type : %s\n", 536*6c5df90dSjeremylt PCMGTypes[pcmgtype], 537*6c5df90dSjeremylt PCMGCycleTypes[pcgmcycletype]); CHKERRQ(ierr); 538*6c5df90dSjeremylt } 539*6c5df90dSjeremylt if (!test_mode) { 540*6c5df90dSjeremylt ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 541*6c5df90dSjeremylt } 542*6c5df90dSjeremylt { 543*6c5df90dSjeremylt PetscReal maxerror; 544*6c5df90dSjeremylt ierr = ComputeErrorMax(userO[numlevels-1], op_error, X[numlevels-1], target, 545*6c5df90dSjeremylt &maxerror); CHKERRQ(ierr); 546*6c5df90dSjeremylt PetscReal tol = 5e-2; 547*6c5df90dSjeremylt if (!test_mode || maxerror > tol) { 548*6c5df90dSjeremylt ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 549*6c5df90dSjeremylt CHKERRQ(ierr); 550*6c5df90dSjeremylt ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 551*6c5df90dSjeremylt CHKERRQ(ierr); 552*6c5df90dSjeremylt ierr = PetscPrintf(comm, 553*6c5df90dSjeremylt " Pointwise Error (max) : %e\n" 554*6c5df90dSjeremylt " CG Solve Time : %g (%g) sec\n", 555*6c5df90dSjeremylt (double)maxerror, rt_max, rt_min); CHKERRQ(ierr); 556*6c5df90dSjeremylt } 557*6c5df90dSjeremylt } 558*6c5df90dSjeremylt if (benchmark_mode && (!test_mode)) { 559*6c5df90dSjeremylt ierr = PetscPrintf(comm, 560*6c5df90dSjeremylt " DoFs/Sec in CG : %g (%g) million\n", 561*6c5df90dSjeremylt 1e-6*gsize[numlevels-1]*its/rt_max, 562*6c5df90dSjeremylt 1e-6*gsize[numlevels-1]*its/rt_min); 563*6c5df90dSjeremylt CHKERRQ(ierr); 564*6c5df90dSjeremylt } 565*6c5df90dSjeremylt } 566*6c5df90dSjeremylt 567*6c5df90dSjeremylt if (write_solution) { 568*6c5df90dSjeremylt PetscViewer vtkviewersoln; 569*6c5df90dSjeremylt 570*6c5df90dSjeremylt ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 571*6c5df90dSjeremylt ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 572*6c5df90dSjeremylt ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 573*6c5df90dSjeremylt ierr = VecView(X[numlevels-1], vtkviewersoln); CHKERRQ(ierr); 574*6c5df90dSjeremylt ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 575*6c5df90dSjeremylt } 576*6c5df90dSjeremylt 577*6c5df90dSjeremylt // Cleanup 578*6c5df90dSjeremylt for (int i=0; i<numlevels; i++) { 579*6c5df90dSjeremylt ierr = VecDestroy(&X[i]); CHKERRQ(ierr); 580*6c5df90dSjeremylt ierr = VecDestroy(&Xloc[i]); CHKERRQ(ierr); 581*6c5df90dSjeremylt ierr = VecDestroy(&mult[i]); CHKERRQ(ierr); 582*6c5df90dSjeremylt ierr = VecDestroy(&userO[i]->Yloc); CHKERRQ(ierr); 583*6c5df90dSjeremylt ierr = VecDestroy(&userO[i]->diag); CHKERRQ(ierr); 584*6c5df90dSjeremylt ierr = MatDestroy(&matO[i]); CHKERRQ(ierr); 585*6c5df90dSjeremylt ierr = PetscFree(userO[i]); CHKERRQ(ierr); 586*6c5df90dSjeremylt if (i > 0) { 587*6c5df90dSjeremylt ierr = MatDestroy(&matI[i]); CHKERRQ(ierr); 588*6c5df90dSjeremylt ierr = PetscFree(userI[i]); CHKERRQ(ierr); 589*6c5df90dSjeremylt ierr = MatDestroy(&matR[i]); CHKERRQ(ierr); 590*6c5df90dSjeremylt ierr = PetscFree(userR[i]); CHKERRQ(ierr); 591*6c5df90dSjeremylt } 592*6c5df90dSjeremylt ierr = CeedDataDestroy(i, ceeddata[i]); CHKERRQ(ierr); 593*6c5df90dSjeremylt ierr = DMDestroy(&dm[i]); CHKERRQ(ierr); 594*6c5df90dSjeremylt } 595*6c5df90dSjeremylt ierr = PetscFree(leveldegrees); CHKERRQ(ierr); 596*6c5df90dSjeremylt ierr = PetscFree(dm); CHKERRQ(ierr); 597*6c5df90dSjeremylt ierr = PetscFree(X); CHKERRQ(ierr); 598*6c5df90dSjeremylt ierr = PetscFree(Xloc); CHKERRQ(ierr); 599*6c5df90dSjeremylt ierr = PetscFree(mult); CHKERRQ(ierr); 600*6c5df90dSjeremylt ierr = PetscFree(matO); CHKERRQ(ierr); 601*6c5df90dSjeremylt ierr = PetscFree(matI); CHKERRQ(ierr); 602*6c5df90dSjeremylt ierr = PetscFree(matR); CHKERRQ(ierr); 603*6c5df90dSjeremylt ierr = PetscFree(ceeddata); CHKERRQ(ierr); 604*6c5df90dSjeremylt ierr = PetscFree(userO); CHKERRQ(ierr); 605*6c5df90dSjeremylt ierr = PetscFree(userI); CHKERRQ(ierr); 606*6c5df90dSjeremylt ierr = PetscFree(userR); CHKERRQ(ierr); 607*6c5df90dSjeremylt ierr = PetscFree(lsize); CHKERRQ(ierr); 608*6c5df90dSjeremylt ierr = PetscFree(xlsize); CHKERRQ(ierr); 609*6c5df90dSjeremylt ierr = PetscFree(gsize); CHKERRQ(ierr); 610*6c5df90dSjeremylt ierr = VecDestroy(&rhs); CHKERRQ(ierr); 611*6c5df90dSjeremylt ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 612*6c5df90dSjeremylt ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 613*6c5df90dSjeremylt ierr = DMDestroy(&dmOrig); CHKERRQ(ierr); 614*6c5df90dSjeremylt CeedVectorDestroy(&target); 615*6c5df90dSjeremylt CeedQFunctionDestroy(&qf_error); 616*6c5df90dSjeremylt CeedQFunctionDestroy(&qf_restrict); 617*6c5df90dSjeremylt CeedQFunctionDestroy(&qf_prolong); 618*6c5df90dSjeremylt CeedOperatorDestroy(&op_error); 619*6c5df90dSjeremylt CeedDestroy(&ceed); 620*6c5df90dSjeremylt return PetscFinalize(); 621*6c5df90dSjeremylt } 622