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