1819eb1b3Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2819eb1b3Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3819eb1b3Sjeremylt // reserved. See files LICENSE and NOTICE for details. 4819eb1b3Sjeremylt // 5819eb1b3Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software 6819eb1b3Sjeremylt // libraries and APIs for efficient high-order finite element and spectral 7819eb1b3Sjeremylt // element discretizations for exascale applications. For more information and 8819eb1b3Sjeremylt // source code availability see http://github.com/ceed. 9819eb1b3Sjeremylt // 10819eb1b3Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11819eb1b3Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office 12819eb1b3Sjeremylt // of Science and the National Nuclear Security Administration) responsible for 13819eb1b3Sjeremylt // the planning and preparation of a capable exascale ecosystem, including 14819eb1b3Sjeremylt // software, applications, hardware, advanced system engineering and early 15819eb1b3Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative. 16819eb1b3Sjeremylt 170c59ef15SJeremy L Thompson // libCEED + PETSc Example: CEED BPs 180c59ef15SJeremy L Thompson // 190c59ef15SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the 200c59ef15SJeremy L Thompson // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps. 210c59ef15SJeremy L Thompson // 220c59ef15SJeremy L Thompson // The code is intentionally "raw", using only low-level communication 230c59ef15SJeremy L Thompson // primitives. 240c59ef15SJeremy L Thompson // 250c59ef15SJeremy L Thompson // Build with: 260c59ef15SJeremy L Thompson // 270c59ef15SJeremy L Thompson // make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 280c59ef15SJeremy L Thompson // 290c59ef15SJeremy L Thompson // Sample runs: 300c59ef15SJeremy L Thompson // 310c59ef15SJeremy L Thompson // bps -problem bp1 320c59ef15SJeremy L Thompson // bps -problem bp2 -ceed /cpu/self 330c59ef15SJeremy L Thompson // bps -problem bp3 -ceed /gpu/occa 340c59ef15SJeremy L Thompson // bps -problem bp4 -ceed /cpu/occa 350c59ef15SJeremy L Thompson // bps -problem bp5 -ceed /omp/occa 360c59ef15SJeremy L Thompson // bps -problem bp6 -ceed /ocl/occa 370c59ef15SJeremy L Thompson // 386d3db1e3Sjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp2 390c59ef15SJeremy L Thompson 400c59ef15SJeremy L Thompson /// @file 410c59ef15SJeremy L Thompson /// CEED BPs example using PETSc 42819eb1b3Sjeremylt /// See bpsdmplex.c for an implementation using DMPlex unstructured grids. 430c59ef15SJeremy L Thompson const char help[] = "Solve CEED BPs using PETSc\n"; 440c59ef15SJeremy L Thompson 450c59ef15SJeremy L Thompson #include <stdbool.h> 460c59ef15SJeremy L Thompson #include <string.h> 474d537eeaSYohann #include <petscksp.h> 484d537eeaSYohann #include <ceed.h> 490c59ef15SJeremy L Thompson #include "common.h" 500c59ef15SJeremy L Thompson #include "bp1.h" 510c59ef15SJeremy L Thompson #include "bp2.h" 520c59ef15SJeremy L Thompson #include "bp3.h" 530c59ef15SJeremy L Thompson #include "bp4.h" 540c59ef15SJeremy L Thompson 550c59ef15SJeremy L Thompson static void Split3(PetscInt size, PetscInt m[3], bool reverse) { 560c59ef15SJeremy L Thompson for (PetscInt d=0,sizeleft=size; d<3; d++) { 570c59ef15SJeremy L Thompson PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d))); 580c59ef15SJeremy L Thompson while (try * (sizeleft / try) != sizeleft) try++; 590c59ef15SJeremy L Thompson m[reverse ? 2-d : d] = try; 600c59ef15SJeremy L Thompson sizeleft /= try; 610c59ef15SJeremy L Thompson } 620c59ef15SJeremy L Thompson } 630c59ef15SJeremy L Thompson 640c59ef15SJeremy L Thompson static PetscInt Max3(const PetscInt a[3]) { 650c59ef15SJeremy L Thompson return PetscMax(a[0], PetscMax(a[1], a[2])); 660c59ef15SJeremy L Thompson } 670c59ef15SJeremy L Thompson static PetscInt Min3(const PetscInt a[3]) { 680c59ef15SJeremy L Thompson return PetscMin(a[0], PetscMin(a[1], a[2])); 690c59ef15SJeremy L Thompson } 70b52ddc88Sjeremylt static void GlobalNodes(const PetscInt p[3], const PetscInt irank[3], 710c59ef15SJeremy L Thompson PetscInt degree, const PetscInt melem[3], 72b52ddc88Sjeremylt PetscInt mnodes[3]) { 730c59ef15SJeremy L Thompson for (int d=0; d<3; d++) 74b52ddc88Sjeremylt mnodes[d] = degree*melem[d] + (irank[d] == p[d]-1); 750c59ef15SJeremy L Thompson } 760c59ef15SJeremy L Thompson static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3], 770c59ef15SJeremy L Thompson PetscInt degree, const PetscInt melem[3]) { 780c59ef15SJeremy L Thompson PetscInt start = 0; 790c59ef15SJeremy L Thompson // Dumb brute-force is easier to read 800c59ef15SJeremy L Thompson for (PetscInt i=0; i<p[0]; i++) { 810c59ef15SJeremy L Thompson for (PetscInt j=0; j<p[1]; j++) { 820c59ef15SJeremy L Thompson for (PetscInt k=0; k<p[2]; k++) { 83b52ddc88Sjeremylt PetscInt mnodes[3], ijkrank[] = {i,j,k}; 840c59ef15SJeremy L Thompson if (i == irank[0] && j == irank[1] && k == irank[2]) return start; 85b52ddc88Sjeremylt GlobalNodes(p, ijkrank, degree, melem, mnodes); 86b52ddc88Sjeremylt start += mnodes[0] * mnodes[1] * mnodes[2]; 870c59ef15SJeremy L Thompson } 880c59ef15SJeremy L Thompson } 890c59ef15SJeremy L Thompson } 900c59ef15SJeremy L Thompson return -1; 910c59ef15SJeremy L Thompson } 920c59ef15SJeremy L Thompson static int CreateRestriction(Ceed ceed, const CeedInt melem[3], 930c59ef15SJeremy L Thompson CeedInt P, CeedInt ncomp, 940c59ef15SJeremy L Thompson CeedElemRestriction *Erestrict) { 95819eb1b3Sjeremylt const PetscInt nelem = melem[0]*melem[1]*melem[2]; 96b52ddc88Sjeremylt PetscInt mnodes[3], *idx, *idxp; 970c59ef15SJeremy L Thompson 98819eb1b3Sjeremylt // Get indicies 99b52ddc88Sjeremylt for (int d=0; d<3; d++) mnodes[d] = melem[d]*(P-1) + 1; 100819eb1b3Sjeremylt idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]); 1010c59ef15SJeremy L Thompson for (CeedInt i=0; i<melem[0]; i++) { 1020c59ef15SJeremy L Thompson for (CeedInt j=0; j<melem[1]; j++) { 1030c59ef15SJeremy L Thompson for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P) { 1040c59ef15SJeremy L Thompson for (CeedInt ii=0; ii<P; ii++) { 1050c59ef15SJeremy L Thompson for (CeedInt jj=0; jj<P; jj++) { 1060c59ef15SJeremy L Thompson for (CeedInt kk=0; kk<P; kk++) { 1070c59ef15SJeremy L Thompson if (0) { // This is the C-style (i,j,k) ordering that I prefer 108b52ddc88Sjeremylt idxp[(ii*P+jj)*P+kk] = (((i*(P-1)+ii)*mnodes[1] 109b52ddc88Sjeremylt + (j*(P-1)+jj))*mnodes[2] 1100c59ef15SJeremy L Thompson + (k*(P-1)+kk)); 1110c59ef15SJeremy L Thompson } else { // (k,j,i) ordering for consistency with MFEM example 112b52ddc88Sjeremylt idxp[ii+P*(jj+P*kk)] = (((i*(P-1)+ii)*mnodes[1] 113b52ddc88Sjeremylt + (j*(P-1)+jj))*mnodes[2] 1140c59ef15SJeremy L Thompson + (k*(P-1)+kk)); 1150c59ef15SJeremy L Thompson } 1160c59ef15SJeremy L Thompson } 1170c59ef15SJeremy L Thompson } 1180c59ef15SJeremy L Thompson } 1190c59ef15SJeremy L Thompson } 1200c59ef15SJeremy L Thompson } 1210c59ef15SJeremy L Thompson } 122819eb1b3Sjeremylt 123819eb1b3Sjeremylt // Setup CEED restriction 124*288c0443SJeremy L Thompson CeedElemRestrictionCreate(ceed, nelem, P*P*P, mnodes[0]*mnodes[1]*mnodes[2], 125*288c0443SJeremy L Thompson ncomp, 1260c59ef15SJeremy L Thompson CEED_MEM_HOST, CEED_OWN_POINTER, idx, Erestrict); 127819eb1b3Sjeremylt 1280c59ef15SJeremy L Thompson PetscFunctionReturn(0); 1290c59ef15SJeremy L Thompson } 1300c59ef15SJeremy L Thompson 1310c59ef15SJeremy L Thompson // Data for PETSc 1320c59ef15SJeremy L Thompson typedef struct User_ *User; 1330c59ef15SJeremy L Thompson struct User_ { 1340c59ef15SJeremy L Thompson MPI_Comm comm; 1350c59ef15SJeremy L Thompson VecScatter ltog; // Scatter for all entries 1360c59ef15SJeremy L Thompson VecScatter ltog0; // Skip Dirichlet values 1370c59ef15SJeremy L Thompson VecScatter gtogD; // global-to-global; only Dirichlet values 1380c59ef15SJeremy L Thompson Vec Xloc, Yloc; 1390c59ef15SJeremy L Thompson CeedVector xceed, yceed; 1400c59ef15SJeremy L Thompson CeedOperator op; 1410c59ef15SJeremy L Thompson CeedVector rho; 1420c59ef15SJeremy L Thompson Ceed ceed; 1430c59ef15SJeremy L Thompson }; 1440c59ef15SJeremy L Thompson 1450c59ef15SJeremy L Thompson // BP Options 1460c59ef15SJeremy L Thompson typedef enum { 1470c59ef15SJeremy L Thompson CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2, 1480c59ef15SJeremy L Thompson CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5 1490c59ef15SJeremy L Thompson } bpType; 1500c59ef15SJeremy L Thompson static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6", 1514d537eeaSYohann "bpType","CEED_BP",0 1524d537eeaSYohann }; 1530c59ef15SJeremy L Thompson 1540c59ef15SJeremy L Thompson // BP specific data 1550c59ef15SJeremy L Thompson typedef struct { 15634d77899SValeria Barra CeedInt ncompu, qdatasize, qextra; 1570c59ef15SJeremy L Thompson CeedQFunctionUser setup, apply, error; 1584d537eeaSYohann const char *setupfname, *applyfname, *errorfname; 1590c59ef15SJeremy L Thompson CeedEvalMode inmode, outmode; 160e07e9ddcSjeremylt CeedQuadMode qmode; 1610c59ef15SJeremy L Thompson } bpData; 1620c59ef15SJeremy L Thompson 1630c59ef15SJeremy L Thompson bpData bpOptions[6] = { 1640c59ef15SJeremy L Thompson [CEED_BP1] = { 16534d77899SValeria Barra .ncompu = 1, 1660c59ef15SJeremy L Thompson .qdatasize = 1, 167194c25f7Sjeremylt .qextra = 1, 1680c59ef15SJeremy L Thompson .setup = SetupMass, 1690c59ef15SJeremy L Thompson .apply = Mass, 1700c59ef15SJeremy L Thompson .error = Error, 1714d537eeaSYohann .setupfname = SetupMass_loc, 1724d537eeaSYohann .applyfname = Mass_loc, 1734d537eeaSYohann .errorfname = Error_loc, 1740c59ef15SJeremy L Thompson .inmode = CEED_EVAL_INTERP, 175e07e9ddcSjeremylt .outmode = CEED_EVAL_INTERP, 176f90c8643Sjeremylt .qmode = CEED_GAUSS 177f90c8643Sjeremylt }, 1780c59ef15SJeremy L Thompson [CEED_BP2] = { 17934d77899SValeria Barra .ncompu = 3, 1800c59ef15SJeremy L Thompson .qdatasize = 1, 181194c25f7Sjeremylt .qextra = 1, 1820c59ef15SJeremy L Thompson .setup = SetupMass3, 1830c59ef15SJeremy L Thompson .apply = Mass3, 1840c59ef15SJeremy L Thompson .error = Error3, 1854d537eeaSYohann .setupfname = SetupMass3_loc, 1864d537eeaSYohann .applyfname = Mass3_loc, 1874d537eeaSYohann .errorfname = Error3_loc, 1880c59ef15SJeremy L Thompson .inmode = CEED_EVAL_INTERP, 189e07e9ddcSjeremylt .outmode = CEED_EVAL_INTERP, 190f90c8643Sjeremylt .qmode = CEED_GAUSS 191f90c8643Sjeremylt }, 1920c59ef15SJeremy L Thompson [CEED_BP3] = { 19334d77899SValeria Barra .ncompu = 1, 1940c59ef15SJeremy L Thompson .qdatasize = 6, 195194c25f7Sjeremylt .qextra = 1, 1960c59ef15SJeremy L Thompson .setup = SetupDiff, 1970c59ef15SJeremy L Thompson .apply = Diff, 1980c59ef15SJeremy L Thompson .error = Error, 1994d537eeaSYohann .setupfname = SetupDiff_loc, 2004d537eeaSYohann .applyfname = Diff_loc, 2014d537eeaSYohann .errorfname = Error_loc, 2020c59ef15SJeremy L Thompson .inmode = CEED_EVAL_GRAD, 203e07e9ddcSjeremylt .outmode = CEED_EVAL_GRAD, 204f90c8643Sjeremylt .qmode = CEED_GAUSS 205f90c8643Sjeremylt }, 2060c59ef15SJeremy L Thompson [CEED_BP4] = { 20734d77899SValeria Barra .ncompu = 3, 2080c59ef15SJeremy L Thompson .qdatasize = 6, 209194c25f7Sjeremylt .qextra = 1, 2100c59ef15SJeremy L Thompson .setup = SetupDiff3, 2110c59ef15SJeremy L Thompson .apply = Diff3, 2124b5b4ec1Sjeremylt .error = Error3, 2134d537eeaSYohann .setupfname = SetupDiff3_loc, 2144d537eeaSYohann .applyfname = Diff3_loc, 2154d537eeaSYohann .errorfname = Error3_loc, 2160c59ef15SJeremy L Thompson .inmode = CEED_EVAL_GRAD, 217e07e9ddcSjeremylt .outmode = CEED_EVAL_GRAD, 218f90c8643Sjeremylt .qmode = CEED_GAUSS 219f90c8643Sjeremylt }, 2200c59ef15SJeremy L Thompson [CEED_BP5] = { 22134d77899SValeria Barra .ncompu = 1, 2220c59ef15SJeremy L Thompson .qdatasize = 6, 223194c25f7Sjeremylt .qextra = 0, 2240c59ef15SJeremy L Thompson .setup = SetupDiff, 2250c59ef15SJeremy L Thompson .apply = Diff, 2260c59ef15SJeremy L Thompson .error = Error, 2274d537eeaSYohann .setupfname = SetupDiff_loc, 2284d537eeaSYohann .applyfname = Diff_loc, 2294d537eeaSYohann .errorfname = Error_loc, 2300c59ef15SJeremy L Thompson .inmode = CEED_EVAL_GRAD, 231e07e9ddcSjeremylt .outmode = CEED_EVAL_GRAD, 232f90c8643Sjeremylt .qmode = CEED_GAUSS_LOBATTO 233f90c8643Sjeremylt }, 2340c59ef15SJeremy L Thompson [CEED_BP6] = { 23534d77899SValeria Barra .ncompu = 3, 2360c59ef15SJeremy L Thompson .qdatasize = 6, 237194c25f7Sjeremylt .qextra = 0, 2380c59ef15SJeremy L Thompson .setup = SetupDiff3, 2390c59ef15SJeremy L Thompson .apply = Diff3, 2404b5b4ec1Sjeremylt .error = Error3, 2414d537eeaSYohann .setupfname = SetupDiff3_loc, 2424d537eeaSYohann .applyfname = Diff3_loc, 2434d537eeaSYohann .errorfname = Error3_loc, 2440c59ef15SJeremy L Thompson .inmode = CEED_EVAL_GRAD, 245e07e9ddcSjeremylt .outmode = CEED_EVAL_GRAD, 246f90c8643Sjeremylt .qmode = CEED_GAUSS_LOBATTO 247f90c8643Sjeremylt } 2480c59ef15SJeremy L Thompson }; 2490c59ef15SJeremy L Thompson 2500c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the mass matrix 2510c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) { 2520c59ef15SJeremy L Thompson PetscErrorCode ierr; 2530c59ef15SJeremy L Thompson User user; 2540c59ef15SJeremy L Thompson PetscScalar *x, *y; 2550c59ef15SJeremy L Thompson 2560c59ef15SJeremy L Thompson PetscFunctionBeginUser; 2570c59ef15SJeremy L Thompson ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 2580c59ef15SJeremy L Thompson ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 2590c59ef15SJeremy L Thompson SCATTER_REVERSE); CHKERRQ(ierr); 2600c59ef15SJeremy L Thompson ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE); 2610c59ef15SJeremy L Thompson CHKERRQ(ierr); 2620c59ef15SJeremy L Thompson ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 2630c59ef15SJeremy L Thompson 2640c59ef15SJeremy L Thompson ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 2650c59ef15SJeremy L Thompson ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 2660c59ef15SJeremy L Thompson CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 2670c59ef15SJeremy L Thompson CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y); 2680c59ef15SJeremy L Thompson 2690c59ef15SJeremy L Thompson CeedOperatorApply(user->op, user->xceed, user->yceed, 2700c59ef15SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 2710c59ef15SJeremy L Thompson ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr); 2720c59ef15SJeremy L Thompson 2730c59ef15SJeremy L Thompson ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 2740c59ef15SJeremy L Thompson ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 2750c59ef15SJeremy L Thompson 2760c59ef15SJeremy L Thompson if (Y) { 2770c59ef15SJeremy L Thompson ierr = VecZeroEntries(Y); CHKERRQ(ierr); 2780c59ef15SJeremy L Thompson ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 2790c59ef15SJeremy L Thompson CHKERRQ(ierr); 2800c59ef15SJeremy L Thompson ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 2810c59ef15SJeremy L Thompson CHKERRQ(ierr); 2820c59ef15SJeremy L Thompson } 2830c59ef15SJeremy L Thompson PetscFunctionReturn(0); 2840c59ef15SJeremy L Thompson } 2850c59ef15SJeremy L Thompson 2860c59ef15SJeremy L Thompson // This function uses libCEED to compute the action of the Laplacian with 2870c59ef15SJeremy L Thompson // Dirichlet boundary conditions 2880c59ef15SJeremy L Thompson static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) { 2890c59ef15SJeremy L Thompson PetscErrorCode ierr; 2900c59ef15SJeremy L Thompson User user; 2910c59ef15SJeremy L Thompson PetscScalar *x, *y; 2920c59ef15SJeremy L Thompson 2930c59ef15SJeremy L Thompson PetscFunctionBeginUser; 2940c59ef15SJeremy L Thompson ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 295819eb1b3Sjeremylt 296819eb1b3Sjeremylt // Global-to-local 2970c59ef15SJeremy L Thompson ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES, 2980c59ef15SJeremy L Thompson SCATTER_REVERSE); CHKERRQ(ierr); 2990c59ef15SJeremy L Thompson ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES, 3000c59ef15SJeremy L Thompson SCATTER_REVERSE); 3010c59ef15SJeremy L Thompson CHKERRQ(ierr); 3020c59ef15SJeremy L Thompson ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 3030c59ef15SJeremy L Thompson 304819eb1b3Sjeremylt // Setup CEED vectors 3050c59ef15SJeremy L Thompson ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 3060c59ef15SJeremy L Thompson ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 3070c59ef15SJeremy L Thompson CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 3080c59ef15SJeremy L Thompson CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y); 3090c59ef15SJeremy L Thompson 310819eb1b3Sjeremylt // Apply CEED operator 3110c59ef15SJeremy L Thompson CeedOperatorApply(user->op, user->xceed, user->yceed, 3120c59ef15SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 3130c59ef15SJeremy L Thompson ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr); 3140c59ef15SJeremy L Thompson 315819eb1b3Sjeremylt // Restore PETSc vectors 3160c59ef15SJeremy L Thompson ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 3170c59ef15SJeremy L Thompson ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 3180c59ef15SJeremy L Thompson 319819eb1b3Sjeremylt // Local-to-global 3200c59ef15SJeremy L Thompson ierr = VecZeroEntries(Y); CHKERRQ(ierr); 3210c59ef15SJeremy L Thompson ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 3220c59ef15SJeremy L Thompson CHKERRQ(ierr); 3230c59ef15SJeremy L Thompson ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 3240c59ef15SJeremy L Thompson CHKERRQ(ierr); 3250c59ef15SJeremy L Thompson ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 3260c59ef15SJeremy L Thompson CHKERRQ(ierr); 3270c59ef15SJeremy L Thompson ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 3280c59ef15SJeremy L Thompson CHKERRQ(ierr); 329819eb1b3Sjeremylt 3300c59ef15SJeremy L Thompson PetscFunctionReturn(0); 3310c59ef15SJeremy L Thompson } 3320c59ef15SJeremy L Thompson 3330c59ef15SJeremy L Thompson // This function calculates the error in the final solution 3340c59ef15SJeremy L Thompson static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X, 3350c59ef15SJeremy L Thompson CeedVector target, PetscReal *maxerror) { 3360c59ef15SJeremy L Thompson PetscErrorCode ierr; 3370c59ef15SJeremy L Thompson PetscScalar *x; 3380c59ef15SJeremy L Thompson CeedVector collocated_error; 3390c59ef15SJeremy L Thompson CeedInt length; 3400c59ef15SJeremy L Thompson 3410c59ef15SJeremy L Thompson PetscFunctionBeginUser; 3420c59ef15SJeremy L Thompson CeedVectorGetLength(target, &length); 3430c59ef15SJeremy L Thompson CeedVectorCreate(user->ceed, length, &collocated_error); 344819eb1b3Sjeremylt 345819eb1b3Sjeremylt // Global-to-local 3460c59ef15SJeremy L Thompson ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 3470c59ef15SJeremy L Thompson SCATTER_REVERSE); CHKERRQ(ierr); 3480c59ef15SJeremy L Thompson ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE); 3490c59ef15SJeremy L Thompson CHKERRQ(ierr); 350819eb1b3Sjeremylt 351819eb1b3Sjeremylt // Setup CEED vector 3520c59ef15SJeremy L Thompson ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 3530c59ef15SJeremy L Thompson CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 354819eb1b3Sjeremylt 355819eb1b3Sjeremylt // Apply CEED operator 3560c59ef15SJeremy L Thompson CeedOperatorApply(op_error, user->xceed, collocated_error, 3570c59ef15SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 358819eb1b3Sjeremylt 359819eb1b3Sjeremylt // Restore PETSc vector 3600c59ef15SJeremy L Thompson VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 3610c59ef15SJeremy L Thompson 362819eb1b3Sjeremylt // Reduce max error 3630c59ef15SJeremy L Thompson *maxerror = 0; 3640c59ef15SJeremy L Thompson const CeedScalar *e; 3650c59ef15SJeremy L Thompson CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e); 3660c59ef15SJeremy L Thompson for (CeedInt i=0; i<length; i++) { 3670c59ef15SJeremy L Thompson *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i])); 3680c59ef15SJeremy L Thompson } 3690c59ef15SJeremy L Thompson CeedVectorRestoreArrayRead(collocated_error, &e); 37000742fa8SJed Brown ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror, 37100742fa8SJed Brown 1, MPIU_REAL, MPIU_MAX, user->comm); CHKERRQ(ierr); 372819eb1b3Sjeremylt 373819eb1b3Sjeremylt // Cleanup 3740c59ef15SJeremy L Thompson CeedVectorDestroy(&collocated_error); 375819eb1b3Sjeremylt 3760c59ef15SJeremy L Thompson PetscFunctionReturn(0); 3770c59ef15SJeremy L Thompson } 3780c59ef15SJeremy L Thompson 3790c59ef15SJeremy L Thompson int main(int argc, char **argv) { 3800c59ef15SJeremy L Thompson PetscInt ierr; 3810c59ef15SJeremy L Thompson MPI_Comm comm; 3820c59ef15SJeremy L Thompson char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 383819eb1b3Sjeremylt double my_rt_start, my_rt, rt_min, rt_max; 384b52ddc88Sjeremylt PetscInt degree, qextra, localnodes, localelem, melem[3], mnodes[3], p[3], 38534d77899SValeria Barra irank[3], lnodes[3], lsize, ncompu = 1; 3860c59ef15SJeremy L Thompson PetscScalar *r; 3876bd9afcaSjeremylt PetscBool test_mode, benchmark_mode, write_solution; 3880c59ef15SJeremy L Thompson PetscMPIInt size, rank; 389819eb1b3Sjeremylt Vec X, Xloc, rhs, rhsloc; 390819eb1b3Sjeremylt Mat mat; 391819eb1b3Sjeremylt KSP ksp; 3924b5b4ec1Sjeremylt VecScatter ltog, ltog0, gtogD; 393819eb1b3Sjeremylt User user; 3940c59ef15SJeremy L Thompson Ceed ceed; 3950c59ef15SJeremy L Thompson CeedBasis basisx, basisu; 3960c59ef15SJeremy L Thompson CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, Erestrictui, 3970c59ef15SJeremy L Thompson Erestrictqdi; 3980c59ef15SJeremy L Thompson CeedQFunction qf_setup, qf_apply, qf_error; 3990c59ef15SJeremy L Thompson CeedOperator op_setup, op_apply, op_error; 4000c59ef15SJeremy L Thompson CeedVector xcoord, rho, rhsceed, target; 4010c59ef15SJeremy L Thompson CeedInt P, Q; 40234d77899SValeria Barra const CeedInt dim = 3, ncompx = 3; 4030c59ef15SJeremy L Thompson bpType bpChoice; 4040c59ef15SJeremy L Thompson 4050c59ef15SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 4060c59ef15SJeremy L Thompson if (ierr) return ierr; 4070c59ef15SJeremy L Thompson comm = PETSC_COMM_WORLD; 4080c59ef15SJeremy L Thompson ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 4090c59ef15SJeremy L Thompson bpChoice = CEED_BP1; 4100c59ef15SJeremy L Thompson ierr = PetscOptionsEnum("-problem", 4110c59ef15SJeremy L Thompson "CEED benchmark problem to solve", NULL, 4120c59ef15SJeremy L Thompson bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice, 4130c59ef15SJeremy L Thompson NULL); CHKERRQ(ierr); 41434d77899SValeria Barra ncompu = bpOptions[bpChoice].ncompu; 4150c59ef15SJeremy L Thompson test_mode = PETSC_FALSE; 4160c59ef15SJeremy L Thompson ierr = PetscOptionsBool("-test", 4170c59ef15SJeremy L Thompson "Testing mode (do not print unless error is large)", 4180c59ef15SJeremy L Thompson NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 4190c59ef15SJeremy L Thompson benchmark_mode = PETSC_FALSE; 4200c59ef15SJeremy L Thompson ierr = PetscOptionsBool("-benchmark", 4210c59ef15SJeremy L Thompson "Benchmarking mode (prints benchmark statistics)", 4220c59ef15SJeremy L Thompson NULL, benchmark_mode, &benchmark_mode, NULL); 4230c59ef15SJeremy L Thompson CHKERRQ(ierr); 4246bd9afcaSjeremylt write_solution = PETSC_FALSE; 4256bd9afcaSjeremylt ierr = PetscOptionsBool("-write_solution", 4266bd9afcaSjeremylt "Write solution for visualization", 4276bd9afcaSjeremylt NULL, write_solution, &write_solution, NULL); 4286bd9afcaSjeremylt CHKERRQ(ierr); 4290c59ef15SJeremy L Thompson degree = test_mode ? 3 : 1; 4300c59ef15SJeremy L Thompson ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 4310c59ef15SJeremy L Thompson NULL, degree, °ree, NULL); CHKERRQ(ierr); 4320c59ef15SJeremy L Thompson qextra = bpOptions[bpChoice].qextra; 4330c59ef15SJeremy L Thompson ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 4340c59ef15SJeremy L Thompson NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 4350c59ef15SJeremy L Thompson ierr = PetscOptionsString("-ceed", "CEED resource specifier", 4360c59ef15SJeremy L Thompson NULL, ceedresource, ceedresource, 4370c59ef15SJeremy L Thompson sizeof(ceedresource), NULL); CHKERRQ(ierr); 438b52ddc88Sjeremylt localnodes = 1000; 4390c59ef15SJeremy L Thompson ierr = PetscOptionsInt("-local", 440432a1099Sjeremylt "Target number of locally owned nodes per process", 441b52ddc88Sjeremylt NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr); 4420c59ef15SJeremy L Thompson ierr = PetscOptionsEnd(); CHKERRQ(ierr); 4433ce86546SJeremy L Thompson P = degree + 1; 4443ce86546SJeremy L Thompson Q = P + qextra; 4450c59ef15SJeremy L Thompson 4460c59ef15SJeremy L Thompson // Determine size of process grid 4470c59ef15SJeremy L Thompson ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr); 4480c59ef15SJeremy L Thompson Split3(size, p, false); 4490c59ef15SJeremy L Thompson 450b52ddc88Sjeremylt // Find a nicely composite number of elements no less than localnodes 451b52ddc88Sjeremylt for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ; 4520c59ef15SJeremy L Thompson localelem++) { 4530c59ef15SJeremy L Thompson Split3(localelem, melem, true); 4540c59ef15SJeremy L Thompson if (Max3(melem) / Min3(melem) <= 2) break; 4550c59ef15SJeremy L Thompson } 4560c59ef15SJeremy L Thompson 4570c59ef15SJeremy L Thompson // Find my location in the process grid 4580c59ef15SJeremy L Thompson ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr); 45934d77899SValeria Barra for (int d=0,rankleft=rank; d<dim; d++) { 4600c59ef15SJeremy L Thompson const int pstride[3] = {p[1] *p[2], p[2], 1}; 4610c59ef15SJeremy L Thompson irank[d] = rankleft / pstride[d]; 4620c59ef15SJeremy L Thompson rankleft -= irank[d] * pstride[d]; 4630c59ef15SJeremy L Thompson } 4640c59ef15SJeremy L Thompson 465b52ddc88Sjeremylt GlobalNodes(p, irank, degree, melem, mnodes); 4660c59ef15SJeremy L Thompson 467819eb1b3Sjeremylt // Setup global vector 4680c59ef15SJeremy L Thompson ierr = VecCreate(comm, &X); CHKERRQ(ierr); 46934d77899SValeria Barra ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE); 4700c59ef15SJeremy L Thompson CHKERRQ(ierr); 4710c59ef15SJeremy L Thompson ierr = VecSetUp(X); CHKERRQ(ierr); 4720c59ef15SJeremy L Thompson 473819eb1b3Sjeremylt // Print summary 4740c59ef15SJeremy L Thompson if (!test_mode) { 4750c59ef15SJeremy L Thompson CeedInt gsize; 4760c59ef15SJeremy L Thompson ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 4773ce86546SJeremy L Thompson ierr = PetscPrintf(comm, 4783ce86546SJeremy L Thompson "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n" 4793ce86546SJeremy L Thompson " libCEED:\n" 4803ce86546SJeremy L Thompson " libCEED Backend : %s\n" 4813ce86546SJeremy L Thompson " Mesh:\n" 4823ce86546SJeremy L Thompson " Number of 1D Basis Nodes (p) : %d\n" 4833ce86546SJeremy L Thompson " Number of 1D Quadrature Points (q) : %d\n" 484b52ddc88Sjeremylt " Global nodes : %D\n" 4853ce86546SJeremy L Thompson " Process Decomposition : %D %D %D\n" 4863ce86546SJeremy L Thompson " Local Elements : %D = %D %D %D\n" 487b52ddc88Sjeremylt " Owned nodes : %D = %D %D %D\n", 48834d77899SValeria Barra bpChoice+1, ceedresource, P, Q, gsize/ncompu, p[0], 4893ce86546SJeremy L Thompson p[1], p[2], localelem, melem[0], melem[1], melem[2], 490b52ddc88Sjeremylt mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1], mnodes[2]); 4913ce86546SJeremy L Thompson CHKERRQ(ierr); 4920c59ef15SJeremy L Thompson } 4930c59ef15SJeremy L Thompson 4940c59ef15SJeremy L Thompson { 4950c59ef15SJeremy L Thompson lsize = 1; 49634d77899SValeria Barra for (int d=0; d<dim; d++) { 497b52ddc88Sjeremylt lnodes[d] = melem[d]*degree + 1; 498b52ddc88Sjeremylt lsize *= lnodes[d]; 4990c59ef15SJeremy L Thompson } 5000c59ef15SJeremy L Thompson ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr); 50134d77899SValeria Barra ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr); 5020c59ef15SJeremy L Thompson ierr = VecSetUp(Xloc); CHKERRQ(ierr); 5030c59ef15SJeremy L Thompson 5040c59ef15SJeremy L Thompson // Create local-to-global scatter 5050c59ef15SJeremy L Thompson PetscInt *ltogind, *ltogind0, *locind, l0count; 5060c59ef15SJeremy L Thompson IS ltogis, ltogis0, locis; 50734d77899SValeria Barra PetscInt gstart[2][2][2], gmnodes[2][2][2][dim]; 5080c59ef15SJeremy L Thompson 5090c59ef15SJeremy L Thompson for (int i=0; i<2; i++) { 5100c59ef15SJeremy L Thompson for (int j=0; j<2; j++) { 5110c59ef15SJeremy L Thompson for (int k=0; k<2; k++) { 5120c59ef15SJeremy L Thompson PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k}; 5130c59ef15SJeremy L Thompson gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem); 514b52ddc88Sjeremylt GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]); 5150c59ef15SJeremy L Thompson } 5160c59ef15SJeremy L Thompson } 5170c59ef15SJeremy L Thompson } 5180c59ef15SJeremy L Thompson 5190c59ef15SJeremy L Thompson ierr = PetscMalloc1(lsize, <ogind); CHKERRQ(ierr); 5200c59ef15SJeremy L Thompson ierr = PetscMalloc1(lsize, <ogind0); CHKERRQ(ierr); 5210c59ef15SJeremy L Thompson ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr); 5220c59ef15SJeremy L Thompson l0count = 0; 523*288c0443SJeremy L Thompson for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++) 524*288c0443SJeremy L Thompson for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++) 525b52ddc88Sjeremylt for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) { 526b52ddc88Sjeremylt PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k; 5270c59ef15SJeremy L Thompson ltogind[here] = 528b52ddc88Sjeremylt gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk; 5290c59ef15SJeremy L Thompson if ((irank[0] == 0 && i == 0) 5300c59ef15SJeremy L Thompson || (irank[1] == 0 && j == 0) 5310c59ef15SJeremy L Thompson || (irank[2] == 0 && k == 0) 532b52ddc88Sjeremylt || (irank[0]+1 == p[0] && i+1 == lnodes[0]) 533b52ddc88Sjeremylt || (irank[1]+1 == p[1] && j+1 == lnodes[1]) 534b52ddc88Sjeremylt || (irank[2]+1 == p[2] && k+1 == lnodes[2])) 5350c59ef15SJeremy L Thompson continue; 5360c59ef15SJeremy L Thompson ltogind0[l0count] = ltogind[here]; 5370c59ef15SJeremy L Thompson locind[l0count++] = here; 5380c59ef15SJeremy L Thompson } 53934d77899SValeria Barra ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER, 5400c59ef15SJeremy L Thompson <ogis); CHKERRQ(ierr); 5410c59ef15SJeremy L Thompson ierr = VecScatterCreate(Xloc, NULL, X, ltogis, <og); CHKERRQ(ierr); 5420c59ef15SJeremy L Thompson CHKERRQ(ierr); 54334d77899SValeria Barra ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER, 5440c59ef15SJeremy L Thompson <ogis0); CHKERRQ(ierr); 54534d77899SValeria Barra ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER, 5460c59ef15SJeremy L Thompson &locis); CHKERRQ(ierr); 5470c59ef15SJeremy L Thompson ierr = VecScatterCreate(Xloc, locis, X, ltogis0, <og0); CHKERRQ(ierr); 5480c59ef15SJeremy L Thompson { 5490c59ef15SJeremy L Thompson // Create global-to-global scatter for Dirichlet values (everything not in 5500c59ef15SJeremy L Thompson // ltogis0, which is the range of ltog0) 5510c59ef15SJeremy L Thompson PetscInt xstart, xend, *indD, countD = 0; 5520c59ef15SJeremy L Thompson IS isD; 5530c59ef15SJeremy L Thompson const PetscScalar *x; 5540c59ef15SJeremy L Thompson ierr = VecZeroEntries(Xloc); CHKERRQ(ierr); 5550c59ef15SJeremy L Thompson ierr = VecSet(X, 1.0); CHKERRQ(ierr); 5560c59ef15SJeremy L Thompson ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 5570c59ef15SJeremy L Thompson CHKERRQ(ierr); 5580c59ef15SJeremy L Thompson ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 5590c59ef15SJeremy L Thompson CHKERRQ(ierr); 5600c59ef15SJeremy L Thompson ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr); 5610c59ef15SJeremy L Thompson ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr); 5620c59ef15SJeremy L Thompson ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr); 5630c59ef15SJeremy L Thompson for (PetscInt i=0; i<xend-xstart; i++) { 5640c59ef15SJeremy L Thompson if (x[i] == 1.) indD[countD++] = xstart + i; 5650c59ef15SJeremy L Thompson } 5660c59ef15SJeremy L Thompson ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr); 5670c59ef15SJeremy L Thompson ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD); 5680c59ef15SJeremy L Thompson CHKERRQ(ierr); 5690c59ef15SJeremy L Thompson ierr = PetscFree(indD); CHKERRQ(ierr); 5700c59ef15SJeremy L Thompson ierr = VecScatterCreate(X, isD, X, isD, >ogD); CHKERRQ(ierr); 5710c59ef15SJeremy L Thompson ierr = ISDestroy(&isD); CHKERRQ(ierr); 5720c59ef15SJeremy L Thompson } 5730c59ef15SJeremy L Thompson ierr = ISDestroy(<ogis); CHKERRQ(ierr); 5740c59ef15SJeremy L Thompson ierr = ISDestroy(<ogis0); CHKERRQ(ierr); 5750c59ef15SJeremy L Thompson ierr = ISDestroy(&locis); CHKERRQ(ierr); 5760c59ef15SJeremy L Thompson } 5770c59ef15SJeremy L Thompson 5780c59ef15SJeremy L Thompson // Set up libCEED 5790c59ef15SJeremy L Thompson CeedInit(ceedresource, &ceed); 58034d77899SValeria Barra 58134d77899SValeria Barra // CEED bases 58234d77899SValeria Barra CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q, 583e07e9ddcSjeremylt bpOptions[bpChoice].qmode, &basisu); 58434d77899SValeria Barra CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q, 585e07e9ddcSjeremylt bpOptions[bpChoice].qmode, &basisx); 5860c59ef15SJeremy L Thompson 58734d77899SValeria Barra // CEED restrictions 58834d77899SValeria Barra CreateRestriction(ceed, melem, P, ncompu, &Erestrictu); 58934d77899SValeria Barra CreateRestriction(ceed, melem, 2, dim, &Erestrictx); 5900c59ef15SJeremy L Thompson CeedInt nelem = melem[0]*melem[1]*melem[2]; 59134d77899SValeria Barra CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, ncompu, 5920c59ef15SJeremy L Thompson &Erestrictui); 5930c59ef15SJeremy L Thompson CeedElemRestrictionCreateIdentity(ceed, nelem, 594ea339264SYohann Q*Q*Q, 595ea339264SYohann nelem*Q*Q*Q, 596ea339264SYohann bpOptions[bpChoice].qdatasize, &Erestrictqdi); 5970c59ef15SJeremy L Thompson CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, 1, 5980c59ef15SJeremy L Thompson &Erestrictxi); 5990c59ef15SJeremy L Thompson { 6000c59ef15SJeremy L Thompson CeedScalar *xloc; 6010c59ef15SJeremy L Thompson CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len = 6020c59ef15SJeremy L Thompson shape[0]*shape[1]*shape[2]; 60334d77899SValeria Barra xloc = malloc(len*ncompx*sizeof xloc[0]); 6040c59ef15SJeremy L Thompson for (CeedInt i=0; i<shape[0]; i++) { 6050c59ef15SJeremy L Thompson for (CeedInt j=0; j<shape[1]; j++) { 6060c59ef15SJeremy L Thompson for (CeedInt k=0; k<shape[2]; k++) { 6070c59ef15SJeremy L Thompson xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) / 6080c59ef15SJeremy L Thompson (p[0]*melem[0]); 6090c59ef15SJeremy L Thompson xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) / 6100c59ef15SJeremy L Thompson (p[1]*melem[1]); 6110c59ef15SJeremy L Thompson xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) / 6120c59ef15SJeremy L Thompson (p[2]*melem[2]); 6130c59ef15SJeremy L Thompson } 6140c59ef15SJeremy L Thompson } 6150c59ef15SJeremy L Thompson } 61634d77899SValeria Barra CeedVectorCreate(ceed, len*ncompx, &xcoord); 6170c59ef15SJeremy L Thompson CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc); 6180c59ef15SJeremy L Thompson } 6190c59ef15SJeremy L Thompson 6200c59ef15SJeremy L Thompson // Create the Q-function that builds the operator (i.e. computes its 6214b5b4ec1Sjeremylt // quadrature data) and set its context data 6220c59ef15SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setup, 6230c59ef15SJeremy L Thompson bpOptions[bpChoice].setupfname, &qf_setup); 62434d77899SValeria Barra CeedQFunctionAddInput(qf_setup, "x", ncompx, CEED_EVAL_INTERP); 62534d77899SValeria Barra CeedQFunctionAddInput(qf_setup, "dx", ncompx*dim, CEED_EVAL_GRAD); 6260c59ef15SJeremy L Thompson CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 6270c59ef15SJeremy L Thompson CeedQFunctionAddOutput(qf_setup, "rho", bpOptions[bpChoice].qdatasize, 6280c59ef15SJeremy L Thompson CEED_EVAL_NONE); 62934d77899SValeria Barra CeedQFunctionAddOutput(qf_setup, "true_soln", ncompu, CEED_EVAL_NONE); 63034d77899SValeria Barra CeedQFunctionAddOutput(qf_setup, "rhs", ncompu, CEED_EVAL_INTERP); 6310c59ef15SJeremy L Thompson 6320c59ef15SJeremy L Thompson // Set up PDE operator 6330c59ef15SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].apply, 6340c59ef15SJeremy L Thompson bpOptions[bpChoice].applyfname, &qf_apply); 6350c59ef15SJeremy L Thompson // Add inputs and outputs 6364d537eeaSYohann CeedInt gradInScale = bpOptions[bpChoice].inmode==CEED_EVAL_GRAD ? 3 : 1; 6374d537eeaSYohann CeedInt gradOutScale = bpOptions[bpChoice].outmode==CEED_EVAL_GRAD ? 3 : 1; 63834d77899SValeria Barra CeedQFunctionAddInput(qf_apply, "u", ncompu*gradInScale, 6394d537eeaSYohann bpOptions[bpChoice].inmode); 6400c59ef15SJeremy L Thompson CeedQFunctionAddInput(qf_apply, "rho", bpOptions[bpChoice].qdatasize, 6410c59ef15SJeremy L Thompson CEED_EVAL_NONE); 64234d77899SValeria Barra CeedQFunctionAddOutput(qf_apply, "v", ncompu*gradOutScale, 6434d537eeaSYohann bpOptions[bpChoice].outmode); 6440c59ef15SJeremy L Thompson 6450c59ef15SJeremy L Thompson // Create the error qfunction 6460c59ef15SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error, 6470c59ef15SJeremy L Thompson bpOptions[bpChoice].errorfname, &qf_error); 64834d77899SValeria Barra CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP); 64934d77899SValeria Barra CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE); 65034d77899SValeria Barra CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE); 6510c59ef15SJeremy L Thompson 6520c59ef15SJeremy L Thompson // Create the persistent vectors that will be needed in setup 653819eb1b3Sjeremylt CeedInt nqpts; 654819eb1b3Sjeremylt CeedBasisGetNumQuadraturePoints(basisu, &nqpts); 655819eb1b3Sjeremylt CeedVectorCreate(ceed, bpOptions[bpChoice].qdatasize*nelem*nqpts, &rho); 65634d77899SValeria Barra CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target); 65734d77899SValeria Barra CeedVectorCreate(ceed, lsize*ncompu, &rhsceed); 6580c59ef15SJeremy L Thompson 6594b5b4ec1Sjeremylt // Create the operator that builds the quadrature data for the ceed operator 6600c59ef15SJeremy L Thompson CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup); 6610c59ef15SJeremy L Thompson CeedOperatorSetField(op_setup, "x", Erestrictx, CEED_NOTRANSPOSE, 6620c59ef15SJeremy L Thompson basisx, CEED_VECTOR_ACTIVE); 6630c59ef15SJeremy L Thompson CeedOperatorSetField(op_setup, "dx", Erestrictx, CEED_NOTRANSPOSE, 6640c59ef15SJeremy L Thompson basisx, CEED_VECTOR_ACTIVE); 6650c59ef15SJeremy L Thompson CeedOperatorSetField(op_setup, "weight", Erestrictxi, CEED_NOTRANSPOSE, 6660c59ef15SJeremy L Thompson basisx, CEED_VECTOR_NONE); 6670c59ef15SJeremy L Thompson CeedOperatorSetField(op_setup, "rho", Erestrictqdi, CEED_NOTRANSPOSE, 6680c59ef15SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 6690c59ef15SJeremy L Thompson CeedOperatorSetField(op_setup, "true_soln", Erestrictui, CEED_NOTRANSPOSE, 6700c59ef15SJeremy L Thompson CEED_BASIS_COLLOCATED, target); 6710c59ef15SJeremy L Thompson CeedOperatorSetField(op_setup, "rhs", Erestrictu, CEED_TRANSPOSE, 6720c59ef15SJeremy L Thompson basisu, rhsceed); 6730c59ef15SJeremy L Thompson 6744b5b4ec1Sjeremylt // Create the mass or diff operator 6750c59ef15SJeremy L Thompson CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply); 6760c59ef15SJeremy L Thompson CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE, 6770c59ef15SJeremy L Thompson basisu, CEED_VECTOR_ACTIVE); 6780c59ef15SJeremy L Thompson CeedOperatorSetField(op_apply, "rho", Erestrictqdi, CEED_NOTRANSPOSE, 6790c59ef15SJeremy L Thompson CEED_BASIS_COLLOCATED, rho); 6800c59ef15SJeremy L Thompson CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE, 6810c59ef15SJeremy L Thompson basisu, CEED_VECTOR_ACTIVE); 6820c59ef15SJeremy L Thompson 6830c59ef15SJeremy L Thompson // Create the error operator 6840c59ef15SJeremy L Thompson CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 6850c59ef15SJeremy L Thompson CeedOperatorSetField(op_error, "u", Erestrictu, CEED_TRANSPOSE, 6860c59ef15SJeremy L Thompson basisu, CEED_VECTOR_ACTIVE); 6870c59ef15SJeremy L Thompson CeedOperatorSetField(op_error, "true_soln", Erestrictui, CEED_NOTRANSPOSE, 6880c59ef15SJeremy L Thompson CEED_BASIS_COLLOCATED, target); 6890c59ef15SJeremy L Thompson CeedOperatorSetField(op_error, "error", Erestrictui, CEED_NOTRANSPOSE, 6900c59ef15SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 6910c59ef15SJeremy L Thompson 6920c59ef15SJeremy L Thompson 6930c59ef15SJeremy L Thompson // Set up Mat 6940c59ef15SJeremy L Thompson ierr = PetscMalloc1(1, &user); CHKERRQ(ierr); 6950c59ef15SJeremy L Thompson user->comm = comm; 6960c59ef15SJeremy L Thompson user->ltog = ltog; 6970c59ef15SJeremy L Thompson if (bpChoice != CEED_BP1 && bpChoice != CEED_BP2) { 6980c59ef15SJeremy L Thompson user->ltog0 = ltog0; 6990c59ef15SJeremy L Thompson user->gtogD = gtogD; 7000c59ef15SJeremy L Thompson } 7010c59ef15SJeremy L Thompson user->Xloc = Xloc; 7020c59ef15SJeremy L Thompson ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr); 70334d77899SValeria Barra CeedVectorCreate(ceed, lsize*ncompu, &user->xceed); 70434d77899SValeria Barra CeedVectorCreate(ceed, lsize*ncompu, &user->yceed); 7050c59ef15SJeremy L Thompson user->op = op_apply; 7060c59ef15SJeremy L Thompson user->rho = rho; 7070c59ef15SJeremy L Thompson user->ceed = ceed; 7080c59ef15SJeremy L Thompson 70934d77899SValeria Barra ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 71034d77899SValeria Barra mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 7110c59ef15SJeremy L Thompson PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr); 7120c59ef15SJeremy L Thompson if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) { 7130c59ef15SJeremy L Thompson ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass); 7140c59ef15SJeremy L Thompson CHKERRQ(ierr); 7150c59ef15SJeremy L Thompson } else { 7160c59ef15SJeremy L Thompson ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff); 7170c59ef15SJeremy L Thompson CHKERRQ(ierr); 7180c59ef15SJeremy L Thompson } 7190c59ef15SJeremy L Thompson ierr = MatCreateVecs(mat, &rhs, NULL); CHKERRQ(ierr); 7200c59ef15SJeremy L Thompson 7210c59ef15SJeremy L Thompson // Get RHS vector 7220c59ef15SJeremy L Thompson ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 7230c59ef15SJeremy L Thompson ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 7240c59ef15SJeremy L Thompson ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 7250c59ef15SJeremy L Thompson CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r); 7260c59ef15SJeremy L Thompson 7270c59ef15SJeremy L Thompson // Setup rho, rhs, and target 7280c59ef15SJeremy L Thompson CeedOperatorApply(op_setup, xcoord, rho, CEED_REQUEST_IMMEDIATE); 7290c59ef15SJeremy L Thompson ierr = CeedVectorSyncArray(rhsceed, CEED_MEM_HOST); CHKERRQ(ierr); 7300c59ef15SJeremy L Thompson CeedVectorDestroy(&xcoord); 7310c59ef15SJeremy L Thompson 7320c59ef15SJeremy L Thompson // Gather RHS 7330c59ef15SJeremy L Thompson ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 7340c59ef15SJeremy L Thompson ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 7350c59ef15SJeremy L Thompson ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 7360c59ef15SJeremy L Thompson CHKERRQ(ierr); 7370c59ef15SJeremy L Thompson ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 7380c59ef15SJeremy L Thompson CHKERRQ(ierr); 7390c59ef15SJeremy L Thompson CeedVectorDestroy(&rhsceed); 7400c59ef15SJeremy L Thompson 7410c59ef15SJeremy L Thompson ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 7420c59ef15SJeremy L Thompson { 7430c59ef15SJeremy L Thompson PC pc; 7440c59ef15SJeremy L Thompson ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 7450c59ef15SJeremy L Thompson if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) { 7460c59ef15SJeremy L Thompson ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 7470c59ef15SJeremy L Thompson ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 7480c59ef15SJeremy L Thompson } else { 7490c59ef15SJeremy L Thompson ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 7500c59ef15SJeremy L Thompson } 7510c59ef15SJeremy L Thompson ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 7520c59ef15SJeremy L Thompson ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 7530c59ef15SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 7540c59ef15SJeremy L Thompson PETSC_DEFAULT); CHKERRQ(ierr); 7550c59ef15SJeremy L Thompson } 7560c59ef15SJeremy L Thompson ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 7570c59ef15SJeremy L Thompson ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr); 7580c59ef15SJeremy L Thompson // First run, if benchmarking 7590c59ef15SJeremy L Thompson if (benchmark_mode) { 7600c59ef15SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 7610c59ef15SJeremy L Thompson CHKERRQ(ierr); 7620c59ef15SJeremy L Thompson my_rt_start = MPI_Wtime(); 7630c59ef15SJeremy L Thompson ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 7640c59ef15SJeremy L Thompson my_rt = MPI_Wtime() - my_rt_start; 7650c59ef15SJeremy L Thompson // Set maxits based on first iteration timing 7660c59ef15SJeremy L Thompson if (my_rt > 0.02) { 7670c59ef15SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 7680c59ef15SJeremy L Thompson CHKERRQ(ierr); 7690c59ef15SJeremy L Thompson } else { 7700c59ef15SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 7710c59ef15SJeremy L Thompson CHKERRQ(ierr); 7720c59ef15SJeremy L Thompson } 7730c59ef15SJeremy L Thompson } 7740c59ef15SJeremy L Thompson // Timed solve 7750c59ef15SJeremy L Thompson my_rt_start = MPI_Wtime(); 7760c59ef15SJeremy L Thompson ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 7770c59ef15SJeremy L Thompson my_rt = MPI_Wtime() - my_rt_start; 7780c59ef15SJeremy L Thompson { 7790c59ef15SJeremy L Thompson KSPType ksptype; 7800c59ef15SJeremy L Thompson KSPConvergedReason reason; 7810c59ef15SJeremy L Thompson PetscReal rnorm; 7820c59ef15SJeremy L Thompson PetscInt its; 7830c59ef15SJeremy L Thompson ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 7840c59ef15SJeremy L Thompson ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 7850c59ef15SJeremy L Thompson ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 7860c59ef15SJeremy L Thompson ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 7870c59ef15SJeremy L Thompson if (!test_mode || reason < 0 || rnorm > 1e-8) { 7883ce86546SJeremy L Thompson ierr = PetscPrintf(comm, 7893ce86546SJeremy L Thompson " KSP:\n" 7903ce86546SJeremy L Thompson " KSP Type : %s\n" 7913ce86546SJeremy L Thompson " KSP Convergence : %s\n" 7923ce86546SJeremy L Thompson " Total KSP Iterations : %D\n" 7933ce86546SJeremy L Thompson " Final rnorm : %e\n", 7943ce86546SJeremy L Thompson ksptype, KSPConvergedReasons[reason], its, 7953ce86546SJeremy L Thompson (double)rnorm); CHKERRQ(ierr); 7960c59ef15SJeremy L Thompson } 7970c59ef15SJeremy L Thompson if (benchmark_mode && (!test_mode)) { 7980c59ef15SJeremy L Thompson CeedInt gsize; 7990c59ef15SJeremy L Thompson ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 8000c59ef15SJeremy L Thompson MPI_Reduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, 0, comm); 8010c59ef15SJeremy L Thompson MPI_Reduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, 0, comm); 8020c59ef15SJeremy L Thompson ierr = PetscPrintf(comm, 8033ce86546SJeremy L Thompson " Performance:\n" 8043ce86546SJeremy L Thompson " CG Solve Time : %g (%g) sec\n" 805b52ddc88Sjeremylt " DoFs/Sec in CG : %g (%g) million\n", 8063ce86546SJeremy L Thompson rt_max, rt_min, 1e-6*gsize*its/rt_max, 8073ce86546SJeremy L Thompson 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 8080c59ef15SJeremy L Thompson } 8090c59ef15SJeremy L Thompson } 8100c59ef15SJeremy L Thompson 8110c59ef15SJeremy L Thompson { 8120c59ef15SJeremy L Thompson PetscReal maxerror; 8130c59ef15SJeremy L Thompson ierr = ComputeErrorMax(user, op_error, X, target, &maxerror); CHKERRQ(ierr); 8140c59ef15SJeremy L Thompson PetscReal tol = (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) ? 5e-3 : 5e-2; 8150c59ef15SJeremy L Thompson if (!test_mode || maxerror > tol) { 8163ce86546SJeremy L Thompson ierr = PetscPrintf(comm, 8173ce86546SJeremy L Thompson " Pointwise Error (max) : %e\n", 8183ce86546SJeremy L Thompson (double)maxerror); CHKERRQ(ierr); 8190c59ef15SJeremy L Thompson } 8200c59ef15SJeremy L Thompson } 8210c59ef15SJeremy L Thompson 8226bd9afcaSjeremylt if (write_solution) { 8236bd9afcaSjeremylt PetscViewer vtkviewersoln; 8246bd9afcaSjeremylt 8256bd9afcaSjeremylt ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 8266bd9afcaSjeremylt ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 8276bd9afcaSjeremylt ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 8286bd9afcaSjeremylt ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 8296bd9afcaSjeremylt ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 8306bd9afcaSjeremylt } 8316bd9afcaSjeremylt 8320c59ef15SJeremy L Thompson ierr = VecDestroy(&rhs); CHKERRQ(ierr); 8330c59ef15SJeremy L Thompson ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 8340c59ef15SJeremy L Thompson ierr = VecDestroy(&X); CHKERRQ(ierr); 8350c59ef15SJeremy L Thompson ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr); 8360c59ef15SJeremy L Thompson ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr); 8370c59ef15SJeremy L Thompson ierr = VecScatterDestroy(<og); CHKERRQ(ierr); 8380c59ef15SJeremy L Thompson ierr = VecScatterDestroy(<og0); CHKERRQ(ierr); 8390c59ef15SJeremy L Thompson ierr = VecScatterDestroy(>ogD); CHKERRQ(ierr); 8400c59ef15SJeremy L Thompson ierr = MatDestroy(&mat); CHKERRQ(ierr); 8410c59ef15SJeremy L Thompson ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 8420c59ef15SJeremy L Thompson 8430c59ef15SJeremy L Thompson CeedVectorDestroy(&user->xceed); 8440c59ef15SJeremy L Thompson CeedVectorDestroy(&user->yceed); 8450c59ef15SJeremy L Thompson CeedVectorDestroy(&user->rho); 8460c59ef15SJeremy L Thompson CeedVectorDestroy(&target); 8470c59ef15SJeremy L Thompson CeedOperatorDestroy(&op_setup); 8480c59ef15SJeremy L Thompson CeedOperatorDestroy(&op_apply); 8490c59ef15SJeremy L Thompson CeedOperatorDestroy(&op_error); 8500c59ef15SJeremy L Thompson CeedElemRestrictionDestroy(&Erestrictu); 8510c59ef15SJeremy L Thompson CeedElemRestrictionDestroy(&Erestrictx); 8520c59ef15SJeremy L Thompson CeedElemRestrictionDestroy(&Erestrictui); 8530c59ef15SJeremy L Thompson CeedElemRestrictionDestroy(&Erestrictxi); 8540c59ef15SJeremy L Thompson CeedElemRestrictionDestroy(&Erestrictqdi); 8550c59ef15SJeremy L Thompson CeedQFunctionDestroy(&qf_setup); 8560c59ef15SJeremy L Thompson CeedQFunctionDestroy(&qf_apply); 8570c59ef15SJeremy L Thompson CeedQFunctionDestroy(&qf_error); 8580c59ef15SJeremy L Thompson CeedBasisDestroy(&basisu); 8590c59ef15SJeremy L Thompson CeedBasisDestroy(&basisx); 8600c59ef15SJeremy L Thompson CeedDestroy(&ceed); 8610c59ef15SJeremy L Thompson ierr = PetscFree(user); CHKERRQ(ierr); 8620c59ef15SJeremy L Thompson return PetscFinalize(); 8630c59ef15SJeremy L Thompson } 864