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