1*cb32e2e7SValeria Barra // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2*cb32e2e7SValeria Barra // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3*cb32e2e7SValeria Barra // reserved. See files LICENSE and NOTICE for details. 4*cb32e2e7SValeria Barra // 5*cb32e2e7SValeria Barra // This file is part of CEED, a collection of benchmarks, miniapps, software 6*cb32e2e7SValeria Barra // libraries and APIs for efficient high-order finite element and spectral 7*cb32e2e7SValeria Barra // element discretizations for exascale applications. For more information and 8*cb32e2e7SValeria Barra // source code availability see http://github.com/ceed. 9*cb32e2e7SValeria Barra // 10*cb32e2e7SValeria Barra // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11*cb32e2e7SValeria Barra // a collaborative effort of two U.S. Department of Energy organizations (Office 12*cb32e2e7SValeria Barra // of Science and the National Nuclear Security Administration) responsible for 13*cb32e2e7SValeria Barra // the planning and preparation of a capable exascale ecosystem, including 14*cb32e2e7SValeria Barra // software, applications, hardware, advanced system engineering and early 15*cb32e2e7SValeria Barra // testbed platforms, in support of the nation's exascale computing imperative. 16*cb32e2e7SValeria Barra 17*cb32e2e7SValeria Barra // libCEED + PETSc Example: Surface Area 18*cb32e2e7SValeria Barra // 19*cb32e2e7SValeria Barra // This example demonstrates a simple usage of libCEED with PETSc to calculate 20*cb32e2e7SValeria Barra // the surface area of a simple closed surface, such as the one of a cube 21*cb32e2e7SValeria Barra // via the mass operator. 22*cb32e2e7SValeria Barra // 23*cb32e2e7SValeria Barra // The code uses higher level communication protocols in DMPlex. 24*cb32e2e7SValeria Barra // 25*cb32e2e7SValeria Barra // Build with: 26*cb32e2e7SValeria Barra // 27*cb32e2e7SValeria Barra // make area [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 28*cb32e2e7SValeria Barra // 29*cb32e2e7SValeria Barra // Sample runs: 30*cb32e2e7SValeria Barra // Sequential: 31*cb32e2e7SValeria Barra // 32*cb32e2e7SValeria Barra // area -petscspace_degree 3 33*cb32e2e7SValeria Barra // 34*cb32e2e7SValeria Barra // In parallel: 35*cb32e2e7SValeria Barra // 36*cb32e2e7SValeria Barra // mpiexec -n 4 area -petscspace_degree 3 37*cb32e2e7SValeria Barra // 38*cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -petscspace_degree 3 39*cb32e2e7SValeria Barra 40*cb32e2e7SValeria Barra /// @file 41*cb32e2e7SValeria Barra /// libCEED example using the mass operator to compute surface area using PETSc with DMPlex 42*cb32e2e7SValeria Barra static const char help[] = "Compute surface area of a cube using DMPlex in PETSc\n"; 43*cb32e2e7SValeria Barra 44*cb32e2e7SValeria Barra #include <string.h> 45*cb32e2e7SValeria Barra #include <petscdmplex.h> 46*cb32e2e7SValeria Barra #include <ceed.h> 47*cb32e2e7SValeria Barra #include "qfunctions/area/area.h" 48*cb32e2e7SValeria Barra 49*cb32e2e7SValeria Barra // Auxiliary function to define CEED restrictions from DMPlex data 50*cb32e2e7SValeria Barra static int CreateRestrictionPlex(Ceed ceed, CeedInt P, CeedInt ncomp, 51*cb32e2e7SValeria Barra CeedElemRestriction *Erestrict, DM dm) { 52*cb32e2e7SValeria Barra PetscInt ierr; 53*cb32e2e7SValeria Barra PetscInt c, cStart, cEnd, nelem, nnodes, *erestrict, eoffset; 54*cb32e2e7SValeria Barra PetscSection section; 55*cb32e2e7SValeria Barra Vec Uloc; 56*cb32e2e7SValeria Barra 57*cb32e2e7SValeria Barra PetscFunctionBegin; 58*cb32e2e7SValeria Barra 59*cb32e2e7SValeria Barra // Get Nelem 60*cb32e2e7SValeria Barra ierr = DMGetSection(dm, §ion); CHKERRQ(ierr); 61*cb32e2e7SValeria Barra ierr = DMPlexGetHeightStratum(dm, 0, &cStart,& cEnd); CHKERRQ(ierr); 62*cb32e2e7SValeria Barra nelem = cEnd - cStart; 63*cb32e2e7SValeria Barra 64*cb32e2e7SValeria Barra // Get indices 65*cb32e2e7SValeria Barra ierr = PetscMalloc1(nelem*P*P, &erestrict); CHKERRQ(ierr); 66*cb32e2e7SValeria Barra for (c=cStart, eoffset = 0; c<cEnd; c++) { 67*cb32e2e7SValeria Barra PetscInt numindices, *indices, i; 68*cb32e2e7SValeria Barra ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices, 69*cb32e2e7SValeria Barra &indices, NULL); CHKERRQ(ierr); 70*cb32e2e7SValeria Barra for (i=0; i<numindices; i+=ncomp) { 71*cb32e2e7SValeria Barra for (PetscInt j=0; j<ncomp; j++) { 72*cb32e2e7SValeria Barra if (indices[i+j] != indices[i] + (PetscInt)(copysign(j, indices[i]))) 73*cb32e2e7SValeria Barra SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 74*cb32e2e7SValeria Barra "Cell %D closure indices not interlaced", c); 75*cb32e2e7SValeria Barra } 76*cb32e2e7SValeria Barra // NO BC on closed surfaces 77*cb32e2e7SValeria Barra PetscInt loc = indices[i]; 78*cb32e2e7SValeria Barra erestrict[eoffset++] = loc/ncomp; 79*cb32e2e7SValeria Barra } 80*cb32e2e7SValeria Barra ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices, 81*cb32e2e7SValeria Barra &indices, NULL); CHKERRQ(ierr); 82*cb32e2e7SValeria Barra } 83*cb32e2e7SValeria Barra 84*cb32e2e7SValeria Barra // Setup CEED restriction 85*cb32e2e7SValeria Barra ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 86*cb32e2e7SValeria Barra ierr = VecGetLocalSize(Uloc, &nnodes); CHKERRQ(ierr); 87*cb32e2e7SValeria Barra 88*cb32e2e7SValeria Barra ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 89*cb32e2e7SValeria Barra CeedElemRestrictionCreate(ceed, nelem, P*P, nnodes/ncomp, ncomp, 90*cb32e2e7SValeria Barra CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, 91*cb32e2e7SValeria Barra Erestrict); 92*cb32e2e7SValeria Barra ierr = PetscFree(erestrict); CHKERRQ(ierr); 93*cb32e2e7SValeria Barra 94*cb32e2e7SValeria Barra PetscFunctionReturn(0); 95*cb32e2e7SValeria Barra } 96*cb32e2e7SValeria Barra 97*cb32e2e7SValeria Barra int main(int argc, char **argv) { 98*cb32e2e7SValeria Barra PetscInt ierr; 99*cb32e2e7SValeria Barra MPI_Comm comm; 100*cb32e2e7SValeria Barra char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self", 101*cb32e2e7SValeria Barra filename[PETSC_MAX_PATH_LEN]; 102*cb32e2e7SValeria Barra PetscInt lsize, gsize, xlsize, 103*cb32e2e7SValeria Barra qextra = 1, // default number of extra quadrature points 104*cb32e2e7SValeria Barra ncompx = 3, // number of components of 3D physical coordinates 105*cb32e2e7SValeria Barra ncompu = 1, // dimension of field to which apply mass operator 106*cb32e2e7SValeria Barra topodim = 2, // topological dimension of manifold 107*cb32e2e7SValeria Barra degree = 3; // default degree for finite element bases 108*cb32e2e7SValeria Barra PetscBool read_mesh = PETSC_FALSE, 109*cb32e2e7SValeria Barra test_mode = PETSC_FALSE; 110*cb32e2e7SValeria Barra PetscSpace sp; 111*cb32e2e7SValeria Barra PetscFE fe; 112*cb32e2e7SValeria Barra Vec X, Xloc, V, Vloc; 113*cb32e2e7SValeria Barra DM dm, dmcoord; 114*cb32e2e7SValeria Barra Ceed ceed; 115*cb32e2e7SValeria Barra CeedInt P, Q; 116*cb32e2e7SValeria Barra 117*cb32e2e7SValeria Barra ierr = PetscInitialize(&argc, &argv, NULL, help); 118*cb32e2e7SValeria Barra if (ierr) return ierr; 119*cb32e2e7SValeria Barra comm = PETSC_COMM_WORLD; 120*cb32e2e7SValeria Barra 121*cb32e2e7SValeria Barra // Read CL options 122*cb32e2e7SValeria Barra ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc", NULL); 123*cb32e2e7SValeria Barra CHKERRQ(ierr); 124*cb32e2e7SValeria Barra ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 125*cb32e2e7SValeria Barra NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 126*cb32e2e7SValeria Barra ierr = PetscOptionsString("-ceed", "CEED resource specifier", 127*cb32e2e7SValeria Barra NULL, ceedresource, ceedresource, 128*cb32e2e7SValeria Barra sizeof(ceedresource), NULL); CHKERRQ(ierr); 129*cb32e2e7SValeria Barra ierr = PetscOptionsBool("-test", 130*cb32e2e7SValeria Barra "Testing mode (do not print unless error is large)", 131*cb32e2e7SValeria Barra NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 132*cb32e2e7SValeria Barra ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 133*cb32e2e7SValeria Barra filename, filename, sizeof(filename), &read_mesh); 134*cb32e2e7SValeria Barra CHKERRQ(ierr); 135*cb32e2e7SValeria Barra ierr = PetscOptionsEnd(); CHKERRQ(ierr); 136*cb32e2e7SValeria Barra 137*cb32e2e7SValeria Barra // Setup DM 138*cb32e2e7SValeria Barra PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube 139*cb32e2e7SValeria Barra if (read_mesh) { 140*cb32e2e7SValeria Barra ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm); 141*cb32e2e7SValeria Barra CHKERRQ(ierr); 142*cb32e2e7SValeria Barra } else { 143*cb32e2e7SValeria Barra // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box. 144*cb32e2e7SValeria Barra PetscBool simplex = PETSC_FALSE; 145*cb32e2e7SValeria Barra ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, &dm); 146*cb32e2e7SValeria Barra CHKERRQ(ierr); 147*cb32e2e7SValeria Barra // Set the object name 148*cb32e2e7SValeria Barra ierr = PetscObjectSetName((PetscObject) dm, "Cube"); CHKERRQ(ierr); 149*cb32e2e7SValeria Barra // Distribute mesh over processes 150*cb32e2e7SValeria Barra { 151*cb32e2e7SValeria Barra DM dmDist = NULL; 152*cb32e2e7SValeria Barra PetscPartitioner part; 153*cb32e2e7SValeria Barra 154*cb32e2e7SValeria Barra ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 155*cb32e2e7SValeria Barra ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 156*cb32e2e7SValeria Barra ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 157*cb32e2e7SValeria Barra if (dmDist) { 158*cb32e2e7SValeria Barra ierr = DMDestroy(&dm); CHKERRQ(ierr); 159*cb32e2e7SValeria Barra dm = dmDist; 160*cb32e2e7SValeria Barra } 161*cb32e2e7SValeria Barra } 162*cb32e2e7SValeria Barra // View DMPlex via runtime option 163*cb32e2e7SValeria Barra ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 164*cb32e2e7SValeria Barra } 165*cb32e2e7SValeria Barra 166*cb32e2e7SValeria Barra // Create FE 167*cb32e2e7SValeria Barra ierr = PetscFECreateDefault(PETSC_COMM_SELF, topodim, ncompu, PETSC_FALSE, NULL, 168*cb32e2e7SValeria Barra PETSC_DETERMINE, &fe); 169*cb32e2e7SValeria Barra CHKERRQ(ierr); 170*cb32e2e7SValeria Barra ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 171*cb32e2e7SValeria Barra ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 172*cb32e2e7SValeria Barra ierr = DMCreateDS(dm); CHKERRQ(ierr); 173*cb32e2e7SValeria Barra ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 174*cb32e2e7SValeria Barra CHKERRQ(ierr); 175*cb32e2e7SValeria Barra 176*cb32e2e7SValeria Barra // Get basis space degree 177*cb32e2e7SValeria Barra ierr = PetscFEGetBasisSpace(fe, &sp); CHKERRQ(ierr); 178*cb32e2e7SValeria Barra ierr = PetscSpaceGetDegree(sp, °ree, NULL); CHKERRQ(ierr); 179*cb32e2e7SValeria Barra ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 180*cb32e2e7SValeria Barra if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 181*cb32e2e7SValeria Barra "-petscspace_degree %D must be at least 1", degree); 182*cb32e2e7SValeria Barra 183*cb32e2e7SValeria Barra // Create vectors 184*cb32e2e7SValeria Barra ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr); 185*cb32e2e7SValeria Barra ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr); 186*cb32e2e7SValeria Barra ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 187*cb32e2e7SValeria Barra ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr); 188*cb32e2e7SValeria Barra ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr); 189*cb32e2e7SValeria Barra ierr = VecDuplicate(X, &V); CHKERRQ(ierr); 190*cb32e2e7SValeria Barra ierr = VecDuplicate(Xloc, &Vloc); CHKERRQ(ierr); 191*cb32e2e7SValeria Barra 192*cb32e2e7SValeria Barra // Set up libCEED 193*cb32e2e7SValeria Barra CeedInit(ceedresource, &ceed); 194*cb32e2e7SValeria Barra 195*cb32e2e7SValeria Barra // Print summary 196*cb32e2e7SValeria Barra P = degree + 1; 197*cb32e2e7SValeria Barra Q = P + qextra; 198*cb32e2e7SValeria Barra const char *usedresource; 199*cb32e2e7SValeria Barra CeedGetResource(ceed, &usedresource); 200*cb32e2e7SValeria Barra if (!test_mode) { 201*cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 202*cb32e2e7SValeria Barra "\n-- libCEED + PETSc Surface Area problem --\n" 203*cb32e2e7SValeria Barra " libCEED:\n" 204*cb32e2e7SValeria Barra " libCEED Backend : %s\n" 205*cb32e2e7SValeria Barra " Mesh:\n" 206*cb32e2e7SValeria Barra " Number of 1D Basis Nodes (p) : %d\n" 207*cb32e2e7SValeria Barra " Number of 1D Quadrature Points (q) : %d\n" 208*cb32e2e7SValeria Barra " Global nodes : %D\n", 209*cb32e2e7SValeria Barra usedresource, P, Q, gsize/ncompu); 210*cb32e2e7SValeria Barra CHKERRQ(ierr); 211*cb32e2e7SValeria Barra } 212*cb32e2e7SValeria Barra 213*cb32e2e7SValeria Barra // Setup libCEED's objects 214*cb32e2e7SValeria Barra // Create CEED operators and API objects they need 215*cb32e2e7SValeria Barra CeedOperator op_setupgeo, op_apply; 216*cb32e2e7SValeria Barra CeedQFunction qf_setupgeo, qf_apply; 217*cb32e2e7SValeria Barra CeedBasis basisx, basisu; 218*cb32e2e7SValeria Barra CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, 219*cb32e2e7SValeria Barra Erestrictqdi; 220*cb32e2e7SValeria Barra 221*cb32e2e7SValeria Barra // Create bases 222*cb32e2e7SValeria Barra CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompu, P, Q, 223*cb32e2e7SValeria Barra CEED_GAUSS, &basisu); 224*cb32e2e7SValeria Barra CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompx, 2, Q, 225*cb32e2e7SValeria Barra CEED_GAUSS, &basisx); 226*cb32e2e7SValeria Barra 227*cb32e2e7SValeria Barra // CEED restrictions 228*cb32e2e7SValeria Barra ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 229*cb32e2e7SValeria Barra ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 230*cb32e2e7SValeria Barra CHKERRQ(ierr); 231*cb32e2e7SValeria Barra 232*cb32e2e7SValeria Barra CreateRestrictionPlex(ceed, 2, ncompx, &Erestrictx, dmcoord); CHKERRQ(ierr); 233*cb32e2e7SValeria Barra CreateRestrictionPlex(ceed, P, ncompu, &Erestrictu, dm); CHKERRQ(ierr); 234*cb32e2e7SValeria Barra 235*cb32e2e7SValeria Barra CeedInt cStart, cEnd; 236*cb32e2e7SValeria Barra ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 237*cb32e2e7SValeria Barra const CeedInt nelem = cEnd - cStart; 238*cb32e2e7SValeria Barra 239*cb32e2e7SValeria Barra // CEED identity restrictions 240*cb32e2e7SValeria Barra const CeedInt qdatasize = 1; 241*cb32e2e7SValeria Barra CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q, nelem*Q*Q, 242*cb32e2e7SValeria Barra qdatasize, &Erestrictqdi); 243*cb32e2e7SValeria Barra CeedElemRestrictionCreateIdentity(ceed, nelem, Q*Q, nelem*Q*Q, 1, 244*cb32e2e7SValeria Barra &Erestrictxi); 245*cb32e2e7SValeria Barra 246*cb32e2e7SValeria Barra // Element coordinates 247*cb32e2e7SValeria Barra Vec coords; 248*cb32e2e7SValeria Barra const PetscScalar *coordArray; 249*cb32e2e7SValeria Barra PetscSection section; 250*cb32e2e7SValeria Barra ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr); 251*cb32e2e7SValeria Barra ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr); 252*cb32e2e7SValeria Barra ierr = DMGetSection(dmcoord, §ion); CHKERRQ(ierr); 253*cb32e2e7SValeria Barra 254*cb32e2e7SValeria Barra CeedVector xcoord; 255*cb32e2e7SValeria Barra CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL); 256*cb32e2e7SValeria Barra CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, 257*cb32e2e7SValeria Barra (PetscScalar *)coordArray); 258*cb32e2e7SValeria Barra ierr = VecRestoreArrayRead(coords, &coordArray); 259*cb32e2e7SValeria Barra 260*cb32e2e7SValeria Barra // Create the vectors that will be needed in setup and apply 261*cb32e2e7SValeria Barra CeedVector uceed, vceed, qdata; 262*cb32e2e7SValeria Barra CeedInt nqpts; 263*cb32e2e7SValeria Barra CeedBasisGetNumQuadraturePoints(basisu, &nqpts); 264*cb32e2e7SValeria Barra CeedVectorCreate(ceed, qdatasize*nelem*nqpts, &qdata); 265*cb32e2e7SValeria Barra CeedVectorCreate(ceed, xlsize, &uceed); 266*cb32e2e7SValeria Barra CeedVectorCreate(ceed, xlsize, &vceed); 267*cb32e2e7SValeria Barra 268*cb32e2e7SValeria Barra /* Create the Q-function that builds the operator for the geomteric factors 269*cb32e2e7SValeria Barra (i.e., the quadrature data) */ 270*cb32e2e7SValeria Barra CeedQFunctionCreateInterior(ceed, 1, SetupMassGeo, 271*cb32e2e7SValeria Barra SetupMassGeo_loc, &qf_setupgeo); 272*cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*topodim, CEED_EVAL_GRAD); 273*cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT); 274*cb32e2e7SValeria Barra CeedQFunctionAddOutput(qf_setupgeo, "qdata", qdatasize, CEED_EVAL_NONE); 275*cb32e2e7SValeria Barra 276*cb32e2e7SValeria Barra // Set up the mass operator 277*cb32e2e7SValeria Barra CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_apply); 278*cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_apply, "u", ncompu, CEED_EVAL_INTERP); 279*cb32e2e7SValeria Barra CeedQFunctionAddInput(qf_apply, "qdata", qdatasize, CEED_EVAL_NONE); 280*cb32e2e7SValeria Barra CeedQFunctionAddOutput(qf_apply, "v", ncompu, CEED_EVAL_INTERP); 281*cb32e2e7SValeria Barra 282*cb32e2e7SValeria Barra // Create the operator that builds the quadrature data for the operator 283*cb32e2e7SValeria Barra CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo); 284*cb32e2e7SValeria Barra CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, CEED_TRANSPOSE, 285*cb32e2e7SValeria Barra basisx, CEED_VECTOR_ACTIVE); 286*cb32e2e7SValeria Barra CeedOperatorSetField(op_setupgeo, "weight", Erestrictxi, CEED_NOTRANSPOSE, 287*cb32e2e7SValeria Barra basisx, CEED_VECTOR_NONE); 288*cb32e2e7SValeria Barra CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi, CEED_NOTRANSPOSE, 289*cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 290*cb32e2e7SValeria Barra 291*cb32e2e7SValeria Barra // Create the mass operator 292*cb32e2e7SValeria Barra CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply); 293*cb32e2e7SValeria Barra CeedOperatorSetField(op_apply, "u", Erestrictu, CEED_TRANSPOSE, 294*cb32e2e7SValeria Barra basisu, CEED_VECTOR_ACTIVE); 295*cb32e2e7SValeria Barra CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_NOTRANSPOSE, 296*cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, qdata); 297*cb32e2e7SValeria Barra CeedOperatorSetField(op_apply, "v", Erestrictu, CEED_TRANSPOSE, 298*cb32e2e7SValeria Barra basisu, CEED_VECTOR_ACTIVE); 299*cb32e2e7SValeria Barra 300*cb32e2e7SValeria Barra // Compute the quadrature data for the mass operator 301*cb32e2e7SValeria Barra CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE); 302*cb32e2e7SValeria Barra 303*cb32e2e7SValeria Barra PetscScalar *v; 304*cb32e2e7SValeria Barra ierr = VecZeroEntries(Vloc); CHKERRQ(ierr); 305*cb32e2e7SValeria Barra ierr = VecGetArray(Vloc, &v); 306*cb32e2e7SValeria Barra CeedVectorSetArray(vceed, CEED_MEM_HOST, CEED_USE_POINTER, v); 307*cb32e2e7SValeria Barra 308*cb32e2e7SValeria Barra // Compute the mesh volume using the mass operator: vol = 1^T \cdot M \cdot 1 309*cb32e2e7SValeria Barra if (!test_mode) { 310*cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 311*cb32e2e7SValeria Barra "Computing the mesh volume using the formula: vol = 1^T M 1\n"); 312*cb32e2e7SValeria Barra CHKERRQ(ierr); 313*cb32e2e7SValeria Barra } 314*cb32e2e7SValeria Barra 315*cb32e2e7SValeria Barra // Initialize u and v with ones 316*cb32e2e7SValeria Barra CeedVectorSetValue(uceed, 1.0); 317*cb32e2e7SValeria Barra CeedVectorSetValue(vceed, 1.0); 318*cb32e2e7SValeria Barra 319*cb32e2e7SValeria Barra // Apply the mass operator: 'u' -> 'v' 320*cb32e2e7SValeria Barra CeedOperatorApply(op_apply, uceed, vceed, CEED_REQUEST_IMMEDIATE); 321*cb32e2e7SValeria Barra CeedVectorSyncArray(vceed, CEED_MEM_HOST); 322*cb32e2e7SValeria Barra 323*cb32e2e7SValeria Barra // Gather output vector 324*cb32e2e7SValeria Barra ierr = VecRestoreArray(Vloc, &v); CHKERRQ(ierr); 325*cb32e2e7SValeria Barra ierr = VecZeroEntries(V); CHKERRQ(ierr); 326*cb32e2e7SValeria Barra ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr); 327*cb32e2e7SValeria Barra ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr); 328*cb32e2e7SValeria Barra 329*cb32e2e7SValeria Barra // Compute and print the sum of the entries of 'v' giving the mesh surface area 330*cb32e2e7SValeria Barra PetscScalar area; 331*cb32e2e7SValeria Barra ierr = VecSum(V, &area); CHKERRQ(ierr); 332*cb32e2e7SValeria Barra 333*cb32e2e7SValeria Barra // Compute the exact surface area and print the result 334*cb32e2e7SValeria Barra CeedScalar exact_surfarea = 6 * (2*l) * (2*l); 335*cb32e2e7SValeria Barra if (!test_mode) { 336*cb32e2e7SValeria Barra ierr = PetscPrintf(comm, "Exact mesh surface area : % .14g\n", 337*cb32e2e7SValeria Barra exact_surfarea); CHKERRQ(ierr); 338*cb32e2e7SValeria Barra ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area); 339*cb32e2e7SValeria Barra CHKERRQ(ierr); 340*cb32e2e7SValeria Barra ierr = PetscPrintf(comm, "Area error : % .14g\n", 341*cb32e2e7SValeria Barra fabs(area - exact_surfarea)); CHKERRQ(ierr); 342*cb32e2e7SValeria Barra } 343*cb32e2e7SValeria Barra 344*cb32e2e7SValeria Barra // PETSc cleanup 345*cb32e2e7SValeria Barra ierr = DMDestroy(&dm); CHKERRQ(ierr); 346*cb32e2e7SValeria Barra ierr = VecDestroy(&X); CHKERRQ(ierr); 347*cb32e2e7SValeria Barra ierr = VecDestroy(&Xloc); CHKERRQ(ierr); 348*cb32e2e7SValeria Barra ierr = VecDestroy(&V); CHKERRQ(ierr); 349*cb32e2e7SValeria Barra ierr = VecDestroy(&Vloc); CHKERRQ(ierr); 350*cb32e2e7SValeria Barra 351*cb32e2e7SValeria Barra // libCEED cleanup 352*cb32e2e7SValeria Barra CeedQFunctionDestroy(&qf_setupgeo); 353*cb32e2e7SValeria Barra CeedOperatorDestroy(&op_setupgeo); 354*cb32e2e7SValeria Barra CeedVectorDestroy(&xcoord); 355*cb32e2e7SValeria Barra CeedVectorDestroy(&uceed); 356*cb32e2e7SValeria Barra CeedVectorDestroy(&vceed); 357*cb32e2e7SValeria Barra CeedVectorDestroy(&qdata); 358*cb32e2e7SValeria Barra CeedBasisDestroy(&basisx); 359*cb32e2e7SValeria Barra CeedBasisDestroy(&basisu); 360*cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictu); 361*cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictx); 362*cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictxi); 363*cb32e2e7SValeria Barra CeedElemRestrictionDestroy(&Erestrictqdi); 364*cb32e2e7SValeria Barra CeedQFunctionDestroy(&qf_apply); 365*cb32e2e7SValeria Barra CeedOperatorDestroy(&op_apply); 366*cb32e2e7SValeria Barra CeedDestroy(&ceed); 367*cb32e2e7SValeria Barra return PetscFinalize(); 368*cb32e2e7SValeria Barra } 369