1 static const char help[] = "Simple libCEED test to calculate surface area using 1^T M 1"; 2 3 /* 4 This is a recreation of libCeed Example 2: https://libceed.readthedocs.io/en/latest/examples/ceed/ 5 */ 6 7 #include <petscdmceed.h> 8 #include <petscdmplexceed.h> 9 #include <petscfeceed.h> 10 #include <petscdmplex.h> 11 #include <petscds.h> 12 13 typedef struct { 14 PetscReal areaExact; 15 CeedQFunctionUser setupgeo, apply; 16 const char *setupgeofname, *applyfname; 17 } AppCtx; 18 19 typedef struct { 20 CeedQFunction qf_apply; 21 CeedOperator op_apply; 22 CeedVector qdata, uceed, vceed; 23 } CeedData; 24 25 static PetscErrorCode CeedDataDestroy(CeedData *data) 26 { 27 PetscFunctionBeginUser; 28 PetscCall(CeedVectorDestroy(&data->qdata)); 29 PetscCall(CeedVectorDestroy(&data->uceed)); 30 PetscCall(CeedVectorDestroy(&data->vceed)); 31 PetscCall(CeedQFunctionDestroy(&data->qf_apply)); 32 PetscCall(CeedOperatorDestroy(&data->op_apply)); 33 PetscFunctionReturn(0); 34 } 35 36 CEED_QFUNCTION(Mass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) 37 { 38 const CeedScalar *u = in[0], *qdata = in[1]; 39 CeedScalar *v = out[0]; 40 41 CeedPragmaSIMD 42 for (CeedInt i = 0; i < Q; ++i) 43 v[i] = qdata[i] * u[i]; 44 45 return 0; 46 } 47 48 /* 49 // Reference (parent) 2D coordinates: X \in [-1, 1]^2 50 // 51 // Global physical coordinates given by the mesh (3D): xx \in [-l, l]^3 52 // 53 // Local physical coordinates on the manifold (2D): x \in [-l, l]^2 54 // 55 // Change of coordinates matrix computed by the library: 56 // (physical 3D coords relative to reference 2D coords) 57 // dxx_j/dX_i (indicial notation) [3 * 2] 58 // 59 // Change of coordinates x (physical 2D) relative to xx (phyisical 3D): 60 // dx_i/dxx_j (indicial notation) [2 * 3] 61 // 62 // Change of coordinates x (physical 2D) relative to X (reference 2D): 63 // (by chain rule) 64 // dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j 65 // 66 // The quadrature data is stored in the array qdata. 67 // 68 // We require the determinant of the Jacobian to properly compute integrals of the form: int(u v) 69 // 70 // Qdata: w * det(dx_i/dX_j) 71 */ 72 CEED_QFUNCTION(SetupMassGeoCube)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) 73 { 74 const CeedScalar *J = in[1], *w = in[2]; 75 CeedScalar *qdata = out[0]; 76 77 CeedPragmaSIMD 78 for (CeedInt i = 0; i < Q; ++i) { 79 // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]] 80 const CeedScalar dxxdX[3][2] = {{J[i+Q*0], J[i+Q*3]}, 81 {J[i+Q*1], J[i+Q*4]}, 82 {J[i+Q*2], J[i+Q*5]}}; 83 // Modulus of dxxdX column vectors 84 const CeedScalar modg1 = PetscSqrtReal(dxxdX[0][0]*dxxdX[0][0] + dxxdX[1][0]*dxxdX[1][0] + dxxdX[2][0]*dxxdX[2][0]); 85 const CeedScalar modg2 = PetscSqrtReal(dxxdX[0][1]*dxxdX[0][1] + dxxdX[1][1]*dxxdX[1][1] + dxxdX[2][1]*dxxdX[2][1]); 86 // Use normalized column vectors of dxxdX as rows of dxdxx 87 const CeedScalar dxdxx[2][3] = {{dxxdX[0][0] / modg1, dxxdX[1][0] / modg1, dxxdX[2][0] / modg1}, 88 {dxxdX[0][1] / modg2, dxxdX[1][1] / modg2, dxxdX[2][1] / modg2}}; 89 90 CeedScalar dxdX[2][2]; 91 for (int j = 0; j < 2; ++j) 92 for (int k = 0; k < 2; ++k) { 93 dxdX[j][k] = 0; 94 for (int l = 0; l < 3; ++l) 95 dxdX[j][k] += dxdxx[j][l]*dxxdX[l][k]; 96 } 97 qdata[i+Q*0] = (dxdX[0][0]*dxdX[1][1] - dxdX[1][0]*dxdX[0][1]) * w[i]; /* det J * weight */ 98 } 99 return 0; 100 } 101 102 /* 103 // Reference (parent) 2D coordinates: X \in [-1, 1]^2 104 // 105 // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3 106 // with R radius of the sphere 107 // 108 // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3 109 // with l half edge of the cube inscribed in the sphere 110 // 111 // Change of coordinates matrix computed by the library: 112 // (physical 3D coords relative to reference 2D coords) 113 // dxx_j/dX_i (indicial notation) [3 * 2] 114 // 115 // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D): 116 // dx_i/dxx_j (indicial notation) [3 * 3] 117 // 118 // Change of coordinates x (on the 2D manifold) relative to X (reference 2D): 119 // (by chain rule) 120 // dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j [3 * 2] 121 // 122 // modJ is given by the magnitude of the cross product of the columns of dx_i/dX_j 123 // 124 // The quadrature data is stored in the array qdata. 125 // 126 // We require the determinant of the Jacobian to properly compute integrals of 127 // the form: int(u v) 128 // 129 // Qdata: modJ * w 130 */ 131 CEED_QFUNCTION(SetupMassGeoSphere)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 132 const CeedScalar *X = in[0], *J = in[1], *w = in[2]; 133 CeedScalar *qdata = out[0]; 134 135 CeedPragmaSIMD 136 for (CeedInt i = 0; i < Q; ++i) { 137 const CeedScalar xx[3][1] = {{X[i+0*Q]}, {X[i+1*Q]}, {X[i+2*Q]}}; 138 // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]] 139 const CeedScalar dxxdX[3][2] = {{J[i+Q*0], J[i+Q*3]}, 140 {J[i+Q*1], J[i+Q*4]}, 141 {J[i+Q*2], J[i+Q*5]}}; 142 // Setup 143 const CeedScalar modxxsq = xx[0][0]*xx[0][0]+xx[1][0]*xx[1][0]+xx[2][0]*xx[2][0]; 144 CeedScalar xxsq[3][3]; 145 for (int j = 0; j < 3; ++j) 146 for (int k = 0; k < 3; ++k) { 147 xxsq[j][k] = 0.; 148 for (int l = 0; l < 1; ++l) 149 xxsq[j][k] += xx[j][l]*xx[k][l] / (sqrt(modxxsq) * modxxsq); 150 } 151 152 const CeedScalar dxdxx[3][3] = {{1./sqrt(modxxsq) - xxsq[0][0], -xxsq[0][1], -xxsq[0][2]}, 153 {-xxsq[1][0], 1./sqrt(modxxsq) - xxsq[1][1], -xxsq[1][2]}, 154 {-xxsq[2][0], -xxsq[2][1], 1./sqrt(modxxsq) - xxsq[2][2]}}; 155 156 CeedScalar dxdX[3][2]; 157 for (int j = 0; j < 3; ++j) 158 for (int k = 0; k < 2; ++k) { 159 dxdX[j][k] = 0.; 160 for (int l = 0; l < 3; ++l) 161 dxdX[j][k] += dxdxx[j][l]*dxxdX[l][k]; 162 } 163 // J is given by the cross product of the columns of dxdX 164 const CeedScalar J[3][1] = {{dxdX[1][0]*dxdX[2][1] - dxdX[2][0]*dxdX[1][1]}, 165 {dxdX[2][0]*dxdX[0][1] - dxdX[0][0]*dxdX[2][1]}, 166 {dxdX[0][0]*dxdX[1][1] - dxdX[1][0]*dxdX[0][1]}}; 167 // Use the magnitude of J as our detJ (volume scaling factor) 168 const CeedScalar modJ = sqrt(J[0][0]*J[0][0]+J[1][0]*J[1][0]+J[2][0]*J[2][0]); 169 qdata[i+Q*0] = modJ * w[i]; 170 } 171 return 0; 172 } 173 174 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *ctx) 175 { 176 DMPlexShape shape = DM_SHAPE_UNKNOWN; 177 178 PetscFunctionBeginUser; 179 PetscOptionsBegin(comm, "", "libCEED Test Options", "DMPLEX"); 180 PetscOptionsEnd(); 181 PetscCall(PetscOptionsGetEnum(NULL, NULL, "-dm_plex_shape", DMPlexShapes, (PetscEnum *) &shape, NULL)); 182 ctx->setupgeo = NULL; 183 ctx->setupgeofname = NULL; 184 ctx->apply = Mass; 185 ctx->applyfname = Mass_loc; 186 ctx->areaExact = 0.0; 187 switch (shape) { 188 case DM_SHAPE_BOX_SURFACE: 189 ctx->setupgeo = SetupMassGeoCube; 190 ctx->setupgeofname = SetupMassGeoCube_loc; 191 ctx->areaExact = 6.0; 192 break; 193 case DM_SHAPE_SPHERE: 194 ctx->setupgeo = SetupMassGeoSphere; 195 ctx->setupgeofname = SetupMassGeoSphere_loc; 196 ctx->areaExact = 4.0*M_PI; 197 break; 198 default: break; 199 } 200 PetscFunctionReturn(0); 201 } 202 203 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 204 { 205 PetscFunctionBegin; 206 PetscCall(DMCreate(comm, dm)); 207 PetscCall(DMSetType(*dm, DMPLEX)); 208 PetscCall(DMSetFromOptions(*dm)); 209 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 210 #ifdef PETSC_HAVE_LIBCEED 211 { 212 Ceed ceed; 213 const char *usedresource; 214 215 PetscCall(DMGetCeed(*dm, &ceed)); 216 PetscCall(CeedGetResource(ceed, &usedresource)); 217 PetscCall(PetscPrintf(PetscObjectComm((PetscObject) *dm), "libCEED Backend: %s\n", usedresource)); 218 } 219 #endif 220 PetscFunctionReturn(0); 221 } 222 223 static PetscErrorCode SetupDiscretization(DM dm) 224 { 225 DM cdm; 226 PetscFE fe, cfe; 227 PetscInt dim, cnc; 228 PetscBool simplex; 229 230 PetscFunctionBeginUser; 231 PetscCall(DMGetDimension(dm, &dim)); 232 PetscCall(DMPlexIsSimplex(dm, &simplex)); 233 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe)); 234 PetscCall(PetscFESetName(fe, "indicator")); 235 PetscCall(DMAddField(dm, NULL, (PetscObject) fe)); 236 PetscCall(PetscFEDestroy(&fe)); 237 PetscCall(DMCreateDS(dm)); 238 PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 239 PetscCall(DMGetCoordinateDim(dm, &cnc)); 240 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, cnc, simplex, NULL, PETSC_DETERMINE, &cfe)); 241 PetscCall(DMProjectCoordinates(dm, cfe)); 242 PetscCall(PetscFEDestroy(&cfe)); 243 PetscCall(DMGetCoordinateDM(dm, &cdm)); 244 PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL)); 245 PetscFunctionReturn(0); 246 } 247 248 static PetscErrorCode LibCeedSetupByDegree(DM dm, AppCtx *ctx, CeedData *data) 249 { 250 PetscDS ds; 251 PetscFE fe, cfe; 252 Ceed ceed; 253 CeedElemRestriction Erestrictx, Erestrictu, Erestrictq; 254 CeedQFunction qf_setupgeo; 255 CeedOperator op_setupgeo; 256 CeedVector xcoord; 257 CeedBasis basisu, basisx; 258 CeedInt Nqdata = 1; 259 CeedInt nqpts, nqptsx; 260 DM cdm; 261 Vec coords; 262 const PetscScalar *coordArray; 263 PetscInt dim, cdim, cStart, cEnd, Ncell; 264 265 PetscFunctionBeginUser; 266 PetscCall(DMGetCeed(dm, &ceed)); 267 PetscCall(DMGetDimension(dm, &dim)); 268 PetscCall(DMGetCoordinateDim(dm, &cdim)); 269 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 270 Ncell = cEnd - cStart; 271 // CEED bases 272 PetscCall(DMGetDS(dm, &ds)); 273 PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *) &fe)); 274 PetscCall(PetscFEGetCeedBasis(fe, &basisu)); 275 PetscCall(DMGetCoordinateDM(dm, &cdm)); 276 PetscCall(DMGetDS(cdm, &ds)); 277 PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *) &cfe)); 278 PetscCall(PetscFEGetCeedBasis(cfe, &basisx)); 279 280 PetscCall(DMPlexGetCeedRestriction(cdm, NULL, 0, 0, 0, &Erestrictx)); 281 PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &Erestrictu)); 282 PetscCall(CeedBasisGetNumQuadraturePoints(basisu, &nqpts)); 283 PetscCall(CeedBasisGetNumQuadraturePoints(basisx, &nqptsx)); 284 PetscCheck(nqptsx == nqpts,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of qpoints for u %" PetscInt_FMT " != %" PetscInt_FMT " Number of qpoints for x", nqpts, nqptsx); 285 PetscCall(CeedElemRestrictionCreateStrided(ceed, Ncell, nqpts, Nqdata, Nqdata*Ncell*nqpts, CEED_STRIDES_BACKEND, &Erestrictq)); 286 287 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 288 PetscCall(VecGetArrayRead(coords, &coordArray)); 289 PetscCall(CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL)); 290 PetscCall(CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *) coordArray)); 291 PetscCall(VecRestoreArrayRead(coords, &coordArray)); 292 293 // Create the vectors that will be needed in setup and apply 294 PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->uceed, NULL)); 295 PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->vceed, NULL)); 296 PetscCall(CeedElemRestrictionCreateVector(Erestrictq, &data->qdata, NULL)); 297 298 // Create the Q-function that builds the operator (i.e. computes its quadrature data) and set its context data 299 PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->setupgeo, ctx->setupgeofname, &qf_setupgeo)); 300 PetscCall(CeedQFunctionAddInput(qf_setupgeo, "x", cdim, CEED_EVAL_INTERP)); 301 PetscCall(CeedQFunctionAddInput(qf_setupgeo, "dx", cdim*dim, CEED_EVAL_GRAD)); 302 PetscCall(CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT)); 303 PetscCall(CeedQFunctionAddOutput(qf_setupgeo, "qdata", Nqdata, CEED_EVAL_NONE)); 304 305 // Set up the mass operator 306 PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->apply, ctx->applyfname, &data->qf_apply)); 307 PetscCall(CeedQFunctionAddInput(data->qf_apply, "u", 1, CEED_EVAL_INTERP)); 308 PetscCall(CeedQFunctionAddInput(data->qf_apply, "qdata", Nqdata, CEED_EVAL_NONE)); 309 PetscCall(CeedQFunctionAddOutput(data->qf_apply, "v", 1, CEED_EVAL_INTERP)); 310 311 // Create the operator that builds the quadrature data for the operator 312 PetscCall(CeedOperatorCreate(ceed, qf_setupgeo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setupgeo)); 313 PetscCall(CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx, CEED_VECTOR_ACTIVE)); 314 PetscCall(CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx, CEED_VECTOR_ACTIVE)); 315 PetscCall(CeedOperatorSetField(op_setupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx, CEED_VECTOR_NONE)); 316 PetscCall(CeedOperatorSetField(op_setupgeo, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 317 318 // Create the mass operator 319 PetscCall(CeedOperatorCreate(ceed, data->qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &data->op_apply)); 320 PetscCall(CeedOperatorSetField(data->op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE)); 321 PetscCall(CeedOperatorSetField(data->op_apply, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, data->qdata)); 322 PetscCall(CeedOperatorSetField(data->op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE)); 323 324 // Setup qdata 325 PetscCall(CeedOperatorApply(op_setupgeo, xcoord, data->qdata, CEED_REQUEST_IMMEDIATE)); 326 327 PetscCall(CeedElemRestrictionDestroy(&Erestrictq)); 328 PetscCall(CeedQFunctionDestroy(&qf_setupgeo)); 329 PetscCall(CeedOperatorDestroy(&op_setupgeo)); 330 PetscCall(CeedVectorDestroy(&xcoord)); 331 PetscFunctionReturn(0); 332 } 333 334 int main(int argc, char **argv) 335 { 336 MPI_Comm comm; 337 DM dm; 338 AppCtx ctx; 339 Vec U, Uloc, V, Vloc; 340 PetscScalar *v; 341 PetscScalar area; 342 CeedData ceeddata; 343 344 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 345 comm = PETSC_COMM_WORLD; 346 PetscCall(ProcessOptions(comm, &ctx)); 347 PetscCall(CreateMesh(comm, &ctx, &dm)); 348 PetscCall(SetupDiscretization(dm)); 349 350 PetscCall(LibCeedSetupByDegree(dm, &ctx, &ceeddata)); 351 352 PetscCall(DMCreateGlobalVector(dm, &U)); 353 PetscCall(DMCreateLocalVector(dm, &Uloc)); 354 PetscCall(VecDuplicate(U, &V)); 355 PetscCall(VecDuplicate(Uloc, &Vloc)); 356 357 /**/ 358 PetscCall(VecSet(Uloc, 1.)); 359 PetscCall(VecZeroEntries(V)); 360 PetscCall(VecZeroEntries(Vloc)); 361 PetscCall(VecGetArray(Vloc, &v)); 362 PetscCall(CeedVectorSetArray(ceeddata.vceed, CEED_MEM_HOST, CEED_USE_POINTER, v)); 363 PetscCall(CeedVectorSetValue(ceeddata.uceed, 1.0)); 364 PetscCall(CeedOperatorApply(ceeddata.op_apply, ceeddata.uceed, ceeddata.vceed, CEED_REQUEST_IMMEDIATE)); 365 PetscCall(CeedVectorTakeArray(ceeddata.vceed, CEED_MEM_HOST, NULL)); 366 PetscCall(VecRestoreArray(Vloc, &v)); 367 PetscCall(DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V)); 368 PetscCall(DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V)); 369 370 PetscCall(VecSum(V, &area)); 371 if (ctx.areaExact > 0.) { 372 PetscReal error = PetscAbsReal(area - ctx.areaExact); 373 PetscReal tol = PETSC_SMALL; 374 375 PetscCall(PetscPrintf(comm, "Exact mesh surface area : % .*f\n", PetscAbsReal(ctx.areaExact - round(ctx.areaExact)) > 1E-15 ? 14 : 1, (double) ctx.areaExact)); 376 PetscCall(PetscPrintf(comm, "Computed mesh surface area : % .*f\n", PetscAbsScalar(area - round(area)) > 1E-15 ? 14 : 1, (double)PetscRealPart(area))); 377 if (error > tol) { 378 PetscCall(PetscPrintf(comm, "Area error : % .14g\n", (double) error)); 379 } else { 380 PetscCall(PetscPrintf(comm, "Area verifies!\n")); 381 } 382 } 383 384 PetscCall(CeedDataDestroy(&ceeddata)); 385 PetscCall(VecDestroy(&U)); 386 PetscCall(VecDestroy(&Uloc)); 387 PetscCall(VecDestroy(&V)); 388 PetscCall(VecDestroy(&Vloc)); 389 PetscCall(DMDestroy(&dm)); 390 return PetscFinalize(); 391 } 392 393 /*TEST 394 395 build: 396 requires: libceed 397 398 testset: 399 args: -dm_plex_simplex 0 -petscspace_degree 3 -dm_view -dm_petscds_view \ 400 -petscfe_default_quadrature_order 4 -coord_dm_default_quadrature_order 4 401 402 test: 403 suffix: cube_3 404 args: -dm_plex_shape box_surface -dm_refine 2 405 406 test: 407 suffix: cube_3_p4 408 nsize: 4 409 args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape box_surface -dm_refine 1 410 411 test: 412 suffix: sphere_3 413 args: -dm_plex_shape sphere -dm_refine 3 414 415 test: 416 suffix: sphere_3_p4 417 nsize: 4 418 args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape sphere -dm_refine 2 419 420 TEST*/ 421