1cb32e2e7SValeria Barra // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2cb32e2e7SValeria Barra // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3cb32e2e7SValeria Barra // reserved. See files LICENSE and NOTICE for details. 4cb32e2e7SValeria Barra // 5cb32e2e7SValeria Barra // This file is part of CEED, a collection of benchmarks, miniapps, software 6cb32e2e7SValeria Barra // libraries and APIs for efficient high-order finite element and spectral 7cb32e2e7SValeria Barra // element discretizations for exascale applications. For more information and 8cb32e2e7SValeria Barra // source code availability see http://github.com/ceed. 9cb32e2e7SValeria Barra // 10cb32e2e7SValeria Barra // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11cb32e2e7SValeria Barra // a collaborative effort of two U.S. Department of Energy organizations (Office 12cb32e2e7SValeria Barra // of Science and the National Nuclear Security Administration) responsible for 13cb32e2e7SValeria Barra // the planning and preparation of a capable exascale ecosystem, including 14cb32e2e7SValeria Barra // software, applications, hardware, advanced system engineering and early 15cb32e2e7SValeria Barra // testbed platforms, in support of the nation's exascale computing imperative. 16cb32e2e7SValeria Barra 17cb32e2e7SValeria Barra // libCEED + PETSc Example: CEED BPs 18cb32e2e7SValeria Barra // 19cb32e2e7SValeria Barra // This example demonstrates a simple usage of libCEED with PETSc to solve the 20cb32e2e7SValeria Barra // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps. 21cb32e2e7SValeria Barra // 22cb32e2e7SValeria Barra // The code is intentionally "raw", using only low-level communication 23cb32e2e7SValeria Barra // primitives. 24cb32e2e7SValeria Barra // 25cb32e2e7SValeria Barra // Build with: 26cb32e2e7SValeria Barra // 27cb32e2e7SValeria Barra // make bpsraw [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 28cb32e2e7SValeria Barra // 29cb32e2e7SValeria Barra // Sample runs: 30cb32e2e7SValeria Barra // 31cb32e2e7SValeria Barra // bpsraw -problem bp1 32cb32e2e7SValeria Barra // bpsraw -problem bp2 -ceed /cpu/self 33cb32e2e7SValeria Barra // bpsraw -problem bp3 -ceed /gpu/occa 34cb32e2e7SValeria Barra // bpsraw -problem bp4 -ceed /cpu/occa 35cb32e2e7SValeria Barra // bpsraw -problem bp5 -ceed /omp/occa 36cb32e2e7SValeria Barra // bpsraw -problem bp6 -ceed /ocl/occa 37cb32e2e7SValeria Barra // 38cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -qextra 5 39cb32e2e7SValeria Barra 40cb32e2e7SValeria Barra /// @file 41cb32e2e7SValeria Barra /// CEED BPs example using PETSc 42cb32e2e7SValeria Barra /// See bps.c for an implementation using DMPlex unstructured grids. 43cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc\n"; 44cb32e2e7SValeria Barra 45cb32e2e7SValeria Barra #include <stdbool.h> 46cb32e2e7SValeria Barra #include <string.h> 47cb32e2e7SValeria Barra #include <petscksp.h> 48cb32e2e7SValeria Barra #include <ceed.h> 49cb32e2e7SValeria Barra #include "qfunctions/bps/common.h" 50cb32e2e7SValeria Barra #include "qfunctions/bps/bp1.h" 51cb32e2e7SValeria Barra #include "qfunctions/bps/bp2.h" 52cb32e2e7SValeria Barra #include "qfunctions/bps/bp3.h" 53cb32e2e7SValeria Barra #include "qfunctions/bps/bp4.h" 54cb32e2e7SValeria Barra 55cb32e2e7SValeria Barra static void Split3(PetscInt size, PetscInt m[3], bool reverse) { 56cb32e2e7SValeria Barra for (PetscInt d=0,sizeleft=size; d<3; d++) { 57cb32e2e7SValeria Barra PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(sizeleft, 1./(3 - d))); 58cb32e2e7SValeria Barra while (try * (sizeleft / try) != sizeleft) try++; 59cb32e2e7SValeria Barra m[reverse ? 2-d : d] = try; 60cb32e2e7SValeria Barra sizeleft /= try; 61cb32e2e7SValeria Barra } 62cb32e2e7SValeria Barra } 63cb32e2e7SValeria Barra 64cb32e2e7SValeria Barra static PetscInt Max3(const PetscInt a[3]) { 65cb32e2e7SValeria Barra return PetscMax(a[0], PetscMax(a[1], a[2])); 66cb32e2e7SValeria Barra } 67cb32e2e7SValeria Barra static PetscInt Min3(const PetscInt a[3]) { 68cb32e2e7SValeria Barra return PetscMin(a[0], PetscMin(a[1], a[2])); 69cb32e2e7SValeria Barra } 70cb32e2e7SValeria Barra static void GlobalNodes(const PetscInt p[3], const PetscInt irank[3], 71cb32e2e7SValeria Barra PetscInt degree, const PetscInt melem[3], 72cb32e2e7SValeria Barra PetscInt mnodes[3]) { 73cb32e2e7SValeria Barra for (int d=0; d<3; d++) 74cb32e2e7SValeria Barra mnodes[d] = degree*melem[d] + (irank[d] == p[d]-1); 75cb32e2e7SValeria Barra } 76cb32e2e7SValeria Barra static PetscInt GlobalStart(const PetscInt p[3], const PetscInt irank[3], 77cb32e2e7SValeria Barra PetscInt degree, const PetscInt melem[3]) { 78cb32e2e7SValeria Barra PetscInt start = 0; 79cb32e2e7SValeria Barra // Dumb brute-force is easier to read 80cb32e2e7SValeria Barra for (PetscInt i=0; i<p[0]; i++) { 81cb32e2e7SValeria Barra for (PetscInt j=0; j<p[1]; j++) { 82cb32e2e7SValeria Barra for (PetscInt k=0; k<p[2]; k++) { 83cb32e2e7SValeria Barra PetscInt mnodes[3], ijkrank[] = {i,j,k}; 84cb32e2e7SValeria Barra if (i == irank[0] && j == irank[1] && k == irank[2]) return start; 85cb32e2e7SValeria Barra GlobalNodes(p, ijkrank, degree, melem, mnodes); 86cb32e2e7SValeria Barra start += mnodes[0] * mnodes[1] * mnodes[2]; 87cb32e2e7SValeria Barra } 88cb32e2e7SValeria Barra } 89cb32e2e7SValeria Barra } 90cb32e2e7SValeria Barra return -1; 91cb32e2e7SValeria Barra } 92cb32e2e7SValeria Barra static int CreateRestriction(Ceed ceed, const CeedInt melem[3], 93cb32e2e7SValeria Barra CeedInt P, CeedInt ncomp, 94cb32e2e7SValeria Barra CeedElemRestriction *Erestrict) { 95cb32e2e7SValeria Barra const PetscInt nelem = melem[0]*melem[1]*melem[2]; 96cb32e2e7SValeria Barra PetscInt mnodes[3], *idx, *idxp; 97cb32e2e7SValeria Barra 98cb32e2e7SValeria Barra // Get indicies 99cb32e2e7SValeria Barra for (int d=0; d<3; d++) mnodes[d] = melem[d]*(P-1) + 1; 100cb32e2e7SValeria Barra idxp = idx = malloc(nelem*P*P*P*sizeof idx[0]); 101cb32e2e7SValeria Barra for (CeedInt i=0; i<melem[0]; i++) 102cb32e2e7SValeria Barra for (CeedInt j=0; j<melem[1]; j++) 103cb32e2e7SValeria Barra for (CeedInt k=0; k<melem[2]; k++,idxp += P*P*P) 104cb32e2e7SValeria Barra for (CeedInt ii=0; ii<P; ii++) 105cb32e2e7SValeria Barra for (CeedInt jj=0; jj<P; jj++) 106cb32e2e7SValeria Barra for (CeedInt kk=0; kk<P; kk++) { 107cb32e2e7SValeria Barra if (0) { // This is the C-style (i,j,k) ordering that I prefer 108cb32e2e7SValeria Barra idxp[(ii*P+jj)*P+kk] = (((i*(P-1)+ii)*mnodes[1] 109cb32e2e7SValeria Barra + (j*(P-1)+jj))*mnodes[2] 110cb32e2e7SValeria Barra + (k*(P-1)+kk)); 111cb32e2e7SValeria Barra } else { // (k,j,i) ordering for consistency with MFEM example 112cb32e2e7SValeria Barra idxp[ii+P*(jj+P*kk)] = (((i*(P-1)+ii)*mnodes[1] 113cb32e2e7SValeria Barra + (j*(P-1)+jj))*mnodes[2] 114cb32e2e7SValeria Barra + (k*(P-1)+kk)); 115cb32e2e7SValeria Barra } 116cb32e2e7SValeria Barra } 117cb32e2e7SValeria Barra 118cb32e2e7SValeria Barra // Setup CEED restriction 119cb32e2e7SValeria Barra CeedElemRestrictionCreate(ceed, nelem, P*P*P, mnodes[0]*mnodes[1]*mnodes[2], 120cb32e2e7SValeria Barra ncomp, CEED_MEM_HOST, CEED_OWN_POINTER, idx, 121cb32e2e7SValeria Barra Erestrict); 122cb32e2e7SValeria Barra 123cb32e2e7SValeria Barra PetscFunctionReturn(0); 124cb32e2e7SValeria Barra } 125cb32e2e7SValeria Barra 126cb32e2e7SValeria Barra // Data for PETSc 127cb32e2e7SValeria Barra typedef struct User_ *User; 128cb32e2e7SValeria Barra struct User_ { 129cb32e2e7SValeria Barra MPI_Comm comm; 130cb32e2e7SValeria Barra VecScatter ltog; // Scatter for all entries 131cb32e2e7SValeria Barra VecScatter ltog0; // Skip Dirichlet values 132cb32e2e7SValeria Barra VecScatter gtogD; // global-to-global; only Dirichlet values 133cb32e2e7SValeria Barra Vec Xloc, Yloc; 134cb32e2e7SValeria Barra CeedVector xceed, yceed; 135cb32e2e7SValeria Barra CeedOperator op; 136cb32e2e7SValeria Barra CeedVector qdata; 137cb32e2e7SValeria Barra Ceed ceed; 138cb32e2e7SValeria Barra }; 139cb32e2e7SValeria Barra 140cb32e2e7SValeria Barra // BP Options 141cb32e2e7SValeria Barra typedef enum { 142cb32e2e7SValeria Barra CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2, 143cb32e2e7SValeria Barra CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5 144cb32e2e7SValeria Barra } bpType; 145cb32e2e7SValeria Barra static const char *const bpTypes[] = {"bp1","bp2","bp3","bp4","bp5","bp6", 146cb32e2e7SValeria Barra "bpType","CEED_BP",0 147cb32e2e7SValeria Barra }; 148cb32e2e7SValeria Barra 149cb32e2e7SValeria Barra // BP specific data 150cb32e2e7SValeria Barra typedef struct { 151cb32e2e7SValeria Barra CeedInt ncompu, qdatasize, qextra; 152cb32e2e7SValeria Barra CeedQFunctionUser setupgeo, setuprhs, apply, error; 153cb32e2e7SValeria Barra const char *setupgeofname, *setuprhsfname, *applyfname, *errorfname; 154cb32e2e7SValeria Barra CeedEvalMode inmode, outmode; 155cb32e2e7SValeria Barra CeedQuadMode qmode; 156cb32e2e7SValeria Barra } bpData; 157cb32e2e7SValeria Barra 158cb32e2e7SValeria Barra bpData bpOptions[6] = { 159cb32e2e7SValeria Barra [CEED_BP1] = { 160cb32e2e7SValeria Barra .ncompu = 1, 161cb32e2e7SValeria Barra .qdatasize = 1, 162cb32e2e7SValeria Barra .qextra = 1, 163cb32e2e7SValeria Barra .setupgeo = SetupMassGeo, 164cb32e2e7SValeria Barra .setuprhs = SetupMassRhs, 165cb32e2e7SValeria Barra .apply = Mass, 166cb32e2e7SValeria Barra .error = Error, 167cb32e2e7SValeria Barra .setupgeofname = SetupMassGeo_loc, 168cb32e2e7SValeria Barra .setuprhsfname = SetupMassRhs_loc, 169cb32e2e7SValeria Barra .applyfname = Mass_loc, 170cb32e2e7SValeria Barra .errorfname = Error_loc, 171cb32e2e7SValeria Barra .inmode = CEED_EVAL_INTERP, 172cb32e2e7SValeria Barra .outmode = CEED_EVAL_INTERP, 173cb32e2e7SValeria Barra .qmode = CEED_GAUSS 174cb32e2e7SValeria Barra }, 175cb32e2e7SValeria Barra [CEED_BP2] = { 176cb32e2e7SValeria Barra .ncompu = 3, 177cb32e2e7SValeria Barra .qdatasize = 1, 178cb32e2e7SValeria Barra .qextra = 1, 179cb32e2e7SValeria Barra .setupgeo = SetupMassGeo, 180cb32e2e7SValeria Barra .setuprhs = SetupMassRhs3, 181cb32e2e7SValeria Barra .apply = Mass3, 182cb32e2e7SValeria Barra .error = Error3, 183cb32e2e7SValeria Barra .setupgeofname = SetupMassGeo_loc, 184cb32e2e7SValeria Barra .setuprhsfname = SetupMassRhs3_loc, 185cb32e2e7SValeria Barra .applyfname = Mass3_loc, 186cb32e2e7SValeria Barra .errorfname = Error3_loc, 187cb32e2e7SValeria Barra .inmode = CEED_EVAL_INTERP, 188cb32e2e7SValeria Barra .outmode = CEED_EVAL_INTERP, 189cb32e2e7SValeria Barra .qmode = CEED_GAUSS 190cb32e2e7SValeria Barra }, 191cb32e2e7SValeria Barra [CEED_BP3] = { 192cb32e2e7SValeria Barra .ncompu = 1, 193cb32e2e7SValeria Barra .qdatasize = 6, 194cb32e2e7SValeria Barra .qextra = 1, 195cb32e2e7SValeria Barra .setupgeo = SetupDiffGeo, 196cb32e2e7SValeria Barra .setuprhs = SetupDiffRhs, 197cb32e2e7SValeria Barra .apply = Diff, 198cb32e2e7SValeria Barra .error = Error, 199cb32e2e7SValeria Barra .setupgeofname = SetupDiffGeo_loc, 200cb32e2e7SValeria Barra .setuprhsfname = SetupDiffRhs_loc, 201cb32e2e7SValeria Barra .applyfname = Diff_loc, 202cb32e2e7SValeria Barra .errorfname = Error_loc, 203cb32e2e7SValeria Barra .inmode = CEED_EVAL_GRAD, 204cb32e2e7SValeria Barra .outmode = CEED_EVAL_GRAD, 205cb32e2e7SValeria Barra .qmode = CEED_GAUSS 206cb32e2e7SValeria Barra }, 207cb32e2e7SValeria Barra [CEED_BP4] = { 208cb32e2e7SValeria Barra .ncompu = 3, 209cb32e2e7SValeria Barra .qdatasize = 6, 210cb32e2e7SValeria Barra .qextra = 1, 211cb32e2e7SValeria Barra .setupgeo = SetupDiffGeo, 212cb32e2e7SValeria Barra .setuprhs = SetupDiffRhs3, 213cb32e2e7SValeria Barra .apply = Diff3, 214cb32e2e7SValeria Barra .error = Error3, 215cb32e2e7SValeria Barra .setupgeofname = SetupDiffGeo_loc, 216cb32e2e7SValeria Barra .setuprhsfname = SetupDiffRhs3_loc, 217cb32e2e7SValeria Barra .applyfname = Diff_loc, 218cb32e2e7SValeria Barra .errorfname = Error3_loc, 219cb32e2e7SValeria Barra .inmode = CEED_EVAL_GRAD, 220cb32e2e7SValeria Barra .outmode = CEED_EVAL_GRAD, 221cb32e2e7SValeria Barra .qmode = CEED_GAUSS 222cb32e2e7SValeria Barra }, 223cb32e2e7SValeria Barra [CEED_BP5] = { 224cb32e2e7SValeria Barra .ncompu = 1, 225cb32e2e7SValeria Barra .qdatasize = 6, 226cb32e2e7SValeria Barra .qextra = 0, 227cb32e2e7SValeria Barra .setupgeo = SetupDiffGeo, 228cb32e2e7SValeria Barra .setuprhs = SetupDiffRhs, 229cb32e2e7SValeria Barra .apply = Diff, 230cb32e2e7SValeria Barra .error = Error, 231cb32e2e7SValeria Barra .setupgeofname = SetupDiffGeo_loc, 232cb32e2e7SValeria Barra .setuprhsfname = SetupDiffRhs_loc, 233cb32e2e7SValeria Barra .applyfname = Diff_loc, 234cb32e2e7SValeria Barra .errorfname = Error_loc, 235cb32e2e7SValeria Barra .inmode = CEED_EVAL_GRAD, 236cb32e2e7SValeria Barra .outmode = CEED_EVAL_GRAD, 237cb32e2e7SValeria Barra .qmode = CEED_GAUSS_LOBATTO 238cb32e2e7SValeria Barra }, 239cb32e2e7SValeria Barra [CEED_BP6] = { 240cb32e2e7SValeria Barra .ncompu = 3, 241cb32e2e7SValeria Barra .qdatasize = 6, 242cb32e2e7SValeria Barra .qextra = 0, 243cb32e2e7SValeria Barra .setupgeo = SetupDiffGeo, 244cb32e2e7SValeria Barra .setuprhs = SetupDiffRhs3, 245cb32e2e7SValeria Barra .apply = Diff3, 246cb32e2e7SValeria Barra .error = Error3, 247cb32e2e7SValeria Barra .setupgeofname = SetupDiffGeo_loc, 248cb32e2e7SValeria Barra .setuprhsfname = SetupDiffRhs3_loc, 249cb32e2e7SValeria Barra .applyfname = Diff_loc, 250cb32e2e7SValeria Barra .errorfname = Error3_loc, 251cb32e2e7SValeria Barra .inmode = CEED_EVAL_GRAD, 252cb32e2e7SValeria Barra .outmode = CEED_EVAL_GRAD, 253cb32e2e7SValeria Barra .qmode = CEED_GAUSS_LOBATTO 254cb32e2e7SValeria Barra } 255cb32e2e7SValeria Barra }; 256cb32e2e7SValeria Barra 257cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the mass matrix 258cb32e2e7SValeria Barra static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) { 259cb32e2e7SValeria Barra PetscErrorCode ierr; 260cb32e2e7SValeria Barra User user; 261cb32e2e7SValeria Barra PetscScalar *x, *y; 262cb32e2e7SValeria Barra 263cb32e2e7SValeria Barra PetscFunctionBeginUser; 264cb32e2e7SValeria Barra ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 265cb32e2e7SValeria Barra ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 266cb32e2e7SValeria Barra SCATTER_REVERSE); CHKERRQ(ierr); 267cb32e2e7SValeria Barra ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE); 268cb32e2e7SValeria Barra CHKERRQ(ierr); 269cb32e2e7SValeria Barra ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 270cb32e2e7SValeria Barra 271cb32e2e7SValeria Barra ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 272cb32e2e7SValeria Barra ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 273cb32e2e7SValeria Barra CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 274cb32e2e7SValeria Barra CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y); 275cb32e2e7SValeria Barra 276cb32e2e7SValeria Barra CeedOperatorApply(user->op, user->xceed, user->yceed, 277cb32e2e7SValeria Barra CEED_REQUEST_IMMEDIATE); 278cb32e2e7SValeria Barra ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr); 279cb32e2e7SValeria Barra 280cb32e2e7SValeria Barra ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 281cb32e2e7SValeria Barra ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 282cb32e2e7SValeria Barra 283cb32e2e7SValeria Barra if (Y) { 284cb32e2e7SValeria Barra ierr = VecZeroEntries(Y); CHKERRQ(ierr); 285cb32e2e7SValeria Barra ierr = VecScatterBegin(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 286cb32e2e7SValeria Barra CHKERRQ(ierr); 287cb32e2e7SValeria Barra ierr = VecScatterEnd(user->ltog, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 288cb32e2e7SValeria Barra CHKERRQ(ierr); 289cb32e2e7SValeria Barra } 290cb32e2e7SValeria Barra PetscFunctionReturn(0); 291cb32e2e7SValeria Barra } 292cb32e2e7SValeria Barra 293cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the Laplacian with 294cb32e2e7SValeria Barra // Dirichlet boundary conditions 295cb32e2e7SValeria Barra static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) { 296cb32e2e7SValeria Barra PetscErrorCode ierr; 297cb32e2e7SValeria Barra User user; 298cb32e2e7SValeria Barra PetscScalar *x, *y; 299cb32e2e7SValeria Barra 300cb32e2e7SValeria Barra PetscFunctionBeginUser; 301cb32e2e7SValeria Barra ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 302cb32e2e7SValeria Barra 303cb32e2e7SValeria Barra // Global-to-local 304cb32e2e7SValeria Barra ierr = VecScatterBegin(user->ltog0, X, user->Xloc, INSERT_VALUES, 305cb32e2e7SValeria Barra SCATTER_REVERSE); CHKERRQ(ierr); 306cb32e2e7SValeria Barra ierr = VecScatterEnd(user->ltog0, X, user->Xloc, INSERT_VALUES, 307cb32e2e7SValeria Barra SCATTER_REVERSE); 308cb32e2e7SValeria Barra CHKERRQ(ierr); 309cb32e2e7SValeria Barra ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr); 310cb32e2e7SValeria Barra 311cb32e2e7SValeria Barra // Setup CEED vectors 312cb32e2e7SValeria Barra ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 313cb32e2e7SValeria Barra ierr = VecGetArray(user->Yloc, &y); CHKERRQ(ierr); 314cb32e2e7SValeria Barra CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 315cb32e2e7SValeria Barra CeedVectorSetArray(user->yceed, CEED_MEM_HOST, CEED_USE_POINTER, y); 316cb32e2e7SValeria Barra 317cb32e2e7SValeria Barra // Apply CEED operator 318cb32e2e7SValeria Barra CeedOperatorApply(user->op, user->xceed, user->yceed, 319cb32e2e7SValeria Barra CEED_REQUEST_IMMEDIATE); 320cb32e2e7SValeria Barra ierr = CeedVectorSyncArray(user->yceed, CEED_MEM_HOST); CHKERRQ(ierr); 321cb32e2e7SValeria Barra 322cb32e2e7SValeria Barra // Restore PETSc vectors 323cb32e2e7SValeria Barra ierr = VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 324cb32e2e7SValeria Barra ierr = VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr); 325cb32e2e7SValeria Barra 326cb32e2e7SValeria Barra // Local-to-global 327cb32e2e7SValeria Barra ierr = VecZeroEntries(Y); CHKERRQ(ierr); 328cb32e2e7SValeria Barra ierr = VecScatterBegin(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 329cb32e2e7SValeria Barra CHKERRQ(ierr); 330cb32e2e7SValeria Barra ierr = VecScatterEnd(user->gtogD, X, Y, INSERT_VALUES, SCATTER_FORWARD); 331cb32e2e7SValeria Barra CHKERRQ(ierr); 332cb32e2e7SValeria Barra ierr = VecScatterBegin(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 333cb32e2e7SValeria Barra CHKERRQ(ierr); 334cb32e2e7SValeria Barra ierr = VecScatterEnd(user->ltog0, user->Yloc, Y, ADD_VALUES, SCATTER_FORWARD); 335cb32e2e7SValeria Barra CHKERRQ(ierr); 336cb32e2e7SValeria Barra 337cb32e2e7SValeria Barra PetscFunctionReturn(0); 338cb32e2e7SValeria Barra } 339cb32e2e7SValeria Barra 340cb32e2e7SValeria Barra // This function calculates the error in the final solution 341cb32e2e7SValeria Barra static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X, 342cb32e2e7SValeria Barra CeedVector target, PetscReal *maxerror) { 343cb32e2e7SValeria Barra PetscErrorCode ierr; 344cb32e2e7SValeria Barra PetscScalar *x; 345cb32e2e7SValeria Barra CeedVector collocated_error; 346cb32e2e7SValeria Barra CeedInt length; 347cb32e2e7SValeria Barra 348cb32e2e7SValeria Barra PetscFunctionBeginUser; 349cb32e2e7SValeria Barra CeedVectorGetLength(target, &length); 350cb32e2e7SValeria Barra CeedVectorCreate(user->ceed, length, &collocated_error); 351cb32e2e7SValeria Barra 352cb32e2e7SValeria Barra // Global-to-local 353cb32e2e7SValeria Barra ierr = VecScatterBegin(user->ltog, X, user->Xloc, INSERT_VALUES, 354cb32e2e7SValeria Barra SCATTER_REVERSE); CHKERRQ(ierr); 355cb32e2e7SValeria Barra ierr = VecScatterEnd(user->ltog, X, user->Xloc, INSERT_VALUES, SCATTER_REVERSE); 356cb32e2e7SValeria Barra CHKERRQ(ierr); 357cb32e2e7SValeria Barra 358cb32e2e7SValeria Barra // Setup CEED vector 359cb32e2e7SValeria Barra ierr = VecGetArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 360cb32e2e7SValeria Barra CeedVectorSetArray(user->xceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 361cb32e2e7SValeria Barra 362cb32e2e7SValeria Barra // Apply CEED operator 363cb32e2e7SValeria Barra CeedOperatorApply(op_error, user->xceed, collocated_error, 364cb32e2e7SValeria Barra CEED_REQUEST_IMMEDIATE); 365cb32e2e7SValeria Barra 366cb32e2e7SValeria Barra // Restore PETSc vector 367cb32e2e7SValeria Barra VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x); CHKERRQ(ierr); 368cb32e2e7SValeria Barra 369cb32e2e7SValeria Barra // Reduce max error 370cb32e2e7SValeria Barra *maxerror = 0; 371cb32e2e7SValeria Barra const CeedScalar *e; 372cb32e2e7SValeria Barra CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e); 373cb32e2e7SValeria Barra for (CeedInt i=0; i<length; i++) { 374cb32e2e7SValeria Barra *maxerror = PetscMax(*maxerror, PetscAbsScalar(e[i])); 375cb32e2e7SValeria Barra } 376cb32e2e7SValeria Barra CeedVectorRestoreArrayRead(collocated_error, &e); 377cb32e2e7SValeria Barra ierr = MPI_Allreduce(MPI_IN_PLACE, maxerror, 378cb32e2e7SValeria Barra 1, MPIU_REAL, MPIU_MAX, user->comm); CHKERRQ(ierr); 379cb32e2e7SValeria Barra 380cb32e2e7SValeria Barra // Cleanup 381cb32e2e7SValeria Barra CeedVectorDestroy(&collocated_error); 382cb32e2e7SValeria Barra 383cb32e2e7SValeria Barra PetscFunctionReturn(0); 384cb32e2e7SValeria Barra } 385cb32e2e7SValeria Barra 386cb32e2e7SValeria Barra int main(int argc, char **argv) { 387cb32e2e7SValeria Barra PetscInt ierr; 388cb32e2e7SValeria Barra MPI_Comm comm; 389cb32e2e7SValeria Barra char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 390cb32e2e7SValeria Barra double my_rt_start, my_rt, rt_min, rt_max; 391cb32e2e7SValeria Barra PetscInt degree, qextra, localnodes, localelem, melem[3], mnodes[3], p[3], 392cb32e2e7SValeria Barra irank[3], lnodes[3], lsize, ncompu = 1; 393cb32e2e7SValeria Barra PetscScalar *r; 394cb32e2e7SValeria Barra PetscBool test_mode, benchmark_mode, write_solution; 395cb32e2e7SValeria Barra PetscMPIInt size, rank; 396cb32e2e7SValeria Barra Vec X, Xloc, rhs, rhsloc; 397cb32e2e7SValeria Barra Mat mat; 398cb32e2e7SValeria Barra KSP ksp; 399cb32e2e7SValeria Barra VecScatter ltog, ltog0, gtogD; 400cb32e2e7SValeria Barra User user; 401cb32e2e7SValeria Barra Ceed ceed; 402cb32e2e7SValeria Barra CeedBasis basisx, basisu; 403cb32e2e7SValeria Barra CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, Erestrictui, 404cb32e2e7SValeria Barra Erestrictqdi; 405cb32e2e7SValeria Barra CeedQFunction qf_setupgeo, qf_setuprhs, qf_apply, qf_error; 406cb32e2e7SValeria Barra CeedOperator op_setupgeo, op_setuprhs, op_apply, op_error; 407cb32e2e7SValeria Barra CeedVector xcoord, qdata, rhsceed, target; 408cb32e2e7SValeria Barra CeedInt P, Q; 409cb32e2e7SValeria Barra const CeedInt dim = 3, ncompx = 3; 410cb32e2e7SValeria Barra bpType bpChoice; 411cb32e2e7SValeria Barra 412cb32e2e7SValeria Barra ierr = PetscInitialize(&argc, &argv, NULL, help); 413cb32e2e7SValeria Barra if (ierr) return ierr; 414cb32e2e7SValeria Barra comm = PETSC_COMM_WORLD; 415cb32e2e7SValeria Barra ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr); 416cb32e2e7SValeria Barra bpChoice = CEED_BP1; 417cb32e2e7SValeria Barra ierr = PetscOptionsEnum("-problem", 418cb32e2e7SValeria Barra "CEED benchmark problem to solve", NULL, 419cb32e2e7SValeria Barra bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice, 420cb32e2e7SValeria Barra NULL); CHKERRQ(ierr); 421cb32e2e7SValeria Barra ncompu = bpOptions[bpChoice].ncompu; 422cb32e2e7SValeria Barra test_mode = PETSC_FALSE; 423cb32e2e7SValeria Barra ierr = PetscOptionsBool("-test", 424cb32e2e7SValeria Barra "Testing mode (do not print unless error is large)", 425cb32e2e7SValeria Barra NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 426cb32e2e7SValeria Barra benchmark_mode = PETSC_FALSE; 427cb32e2e7SValeria Barra ierr = PetscOptionsBool("-benchmark", 428cb32e2e7SValeria Barra "Benchmarking mode (prints benchmark statistics)", 429cb32e2e7SValeria Barra NULL, benchmark_mode, &benchmark_mode, NULL); 430cb32e2e7SValeria Barra CHKERRQ(ierr); 431cb32e2e7SValeria Barra write_solution = PETSC_FALSE; 432cb32e2e7SValeria Barra ierr = PetscOptionsBool("-write_solution", 433cb32e2e7SValeria Barra "Write solution for visualization", 434cb32e2e7SValeria Barra NULL, write_solution, &write_solution, NULL); 435cb32e2e7SValeria Barra CHKERRQ(ierr); 436cb32e2e7SValeria Barra degree = test_mode ? 3 : 1; 437cb32e2e7SValeria Barra ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 438cb32e2e7SValeria Barra NULL, degree, °ree, NULL); CHKERRQ(ierr); 439cb32e2e7SValeria Barra qextra = bpOptions[bpChoice].qextra; 440cb32e2e7SValeria Barra ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 441cb32e2e7SValeria Barra NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 442cb32e2e7SValeria Barra ierr = PetscOptionsString("-ceed", "CEED resource specifier", 443cb32e2e7SValeria Barra NULL, ceedresource, ceedresource, 444cb32e2e7SValeria Barra sizeof(ceedresource), NULL); CHKERRQ(ierr); 445cb32e2e7SValeria Barra localnodes = 1000; 446cb32e2e7SValeria Barra ierr = PetscOptionsInt("-local", 447cb32e2e7SValeria Barra "Target number of locally owned nodes per process", 448cb32e2e7SValeria Barra NULL, localnodes, &localnodes, NULL); CHKERRQ(ierr); 449cb32e2e7SValeria Barra ierr = PetscOptionsEnd(); CHKERRQ(ierr); 450cb32e2e7SValeria Barra P = degree + 1; 451cb32e2e7SValeria Barra Q = P + qextra; 452cb32e2e7SValeria Barra 453cb32e2e7SValeria Barra // Determine size of process grid 454cb32e2e7SValeria Barra ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr); 455cb32e2e7SValeria Barra Split3(size, p, false); 456cb32e2e7SValeria Barra 457cb32e2e7SValeria Barra // Find a nicely composite number of elements no less than localnodes 458cb32e2e7SValeria Barra for (localelem = PetscMax(1, localnodes / (degree*degree*degree)); ; 459cb32e2e7SValeria Barra localelem++) { 460cb32e2e7SValeria Barra Split3(localelem, melem, true); 461cb32e2e7SValeria Barra if (Max3(melem) / Min3(melem) <= 2) break; 462cb32e2e7SValeria Barra } 463cb32e2e7SValeria Barra 464cb32e2e7SValeria Barra // Find my location in the process grid 465cb32e2e7SValeria Barra ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr); 466cb32e2e7SValeria Barra for (int d=0,rankleft=rank; d<dim; d++) { 467cb32e2e7SValeria Barra const int pstride[3] = {p[1] *p[2], p[2], 1}; 468cb32e2e7SValeria Barra irank[d] = rankleft / pstride[d]; 469cb32e2e7SValeria Barra rankleft -= irank[d] * pstride[d]; 470cb32e2e7SValeria Barra } 471cb32e2e7SValeria Barra 472cb32e2e7SValeria Barra GlobalNodes(p, irank, degree, melem, mnodes); 473cb32e2e7SValeria Barra 474cb32e2e7SValeria Barra // Setup global vector 475cb32e2e7SValeria Barra ierr = VecCreate(comm, &X); CHKERRQ(ierr); 476cb32e2e7SValeria Barra ierr = VecSetSizes(X, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, PETSC_DECIDE); 477cb32e2e7SValeria Barra CHKERRQ(ierr); 478cb32e2e7SValeria Barra ierr = VecSetUp(X); CHKERRQ(ierr); 479cb32e2e7SValeria Barra 480cb32e2e7SValeria Barra // Set up libCEED 481cb32e2e7SValeria Barra CeedInit(ceedresource, &ceed); 482cb32e2e7SValeria Barra 483cb32e2e7SValeria Barra // Print summary 484cb32e2e7SValeria Barra CeedInt gsize; 485cb32e2e7SValeria Barra ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 486dc7d240cSValeria Barra if (!test_mode) { 487cb32e2e7SValeria Barra const char *usedresource; 488cb32e2e7SValeria Barra CeedGetResource(ceed, &usedresource); 489cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 490cb32e2e7SValeria Barra "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n" 491cb32e2e7SValeria Barra " libCEED:\n" 492cb32e2e7SValeria Barra " libCEED Backend : %s\n" 493cb32e2e7SValeria Barra " Mesh:\n" 494cb32e2e7SValeria Barra " Number of 1D Basis Nodes (p) : %d\n" 495cb32e2e7SValeria Barra " Number of 1D Quadrature Points (q) : %d\n" 496cb32e2e7SValeria Barra " Global nodes : %D\n" 497cb32e2e7SValeria Barra " Process Decomposition : %D %D %D\n" 498cb32e2e7SValeria Barra " Local Elements : %D = %D %D %D\n" 499cb32e2e7SValeria Barra " Owned nodes : %D = %D %D %D\n", 500cb32e2e7SValeria Barra bpChoice+1, usedresource, P, Q, gsize/ncompu, p[0], 501cb32e2e7SValeria Barra p[1], p[2], localelem, melem[0], melem[1], melem[2], 502cb32e2e7SValeria Barra mnodes[0]*mnodes[1]*mnodes[2], mnodes[0], mnodes[1], 503cb32e2e7SValeria Barra mnodes[2]); CHKERRQ(ierr); 504cb32e2e7SValeria Barra } 505cb32e2e7SValeria Barra 506cb32e2e7SValeria Barra { 507cb32e2e7SValeria Barra lsize = 1; 508cb32e2e7SValeria Barra for (int d=0; d<dim; d++) { 509cb32e2e7SValeria Barra lnodes[d] = melem[d]*degree + 1; 510cb32e2e7SValeria Barra lsize *= lnodes[d]; 511cb32e2e7SValeria Barra } 512cb32e2e7SValeria Barra ierr = VecCreate(PETSC_COMM_SELF, &Xloc); CHKERRQ(ierr); 513cb32e2e7SValeria Barra ierr = VecSetSizes(Xloc, lsize*ncompu, PETSC_DECIDE); CHKERRQ(ierr); 514cb32e2e7SValeria Barra ierr = VecSetUp(Xloc); CHKERRQ(ierr); 515cb32e2e7SValeria Barra 516cb32e2e7SValeria Barra // Create local-to-global scatter 517cb32e2e7SValeria Barra PetscInt *ltogind, *ltogind0, *locind, l0count; 518cb32e2e7SValeria Barra IS ltogis, ltogis0, locis; 519cb32e2e7SValeria Barra PetscInt gstart[2][2][2], gmnodes[2][2][2][dim]; 520cb32e2e7SValeria Barra 521cb32e2e7SValeria Barra for (int i=0; i<2; i++) { 522cb32e2e7SValeria Barra for (int j=0; j<2; j++) { 523cb32e2e7SValeria Barra for (int k=0; k<2; k++) { 524cb32e2e7SValeria Barra PetscInt ijkrank[3] = {irank[0]+i, irank[1]+j, irank[2]+k}; 525cb32e2e7SValeria Barra gstart[i][j][k] = GlobalStart(p, ijkrank, degree, melem); 526cb32e2e7SValeria Barra GlobalNodes(p, ijkrank, degree, melem, gmnodes[i][j][k]); 527cb32e2e7SValeria Barra } 528cb32e2e7SValeria Barra } 529cb32e2e7SValeria Barra } 530cb32e2e7SValeria Barra 531cb32e2e7SValeria Barra ierr = PetscMalloc1(lsize, <ogind); CHKERRQ(ierr); 532cb32e2e7SValeria Barra ierr = PetscMalloc1(lsize, <ogind0); CHKERRQ(ierr); 533cb32e2e7SValeria Barra ierr = PetscMalloc1(lsize, &locind); CHKERRQ(ierr); 534cb32e2e7SValeria Barra l0count = 0; 535cb32e2e7SValeria Barra for (PetscInt i=0,ir,ii; ir=i>=mnodes[0], ii=i-ir*mnodes[0], i<lnodes[0]; i++) 536cb32e2e7SValeria Barra for (PetscInt j=0,jr,jj; jr=j>=mnodes[1], jj=j-jr*mnodes[1], j<lnodes[1]; j++) 537cb32e2e7SValeria Barra for (PetscInt k=0,kr,kk; kr=k>=mnodes[2], kk=k-kr*mnodes[2], k<lnodes[2]; k++) { 538cb32e2e7SValeria Barra PetscInt here = (i*lnodes[1]+j)*lnodes[2]+k; 539cb32e2e7SValeria Barra ltogind[here] = 540cb32e2e7SValeria Barra gstart[ir][jr][kr] + (ii*gmnodes[ir][jr][kr][1]+jj)*gmnodes[ir][jr][kr][2]+kk; 541cb32e2e7SValeria Barra if ((irank[0] == 0 && i == 0) 542cb32e2e7SValeria Barra || (irank[1] == 0 && j == 0) 543cb32e2e7SValeria Barra || (irank[2] == 0 && k == 0) 544cb32e2e7SValeria Barra || (irank[0]+1 == p[0] && i+1 == lnodes[0]) 545cb32e2e7SValeria Barra || (irank[1]+1 == p[1] && j+1 == lnodes[1]) 546cb32e2e7SValeria Barra || (irank[2]+1 == p[2] && k+1 == lnodes[2])) 547cb32e2e7SValeria Barra continue; 548cb32e2e7SValeria Barra ltogind0[l0count] = ltogind[here]; 549cb32e2e7SValeria Barra locind[l0count++] = here; 550cb32e2e7SValeria Barra } 551cb32e2e7SValeria Barra ierr = ISCreateBlock(comm, ncompu, lsize, ltogind, PETSC_OWN_POINTER, 552cb32e2e7SValeria Barra <ogis); CHKERRQ(ierr); 553cb32e2e7SValeria Barra ierr = VecScatterCreate(Xloc, NULL, X, ltogis, <og); CHKERRQ(ierr); 554cb32e2e7SValeria Barra CHKERRQ(ierr); 555cb32e2e7SValeria Barra ierr = ISCreateBlock(comm, ncompu, l0count, ltogind0, PETSC_OWN_POINTER, 556cb32e2e7SValeria Barra <ogis0); CHKERRQ(ierr); 557cb32e2e7SValeria Barra ierr = ISCreateBlock(comm, ncompu, l0count, locind, PETSC_OWN_POINTER, 558cb32e2e7SValeria Barra &locis); CHKERRQ(ierr); 559cb32e2e7SValeria Barra ierr = VecScatterCreate(Xloc, locis, X, ltogis0, <og0); CHKERRQ(ierr); 560cb32e2e7SValeria Barra { 561cb32e2e7SValeria Barra // Create global-to-global scatter for Dirichlet values (everything not in 562cb32e2e7SValeria Barra // ltogis0, which is the range of ltog0) 563cb32e2e7SValeria Barra PetscInt xstart, xend, *indD, countD = 0; 564cb32e2e7SValeria Barra IS isD; 565cb32e2e7SValeria Barra const PetscScalar *x; 566cb32e2e7SValeria Barra ierr = VecZeroEntries(Xloc); CHKERRQ(ierr); 567cb32e2e7SValeria Barra ierr = VecSet(X, 1.0); CHKERRQ(ierr); 568cb32e2e7SValeria Barra ierr = VecScatterBegin(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 569cb32e2e7SValeria Barra CHKERRQ(ierr); 570cb32e2e7SValeria Barra ierr = VecScatterEnd(ltog0, Xloc, X, INSERT_VALUES, SCATTER_FORWARD); 571cb32e2e7SValeria Barra CHKERRQ(ierr); 572cb32e2e7SValeria Barra ierr = VecGetOwnershipRange(X, &xstart, &xend); CHKERRQ(ierr); 573cb32e2e7SValeria Barra ierr = PetscMalloc1(xend-xstart, &indD); CHKERRQ(ierr); 574cb32e2e7SValeria Barra ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr); 575cb32e2e7SValeria Barra for (PetscInt i=0; i<xend-xstart; i++) { 576cb32e2e7SValeria Barra if (x[i] == 1.) indD[countD++] = xstart + i; 577cb32e2e7SValeria Barra } 578cb32e2e7SValeria Barra ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr); 579cb32e2e7SValeria Barra ierr = ISCreateGeneral(comm, countD, indD, PETSC_COPY_VALUES, &isD); 580cb32e2e7SValeria Barra CHKERRQ(ierr); 581cb32e2e7SValeria Barra ierr = PetscFree(indD); CHKERRQ(ierr); 582cb32e2e7SValeria Barra ierr = VecScatterCreate(X, isD, X, isD, >ogD); CHKERRQ(ierr); 583cb32e2e7SValeria Barra ierr = ISDestroy(&isD); CHKERRQ(ierr); 584cb32e2e7SValeria Barra } 585cb32e2e7SValeria Barra ierr = ISDestroy(<ogis); CHKERRQ(ierr); 586cb32e2e7SValeria Barra ierr = ISDestroy(<ogis0); CHKERRQ(ierr); 587cb32e2e7SValeria Barra ierr = ISDestroy(&locis); CHKERRQ(ierr); 588cb32e2e7SValeria Barra } 589cb32e2e7SValeria Barra 590cb32e2e7SValeria Barra // CEED bases 591cb32e2e7SValeria Barra CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompu, P, Q, 592cb32e2e7SValeria Barra bpOptions[bpChoice].qmode, &basisu); 593cb32e2e7SValeria Barra CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, Q, 594cb32e2e7SValeria Barra bpOptions[bpChoice].qmode, &basisx); 595cb32e2e7SValeria Barra 596cb32e2e7SValeria Barra // CEED restrictions 597cb32e2e7SValeria Barra CreateRestriction(ceed, melem, P, ncompu, &Erestrictu); 598cb32e2e7SValeria Barra CreateRestriction(ceed, melem, 2, dim, &Erestrictx); 599cb32e2e7SValeria Barra CeedInt nelem = melem[0]*melem[1]*melem[2]; 600cb32e2e7SValeria Barra CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, ncompu, 601cb32e2e7SValeria Barra &Erestrictui); 602cb32e2e7SValeria Barra CeedElemRestrictionCreateIdentity(ceed, nelem, 603cb32e2e7SValeria Barra Q*Q*Q, 604cb32e2e7SValeria Barra nelem*Q*Q*Q, 605cb32e2e7SValeria Barra bpOptions[bpChoice].qdatasize, &Erestrictqdi); 606cb32e2e7SValeria Barra CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q*Q, nelem*Q*Q*Q, 1, 607cb32e2e7SValeria Barra &Erestrictxi); 608cb32e2e7SValeria Barra { 609cb32e2e7SValeria Barra CeedScalar *xloc; 610cb32e2e7SValeria Barra CeedInt shape[3] = {melem[0]+1, melem[1]+1, melem[2]+1}, len = 611cb32e2e7SValeria Barra shape[0]*shape[1]*shape[2]; 612cb32e2e7SValeria Barra xloc = malloc(len*ncompx*sizeof xloc[0]); 613cb32e2e7SValeria Barra for (CeedInt i=0; i<shape[0]; i++) { 614cb32e2e7SValeria Barra for (CeedInt j=0; j<shape[1]; j++) { 615cb32e2e7SValeria Barra for (CeedInt k=0; k<shape[2]; k++) { 616cb32e2e7SValeria Barra xloc[((i*shape[1]+j)*shape[2]+k) + 0*len] = 1.*(irank[0]*melem[0]+i) / 617cb32e2e7SValeria Barra (p[0]*melem[0]); 618cb32e2e7SValeria Barra xloc[((i*shape[1]+j)*shape[2]+k) + 1*len] = 1.*(irank[1]*melem[1]+j) / 619cb32e2e7SValeria Barra (p[1]*melem[1]); 620cb32e2e7SValeria Barra xloc[((i*shape[1]+j)*shape[2]+k) + 2*len] = 1.*(irank[2]*melem[2]+k) / 621cb32e2e7SValeria Barra (p[2]*melem[2]); 622cb32e2e7SValeria Barra } 623cb32e2e7SValeria Barra } 624cb32e2e7SValeria Barra } 625cb32e2e7SValeria Barra CeedVectorCreate(ceed, len*ncompx, &xcoord); 626cb32e2e7SValeria Barra CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_OWN_POINTER, xloc); 627cb32e2e7SValeria Barra } 628cb32e2e7SValeria Barra 629cb32e2e7SValeria Barra // Create the Qfunction that builds the operator quadrature data 630cb32e2e7SValeria Barra CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setupgeo, 631cb32e2e7SValeria Barra bpOptions[bpChoice].setupgeofname, &qf_setupgeo); 632cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*dim, CEED_EVAL_GRAD); 633cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT); 634cb32e2e7SValeria Barra CeedQFunctionAddOutput(qf_setupgeo, "qdata", bpOptions[bpChoice].qdatasize, 635cb32e2e7SValeria Barra CEED_EVAL_NONE); 636cb32e2e7SValeria Barra 637cb32e2e7SValeria Barra // Create the Qfunction that sets up the RHS and true solution 638cb32e2e7SValeria Barra CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].setuprhs, 639cb32e2e7SValeria Barra bpOptions[bpChoice].setuprhsfname, &qf_setuprhs); 640cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_setuprhs, "x", ncompx, CEED_EVAL_INTERP); 641cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_setuprhs, "dx", ncompx*dim, CEED_EVAL_GRAD); 642cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_setuprhs, "weight", 1, CEED_EVAL_WEIGHT); 643cb32e2e7SValeria Barra CeedQFunctionAddOutput(qf_setuprhs, "true_soln", ncompu, CEED_EVAL_NONE); 644cb32e2e7SValeria Barra CeedQFunctionAddOutput(qf_setuprhs, "rhs", ncompu, CEED_EVAL_INTERP); 645cb32e2e7SValeria Barra 646cb32e2e7SValeria Barra // Set up PDE operator 647cb32e2e7SValeria Barra CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].apply, 648cb32e2e7SValeria Barra bpOptions[bpChoice].applyfname, &qf_apply); 649cb32e2e7SValeria Barra // Add inputs and outputs 650cb32e2e7SValeria Barra CeedInt inscale = bpOptions[bpChoice].inmode==CEED_EVAL_GRAD ? 3 : 1; 651cb32e2e7SValeria Barra CeedInt outscale = bpOptions[bpChoice].outmode==CEED_EVAL_GRAD ? 3 : 1; 652cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_apply, "u", ncompu*inscale, 653cb32e2e7SValeria Barra bpOptions[bpChoice].inmode); 654cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_apply, "qdata", bpOptions[bpChoice].qdatasize, 655cb32e2e7SValeria Barra CEED_EVAL_NONE); 656cb32e2e7SValeria Barra CeedQFunctionAddOutput(qf_apply, "v", ncompu*outscale, 657cb32e2e7SValeria Barra bpOptions[bpChoice].outmode); 658cb32e2e7SValeria Barra 659cb32e2e7SValeria Barra // Create the error qfunction 660cb32e2e7SValeria Barra CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error, 661cb32e2e7SValeria Barra bpOptions[bpChoice].errorfname, &qf_error); 662cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP); 663cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE); 664cb32e2e7SValeria Barra CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE); 665cb32e2e7SValeria Barra 666cb32e2e7SValeria Barra // Create the persistent vectors that will be needed in setup 667cb32e2e7SValeria Barra CeedInt nqpts; 668cb32e2e7SValeria Barra CeedBasisGetNumQuadraturePoints(basisu, &nqpts); 669cb32e2e7SValeria Barra CeedVectorCreate(ceed, bpOptions[bpChoice].qdatasize*nelem*nqpts, &qdata); 670cb32e2e7SValeria Barra CeedVectorCreate(ceed, nelem*nqpts*ncompu, &target); 671cb32e2e7SValeria Barra CeedVectorCreate(ceed, lsize*ncompu, &rhsceed); 672cb32e2e7SValeria Barra 673cb32e2e7SValeria Barra // Create the operator that builds the quadrature data for the ceed operator 674*442e7f0bSjeremylt CeedOperatorCreate(ceed, qf_setupgeo, CEED_QFUNCTION_NONE, 675*442e7f0bSjeremylt CEED_QFUNCTION_NONE, &op_setupgeo); 676cb32e2e7SValeria Barra CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, CEED_NOTRANSPOSE, 677cb32e2e7SValeria Barra basisx, CEED_VECTOR_ACTIVE); 678cb32e2e7SValeria Barra CeedOperatorSetField(op_setupgeo, "weight", Erestrictxi, CEED_NOTRANSPOSE, 679cb32e2e7SValeria Barra basisx, CEED_VECTOR_NONE); 680cb32e2e7SValeria Barra CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi, CEED_NOTRANSPOSE, 681cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 682cb32e2e7SValeria Barra 683cb32e2e7SValeria Barra // Create the operator that builds the RHS and true solution 684*442e7f0bSjeremylt CeedOperatorCreate(ceed, qf_setuprhs, CEED_QFUNCTION_NONE, 685*442e7f0bSjeremylt CEED_QFUNCTION_NONE, &op_setuprhs); 686cb32e2e7SValeria Barra CeedOperatorSetField(op_setuprhs, "x", Erestrictx, CEED_NOTRANSPOSE, 687cb32e2e7SValeria Barra basisx, CEED_VECTOR_ACTIVE); 688cb32e2e7SValeria Barra CeedOperatorSetField(op_setuprhs, "dx", Erestrictx, CEED_NOTRANSPOSE, 689cb32e2e7SValeria Barra basisx, CEED_VECTOR_ACTIVE); 690cb32e2e7SValeria Barra CeedOperatorSetField(op_setuprhs, "weight", Erestrictxi, CEED_NOTRANSPOSE, 691cb32e2e7SValeria Barra basisx, CEED_VECTOR_NONE); 692cb32e2e7SValeria Barra CeedOperatorSetField(op_setuprhs, "true_soln", Erestrictui, CEED_NOTRANSPOSE, 693cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, target); 694cb32e2e7SValeria Barra CeedOperatorSetField(op_setuprhs, "rhs", Erestrictu, CEED_TRANSPOSE, 695cb32e2e7SValeria Barra basisu, CEED_VECTOR_ACTIVE); 696cb32e2e7SValeria Barra 697cb32e2e7SValeria Barra // Create the mass or diff operator 698*442e7f0bSjeremylt CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 699*442e7f0bSjeremylt &op_apply); 700cb32e2e7SValeria Barra CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE, 701cb32e2e7SValeria Barra basisu, CEED_VECTOR_ACTIVE); 702cb32e2e7SValeria Barra CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_NOTRANSPOSE, 703cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, qdata); 704cb32e2e7SValeria Barra CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE, 705cb32e2e7SValeria Barra basisu, CEED_VECTOR_ACTIVE); 706cb32e2e7SValeria Barra 707cb32e2e7SValeria Barra // Create the error operator 708*442e7f0bSjeremylt CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 709*442e7f0bSjeremylt &op_error); 710cb32e2e7SValeria Barra CeedOperatorSetField(op_error, "u", Erestrictu, CEED_TRANSPOSE, 711cb32e2e7SValeria Barra basisu, CEED_VECTOR_ACTIVE); 712cb32e2e7SValeria Barra CeedOperatorSetField(op_error, "true_soln", Erestrictui, CEED_NOTRANSPOSE, 713cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, target); 714cb32e2e7SValeria Barra CeedOperatorSetField(op_error, "error", Erestrictui, CEED_NOTRANSPOSE, 715cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 716cb32e2e7SValeria Barra 717cb32e2e7SValeria Barra // Set up Mat 718cb32e2e7SValeria Barra ierr = PetscMalloc1(1, &user); CHKERRQ(ierr); 719cb32e2e7SValeria Barra user->comm = comm; 720cb32e2e7SValeria Barra user->ltog = ltog; 721cb32e2e7SValeria Barra if (bpChoice != CEED_BP1 && bpChoice != CEED_BP2) { 722cb32e2e7SValeria Barra user->ltog0 = ltog0; 723cb32e2e7SValeria Barra user->gtogD = gtogD; 724cb32e2e7SValeria Barra } 725cb32e2e7SValeria Barra user->Xloc = Xloc; 726cb32e2e7SValeria Barra ierr = VecDuplicate(Xloc, &user->Yloc); CHKERRQ(ierr); 727cb32e2e7SValeria Barra CeedVectorCreate(ceed, lsize*ncompu, &user->xceed); 728cb32e2e7SValeria Barra CeedVectorCreate(ceed, lsize*ncompu, &user->yceed); 729cb32e2e7SValeria Barra user->op = op_apply; 730cb32e2e7SValeria Barra user->qdata = qdata; 731cb32e2e7SValeria Barra user->ceed = ceed; 732cb32e2e7SValeria Barra 733cb32e2e7SValeria Barra ierr = MatCreateShell(comm, mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 734cb32e2e7SValeria Barra mnodes[0]*mnodes[1]*mnodes[2]*ncompu, 735cb32e2e7SValeria Barra PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr); 736cb32e2e7SValeria Barra if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) { 737cb32e2e7SValeria Barra ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass); 738cb32e2e7SValeria Barra CHKERRQ(ierr); 739cb32e2e7SValeria Barra } else { 740cb32e2e7SValeria Barra ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff); 741cb32e2e7SValeria Barra CHKERRQ(ierr); 742cb32e2e7SValeria Barra } 743cb32e2e7SValeria Barra ierr = MatCreateVecs(mat, &rhs, NULL); CHKERRQ(ierr); 744cb32e2e7SValeria Barra 745cb32e2e7SValeria Barra // Get RHS vector 746cb32e2e7SValeria Barra ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr); 747cb32e2e7SValeria Barra ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr); 748cb32e2e7SValeria Barra ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr); 749cb32e2e7SValeria Barra CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r); 750cb32e2e7SValeria Barra 751cb32e2e7SValeria Barra // Setup qdata, rhs, and target 752cb32e2e7SValeria Barra CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE); 753cb32e2e7SValeria Barra CeedOperatorApply(op_setuprhs, xcoord, rhsceed, CEED_REQUEST_IMMEDIATE); 754cb32e2e7SValeria Barra ierr = CeedVectorSyncArray(rhsceed, CEED_MEM_HOST); CHKERRQ(ierr); 755cb32e2e7SValeria Barra CeedVectorDestroy(&xcoord); 756cb32e2e7SValeria Barra 757cb32e2e7SValeria Barra // Gather RHS 758cb32e2e7SValeria Barra ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr); 759cb32e2e7SValeria Barra ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 760cb32e2e7SValeria Barra ierr = VecScatterBegin(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 761cb32e2e7SValeria Barra CHKERRQ(ierr); 762cb32e2e7SValeria Barra ierr = VecScatterEnd(ltog, rhsloc, rhs, ADD_VALUES, SCATTER_FORWARD); 763cb32e2e7SValeria Barra CHKERRQ(ierr); 764cb32e2e7SValeria Barra CeedVectorDestroy(&rhsceed); 765cb32e2e7SValeria Barra 766cb32e2e7SValeria Barra ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 767cb32e2e7SValeria Barra { 768cb32e2e7SValeria Barra PC pc; 769cb32e2e7SValeria Barra ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 770cb32e2e7SValeria Barra if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) { 771cb32e2e7SValeria Barra ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 772cb32e2e7SValeria Barra ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 773cb32e2e7SValeria Barra } else { 774cb32e2e7SValeria Barra ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 775cb32e2e7SValeria Barra } 776cb32e2e7SValeria Barra ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 777cb32e2e7SValeria Barra ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 778cb32e2e7SValeria Barra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 779cb32e2e7SValeria Barra PETSC_DEFAULT); CHKERRQ(ierr); 780cb32e2e7SValeria Barra } 781cb32e2e7SValeria Barra ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 782cb32e2e7SValeria Barra ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr); 783cb32e2e7SValeria Barra // First run, if benchmarking 784cb32e2e7SValeria Barra if (benchmark_mode) { 785cb32e2e7SValeria Barra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 786cb32e2e7SValeria Barra CHKERRQ(ierr); 787cb32e2e7SValeria Barra my_rt_start = MPI_Wtime(); 788cb32e2e7SValeria Barra ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 789cb32e2e7SValeria Barra my_rt = MPI_Wtime() - my_rt_start; 790cb32e2e7SValeria Barra ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 791cb32e2e7SValeria Barra CHKERRQ(ierr); 792cb32e2e7SValeria Barra // Set maxits based on first iteration timing 793cb32e2e7SValeria Barra if (my_rt > 0.02) { 794cb32e2e7SValeria Barra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5); 795cb32e2e7SValeria Barra CHKERRQ(ierr); 796cb32e2e7SValeria Barra } else { 797cb32e2e7SValeria Barra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20); 798cb32e2e7SValeria Barra CHKERRQ(ierr); 799cb32e2e7SValeria Barra } 800cb32e2e7SValeria Barra } 801cb32e2e7SValeria Barra // Timed solve 802cb32e2e7SValeria Barra ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 803cb32e2e7SValeria Barra my_rt_start = MPI_Wtime(); 804cb32e2e7SValeria Barra ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 805cb32e2e7SValeria Barra my_rt = MPI_Wtime() - my_rt_start; 806cb32e2e7SValeria Barra { 807cb32e2e7SValeria Barra KSPType ksptype; 808cb32e2e7SValeria Barra KSPConvergedReason reason; 809cb32e2e7SValeria Barra PetscReal rnorm; 810cb32e2e7SValeria Barra PetscInt its; 811cb32e2e7SValeria Barra ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr); 812cb32e2e7SValeria Barra ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 813cb32e2e7SValeria Barra ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 814cb32e2e7SValeria Barra ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 815cb32e2e7SValeria Barra if (!test_mode || reason < 0 || rnorm > 1e-8) { 816dc7d240cSValeria Barra ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 817dc7d240cSValeria Barra CHKERRQ(ierr); 818dc7d240cSValeria Barra ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 819dc7d240cSValeria Barra CHKERRQ(ierr); 820cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 821cb32e2e7SValeria Barra " KSP:\n" 822cb32e2e7SValeria Barra " KSP Type : %s\n" 823cb32e2e7SValeria Barra " KSP Convergence : %s\n" 824cb32e2e7SValeria Barra " Total KSP Iterations : %D\n" 825cb32e2e7SValeria Barra " Final rnorm : %e\n", 826cb32e2e7SValeria Barra ksptype, KSPConvergedReasons[reason], its, 827cb32e2e7SValeria Barra (double)rnorm); CHKERRQ(ierr); 828dc7d240cSValeria Barra ierr = PetscPrintf(comm, 829dc7d240cSValeria Barra " Performance:\n" 830dc7d240cSValeria Barra " CG Solve Time : %g (%g) sec\n" 831dc7d240cSValeria Barra " DoFs/Sec in CG : %g (%g) million\n", 832dc7d240cSValeria Barra rt_max, rt_min, 1e-6*gsize*its/rt_max, 833dc7d240cSValeria Barra 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 834cb32e2e7SValeria Barra } 835cb32e2e7SValeria Barra if (benchmark_mode && (!test_mode)) { 836cb32e2e7SValeria Barra ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 837cb32e2e7SValeria Barra CHKERRQ(ierr); 838cb32e2e7SValeria Barra ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 839cb32e2e7SValeria Barra CHKERRQ(ierr); 840cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 841cb32e2e7SValeria Barra " Performance:\n" 842cb32e2e7SValeria Barra " CG Solve Time : %g (%g) sec\n" 843cb32e2e7SValeria Barra " DoFs/Sec in CG : %g (%g) million\n", 844cb32e2e7SValeria Barra rt_max, rt_min, 1e-6*gsize*its/rt_max, 845cb32e2e7SValeria Barra 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 846cb32e2e7SValeria Barra } 847cb32e2e7SValeria Barra } 848cb32e2e7SValeria Barra 849cb32e2e7SValeria Barra { 850cb32e2e7SValeria Barra PetscReal maxerror; 851cb32e2e7SValeria Barra ierr = ComputeErrorMax(user, op_error, X, target, &maxerror); CHKERRQ(ierr); 852cb32e2e7SValeria Barra PetscReal tol = (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) ? 5e-3 : 5e-2; 853cb32e2e7SValeria Barra if (!test_mode || maxerror > tol) { 854cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 855cb32e2e7SValeria Barra " Pointwise Error (max) : %e\n", 856cb32e2e7SValeria Barra (double)maxerror); CHKERRQ(ierr); 857cb32e2e7SValeria Barra } 858cb32e2e7SValeria Barra } 859cb32e2e7SValeria Barra 860cb32e2e7SValeria Barra if (write_solution) { 861cb32e2e7SValeria Barra PetscViewer vtkviewersoln; 862cb32e2e7SValeria Barra 863cb32e2e7SValeria Barra ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr); 864cb32e2e7SValeria Barra ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr); 865cb32e2e7SValeria Barra ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr); 866cb32e2e7SValeria Barra ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr); 867cb32e2e7SValeria Barra ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr); 868cb32e2e7SValeria Barra } 869cb32e2e7SValeria Barra 870cb32e2e7SValeria Barra ierr = VecDestroy(&rhs); CHKERRQ(ierr); 871cb32e2e7SValeria Barra ierr = VecDestroy(&rhsloc); CHKERRQ(ierr); 872cb32e2e7SValeria Barra ierr = VecDestroy(&X); CHKERRQ(ierr); 873cb32e2e7SValeria Barra ierr = VecDestroy(&user->Xloc); CHKERRQ(ierr); 874cb32e2e7SValeria Barra ierr = VecDestroy(&user->Yloc); CHKERRQ(ierr); 875cb32e2e7SValeria Barra ierr = VecScatterDestroy(<og); CHKERRQ(ierr); 876cb32e2e7SValeria Barra ierr = VecScatterDestroy(<og0); CHKERRQ(ierr); 877cb32e2e7SValeria Barra ierr = VecScatterDestroy(>ogD); CHKERRQ(ierr); 878cb32e2e7SValeria Barra ierr = MatDestroy(&mat); CHKERRQ(ierr); 879cb32e2e7SValeria Barra ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 880cb32e2e7SValeria Barra 881cb32e2e7SValeria Barra CeedVectorDestroy(&user->xceed); 882cb32e2e7SValeria Barra CeedVectorDestroy(&user->yceed); 883cb32e2e7SValeria Barra CeedVectorDestroy(&user->qdata); 884cb32e2e7SValeria Barra CeedVectorDestroy(&target); 885cb32e2e7SValeria Barra CeedOperatorDestroy(&op_setupgeo); 886cb32e2e7SValeria Barra CeedOperatorDestroy(&op_setuprhs); 887cb32e2e7SValeria Barra CeedOperatorDestroy(&op_apply); 888cb32e2e7SValeria Barra CeedOperatorDestroy(&op_error); 889cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictu); 890cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictx); 891cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictui); 892cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictxi); 893cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictqdi); 894cb32e2e7SValeria Barra CeedQFunctionDestroy(&qf_setupgeo); 895cb32e2e7SValeria Barra CeedQFunctionDestroy(&qf_setuprhs); 896cb32e2e7SValeria Barra CeedQFunctionDestroy(&qf_apply); 897cb32e2e7SValeria Barra CeedQFunctionDestroy(&qf_error); 898cb32e2e7SValeria Barra CeedBasisDestroy(&basisu); 899cb32e2e7SValeria Barra CeedBasisDestroy(&basisx); 900cb32e2e7SValeria Barra CeedDestroy(&ceed); 901cb32e2e7SValeria Barra ierr = PetscFree(user); CHKERRQ(ierr); 902cb32e2e7SValeria Barra return PetscFinalize(); 903cb32e2e7SValeria Barra } 904