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, CeedTransposeMode lmode, CeedInt P, 64 CeedInt ncomp, CeedElemRestriction *Erestrict, 65 DM dm) { 66 PetscInt ierr; 67 PetscInt c, cStart, cEnd, nelem, nnodes, *erestrict, eoffset; 68 PetscSection section; 69 Vec Uloc; 70 71 PetscFunctionBegin; 72 73 // Get Nelem 74 ierr = DMGetSection(dm, §ion); CHKERRQ(ierr); 75 ierr = DMPlexGetHeightStratum(dm, 0, &cStart,& cEnd); CHKERRQ(ierr); 76 nelem = cEnd - cStart; 77 78 // Get indices 79 ierr = PetscMalloc1(nelem*P*P, &erestrict); CHKERRQ(ierr); 80 for (c=cStart, eoffset = 0; c<cEnd; c++) { 81 PetscInt numindices, *indices, i; 82 ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices, 83 &indices, NULL); CHKERRQ(ierr); 84 for (i=0; i<numindices; i+=ncomp) { 85 for (PetscInt j=0; j<ncomp; j++) { 86 if (indices[i+j] != indices[i] + (PetscInt)(copysign(j, indices[i]))) 87 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 88 "Cell %D closure indices not interlaced", c); 89 } 90 // NO BC on closed surfaces 91 PetscInt loc = indices[i]; 92 erestrict[eoffset++] = loc/ncomp; 93 } 94 ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices, 95 &indices, NULL); CHKERRQ(ierr); 96 } 97 98 // Setup CEED restriction 99 ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 100 ierr = VecGetLocalSize(Uloc, &nnodes); CHKERRQ(ierr); 101 102 ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 103 CeedElemRestrictionCreate(ceed, lmode, nelem, P*P, nnodes/ncomp, ncomp, 104 CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, 105 Erestrict); 106 ierr = PetscFree(erestrict); CHKERRQ(ierr); 107 108 PetscFunctionReturn(0); 109 } 110 111 // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c 112 static PetscErrorCode ProjectToUnitSphere(DM dm) { 113 Vec coordinates; 114 PetscScalar *coords; 115 PetscInt Nv, v, dim, d; 116 PetscErrorCode ierr; 117 118 PetscFunctionBeginUser; 119 ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr); 120 ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr); 121 ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr); 122 Nv /= dim; 123 ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr); 124 for (v = 0; v < Nv; ++v) { 125 PetscReal r = 0.0; 126 127 for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d])); 128 r = PetscSqrtReal(r); 129 for (d = 0; d < dim; ++d) coords[v*dim+d] /= r; 130 } 131 ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr); 132 PetscFunctionReturn(0); 133 } 134 135 int main(int argc, char **argv) { 136 PetscInt ierr; 137 MPI_Comm comm; 138 char filename[PETSC_MAX_PATH_LEN], 139 ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 140 PetscInt lsize, gsize, xlsize, 141 qextra = 1, // default number of extra quadrature points 142 ncompx = 3, // number of components of 3D physical coordinates 143 ncompu = 1, // dimension of field to which apply mass operator 144 topodim = 2, // topological dimension of manifold 145 degree = 3; // default degree for finite element bases 146 PetscBool read_mesh = PETSC_FALSE, 147 test_mode = PETSC_FALSE; 148 PetscSpace sp; 149 PetscFE fe; 150 Vec X, Xloc, V, Vloc; 151 DM dm, dmcoord; 152 Ceed ceed; 153 CeedInt P, Q; 154 CeedOperator op_setupgeo, op_apply; 155 CeedQFunction qf_setupgeo, qf_apply; 156 CeedBasis basisx, basisu; 157 CeedTransposeMode lmodeceed = CEED_NOTRANSPOSE, lmodepetsc = CEED_TRANSPOSE; 158 CeedElemRestriction Erestrictx, Erestrictu, Erestrictxi, 159 Erestrictqdi; 160 PetscFunctionList geomfactorlist = NULL; 161 char problemtype[PETSC_MAX_PATH_LEN] = "sphere"; 162 163 ierr = PetscInitialize(&argc, &argv, NULL, help); 164 if (ierr) return ierr; 165 comm = PETSC_COMM_WORLD; 166 167 // Set up problem type command line option 168 ierr = PetscFunctionListAdd(&geomfactorlist, "cube", &SetupMassGeoCube); 169 CHKERRQ(ierr); 170 ierr = PetscFunctionListAdd(&geomfactorlist, "sphere", &SetupMassGeoSphere); 171 CHKERRQ(ierr); 172 173 // Read command line options 174 ierr = PetscOptionsBegin(comm, NULL, "CEED surface area problem with PETSc", 175 NULL); 176 CHKERRQ(ierr); 177 ierr = PetscOptionsFList("-problem", "Problem to solve", NULL, geomfactorlist, 178 problemtype, problemtype, sizeof problemtype, NULL); 179 CHKERRQ(ierr); 180 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 181 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 182 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 183 NULL, ceedresource, ceedresource, 184 sizeof(ceedresource), NULL); CHKERRQ(ierr); 185 ierr = PetscOptionsBool("-test", 186 "Testing mode (do not print unless error is large)", 187 NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 188 ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL, 189 filename, filename, sizeof(filename), &read_mesh); 190 CHKERRQ(ierr); 191 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 192 193 // Setup function pointer for geometric factors 194 int (*geomfp)(void *ctx, const CeedInt Q, const CeedScalar *const *in, 195 CeedScalar *const *out); 196 ierr = PetscFunctionListFind(geomfactorlist, problemtype, 197 (void(* *)(void))&geomfp); CHKERRQ(ierr); 198 const char *str; 199 if (geomfp == SetupMassGeoCube) 200 str = SetupMassGeoCube_loc; 201 else if (geomfp == SetupMassGeoSphere) 202 str = SetupMassGeoSphere_loc; 203 else 204 return CeedError(ceed, 1, "Function not found in the list"); 205 206 // Setup DM 207 if (read_mesh) { 208 ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm); 209 CHKERRQ(ierr); 210 } else { 211 // Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box. 212 PetscBool simplex = PETSC_FALSE; 213 ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, &dm); 214 CHKERRQ(ierr); 215 // Set the object name 216 ierr = PetscObjectSetName((PetscObject)dm, problemtype); CHKERRQ(ierr); 217 // Distribute mesh over processes 218 { 219 DM dmDist = NULL; 220 PetscPartitioner part; 221 222 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 223 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 224 ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 225 if (dmDist) { 226 ierr = DMDestroy(&dm); CHKERRQ(ierr); 227 dm = dmDist; 228 } 229 } 230 // Refine DMPlex with uniform refinement using runtime option -dm_refine 231 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 232 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 233 if (!strcmp(problemtype, "sphere")) 234 ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr); 235 // View DMPlex via runtime option 236 ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 237 } 238 239 // Create FE 240 ierr = PetscFECreateDefault(PETSC_COMM_SELF, topodim, ncompu, PETSC_FALSE, NULL, 241 PETSC_DETERMINE, &fe); 242 CHKERRQ(ierr); 243 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 244 ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr); 245 ierr = DMCreateDS(dm); CHKERRQ(ierr); 246 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 247 CHKERRQ(ierr); 248 249 // Get basis space degree 250 ierr = PetscFEGetBasisSpace(fe, &sp); CHKERRQ(ierr); 251 ierr = PetscSpaceGetDegree(sp, °ree, NULL); CHKERRQ(ierr); 252 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 253 if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 254 "-petscspace_degree %D must be at least 1", degree); 255 256 // Create vectors 257 ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr); 258 ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr); 259 ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 260 ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr); 261 ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr); 262 ierr = VecDuplicate(X, &V); CHKERRQ(ierr); 263 ierr = VecDuplicate(Xloc, &Vloc); CHKERRQ(ierr); 264 265 // Set up libCEED 266 CeedInit(ceedresource, &ceed); 267 268 // Print summary 269 P = degree + 1; 270 Q = P + qextra; 271 const char *usedresource; 272 CeedGetResource(ceed, &usedresource); 273 if (!test_mode) { 274 ierr = PetscPrintf(comm, 275 "\n-- libCEED + PETSc Surface Area of a Manifold --\n" 276 " libCEED:\n" 277 " libCEED Backend : %s\n" 278 " Mesh:\n" 279 " Number of 1D Basis Nodes (p) : %d\n" 280 " Number of 1D Quadrature Points (q) : %d\n" 281 " Global nodes : %D\n" 282 " DoF per node : %D\n", 283 usedresource, P, Q, gsize/ncompu, ncompu); 284 CHKERRQ(ierr); 285 } 286 287 // Setup libCEED's objects: 288 // Create bases 289 CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompu, P, Q, 290 CEED_GAUSS, &basisu); 291 CeedBasisCreateTensorH1Lagrange(ceed, topodim, ncompx, 2, Q, 292 CEED_GAUSS, &basisx); 293 294 // CEED restrictions 295 ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 296 ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 297 CHKERRQ(ierr); 298 299 CreateRestrictionPlex(ceed, lmodepetsc, 2, ncompx, &Erestrictx, dmcoord); 300 CHKERRQ(ierr); 301 CreateRestrictionPlex(ceed, lmodepetsc, P, ncompu, &Erestrictu, dm); 302 CHKERRQ(ierr); 303 304 CeedInt cStart, cEnd; 305 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 306 const CeedInt nelem = cEnd - cStart; 307 308 // CEED identity restrictions 309 const CeedInt qdatasize = 1; 310 CeedElemRestrictionCreateIdentity(ceed, lmodeceed, nelem, Q*Q, nelem*Q*Q, 311 qdatasize, &Erestrictqdi); 312 CeedElemRestrictionCreateIdentity(ceed, lmodeceed, nelem, Q*Q, nelem*Q*Q, 1, 313 &Erestrictxi); 314 315 // Element coordinates 316 Vec coords; 317 const PetscScalar *coordArray; 318 PetscSection section; 319 ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr); 320 ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr); 321 ierr = DMGetSection(dmcoord, §ion); CHKERRQ(ierr); 322 323 CeedVector xcoord; 324 CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL); 325 CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, 326 (PetscScalar *)coordArray); 327 ierr = VecRestoreArrayRead(coords, &coordArray); 328 329 // Create the vectors that will be needed in setup and apply 330 CeedVector uceed, vceed, qdata; 331 CeedInt nqpts; 332 CeedBasisGetNumQuadraturePoints(basisu, &nqpts); 333 CeedVectorCreate(ceed, qdatasize*nelem*nqpts, &qdata); 334 CeedVectorCreate(ceed, xlsize, &uceed); 335 CeedVectorCreate(ceed, xlsize, &vceed); 336 337 // Create the Q-function that builds the operator for the geomteric factors 338 // (i.e., the quadrature data) 339 CeedQFunctionCreateInterior(ceed, 1, geomfp, str, &qf_setupgeo); 340 CeedQFunctionAddInput(qf_setupgeo, "x", ncompx, CEED_EVAL_INTERP); 341 CeedQFunctionAddInput(qf_setupgeo, "dx", ncompx*topodim, CEED_EVAL_GRAD); 342 CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT); 343 CeedQFunctionAddOutput(qf_setupgeo, "qdata", qdatasize, CEED_EVAL_NONE); 344 345 // Set up the mass operator 346 CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_apply); 347 CeedQFunctionAddInput(qf_apply, "u", ncompu, CEED_EVAL_INTERP); 348 CeedQFunctionAddInput(qf_apply, "qdata", qdatasize, CEED_EVAL_NONE); 349 CeedQFunctionAddOutput(qf_apply, "v", ncompu, CEED_EVAL_INTERP); 350 351 // Create the operator that builds the quadrature data for the operator 352 CeedOperatorCreate(ceed, qf_setupgeo, NULL, NULL, &op_setupgeo); 353 CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx, 354 CEED_VECTOR_ACTIVE); 355 CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx, 356 CEED_VECTOR_ACTIVE); 357 CeedOperatorSetField(op_setupgeo, "weight", Erestrictxi, basisx, 358 CEED_VECTOR_NONE); 359 CeedOperatorSetField(op_setupgeo, "qdata", Erestrictqdi, 360 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 361 362 // Create the mass operator 363 CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply); 364 CeedOperatorSetField(op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE); 365 CeedOperatorSetField(op_apply, "qdata", Erestrictqdi, CEED_BASIS_COLLOCATED, 366 qdata); 367 CeedOperatorSetField(op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE); 368 369 // Compute the quadrature data for the mass operator 370 CeedOperatorApply(op_setupgeo, xcoord, qdata, CEED_REQUEST_IMMEDIATE); 371 372 PetscScalar *v; 373 ierr = VecZeroEntries(Vloc); CHKERRQ(ierr); 374 ierr = VecGetArray(Vloc, &v); 375 CeedVectorSetArray(vceed, CEED_MEM_HOST, CEED_USE_POINTER, v); 376 377 // Compute the mesh volume using the mass operator: vol = 1^T \cdot M \cdot 1 378 if (!test_mode) { 379 ierr = PetscPrintf(comm, 380 "Computing the mesh area using the formula: area = 1^T M 1\n"); 381 CHKERRQ(ierr); 382 } 383 384 // Initialize u and v with ones 385 CeedVectorSetValue(uceed, 1.0); 386 387 // Apply the mass operator: 'u' -> 'v' 388 CeedOperatorApply(op_apply, uceed, vceed, CEED_REQUEST_IMMEDIATE); 389 CeedVectorSyncArray(vceed, CEED_MEM_HOST); 390 391 // Gather output vector 392 ierr = VecRestoreArray(Vloc, &v); CHKERRQ(ierr); 393 ierr = VecZeroEntries(V); CHKERRQ(ierr); 394 ierr = DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr); 395 ierr = DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V); CHKERRQ(ierr); 396 397 // Compute and print the sum of the entries of 'v' giving the mesh surface area 398 PetscScalar area; 399 ierr = VecSum(V, &area); CHKERRQ(ierr); 400 401 // Compute the exact surface area and print the result 402 CeedScalar exact_surfarea = 4 * M_PI; 403 if (!strcmp(problemtype, "cube")) { 404 PetscScalar l = 1.0/PetscSqrtReal(3.0); // half edge of the cube 405 exact_surfarea = 6 * (2*l) * (2*l); 406 } 407 408 if (!test_mode) { 409 ierr = PetscPrintf(comm, "Exact mesh surface area : % .14g\n", 410 exact_surfarea); CHKERRQ(ierr); 411 ierr = PetscPrintf(comm, "Computed mesh surface area : % .14g\n", area); 412 CHKERRQ(ierr); 413 ierr = PetscPrintf(comm, "Area error : % .14g\n", 414 fabs(area - exact_surfarea)); CHKERRQ(ierr); 415 } 416 417 // PETSc cleanup 418 ierr = DMDestroy(&dm); CHKERRQ(ierr); 419 ierr = VecDestroy(&X); CHKERRQ(ierr); 420 ierr = VecDestroy(&Xloc); CHKERRQ(ierr); 421 ierr = VecDestroy(&V); CHKERRQ(ierr); 422 ierr = VecDestroy(&Vloc); CHKERRQ(ierr); 423 424 // libCEED cleanup 425 CeedQFunctionDestroy(&qf_setupgeo); 426 CeedOperatorDestroy(&op_setupgeo); 427 CeedVectorDestroy(&xcoord); 428 CeedVectorDestroy(&uceed); 429 CeedVectorDestroy(&vceed); 430 CeedVectorDestroy(&qdata); 431 CeedBasisDestroy(&basisx); 432 CeedBasisDestroy(&basisu); 433 CeedElemRestrictionDestroy(&Erestrictu); 434 CeedElemRestrictionDestroy(&Erestrictx); 435 CeedElemRestrictionDestroy(&Erestrictxi); 436 CeedElemRestrictionDestroy(&Erestrictqdi); 437 CeedQFunctionDestroy(&qf_apply); 438 CeedOperatorDestroy(&op_apply); 439 CeedDestroy(&ceed); 440 return PetscFinalize(); 441 } 442