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 PetscFunctionBeginUser; 289566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&data->qdata)); 299566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&data->uceed)); 309566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&data->vceed)); 319566063dSJacob Faibussowitsch PetscCall(CeedQFunctionDestroy(&data->qf_apply)); 329566063dSJacob Faibussowitsch PetscCall(CeedOperatorDestroy(&data->op_apply)); 33f918ec44SMatthew G. Knepley PetscFunctionReturn(0); 34f918ec44SMatthew G. Knepley } 35f918ec44SMatthew G. Knepley 36f918ec44SMatthew G. Knepley CEED_QFUNCTION(Mass)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) 37f918ec44SMatthew G. Knepley { 38f918ec44SMatthew G. Knepley const CeedScalar *u = in[0], *qdata = in[1]; 39f918ec44SMatthew G. Knepley CeedScalar *v = out[0]; 40f918ec44SMatthew G. Knepley 41f918ec44SMatthew G. Knepley CeedPragmaSIMD 42f918ec44SMatthew G. Knepley for (CeedInt i = 0; i < Q; ++i) 43f918ec44SMatthew G. Knepley v[i] = qdata[i] * u[i]; 44f918ec44SMatthew G. Knepley 45f918ec44SMatthew G. Knepley return 0; 46f918ec44SMatthew G. Knepley } 47f918ec44SMatthew G. Knepley 48f918ec44SMatthew G. Knepley /* 49f918ec44SMatthew G. Knepley // Reference (parent) 2D coordinates: X \in [-1, 1]^2 50f918ec44SMatthew G. Knepley // 51f918ec44SMatthew G. Knepley // Global physical coordinates given by the mesh (3D): xx \in [-l, l]^3 52f918ec44SMatthew G. Knepley // 53f918ec44SMatthew G. Knepley // Local physical coordinates on the manifold (2D): x \in [-l, l]^2 54f918ec44SMatthew G. Knepley // 55f918ec44SMatthew G. Knepley // Change of coordinates matrix computed by the library: 56f918ec44SMatthew G. Knepley // (physical 3D coords relative to reference 2D coords) 57f918ec44SMatthew G. Knepley // dxx_j/dX_i (indicial notation) [3 * 2] 58f918ec44SMatthew G. Knepley // 59f918ec44SMatthew G. Knepley // Change of coordinates x (physical 2D) relative to xx (phyisical 3D): 60f918ec44SMatthew G. Knepley // dx_i/dxx_j (indicial notation) [2 * 3] 61f918ec44SMatthew G. Knepley // 62f918ec44SMatthew G. Knepley // Change of coordinates x (physical 2D) relative to X (reference 2D): 63f918ec44SMatthew G. Knepley // (by chain rule) 64f918ec44SMatthew G. Knepley // dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j 65f918ec44SMatthew G. Knepley // 66f918ec44SMatthew G. Knepley // The quadrature data is stored in the array qdata. 67f918ec44SMatthew G. Knepley // 68f918ec44SMatthew G. Knepley // We require the determinant of the Jacobian to properly compute integrals of the form: int(u v) 69f918ec44SMatthew G. Knepley // 70f918ec44SMatthew G. Knepley // Qdata: w * det(dx_i/dX_j) 71f918ec44SMatthew G. Knepley */ 72f918ec44SMatthew G. Knepley CEED_QFUNCTION(SetupMassGeoCube)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) 73f918ec44SMatthew G. Knepley { 74f918ec44SMatthew G. Knepley const CeedScalar *J = in[1], *w = in[2]; 75f918ec44SMatthew G. Knepley CeedScalar *qdata = out[0]; 76f918ec44SMatthew G. Knepley 77f918ec44SMatthew G. Knepley CeedPragmaSIMD 78f918ec44SMatthew G. Knepley for (CeedInt i = 0; i < Q; ++i) { 79f918ec44SMatthew G. Knepley // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]] 80f918ec44SMatthew G. Knepley const CeedScalar dxxdX[3][2] = {{J[i+Q*0], J[i+Q*3]}, 81f918ec44SMatthew G. Knepley {J[i+Q*1], J[i+Q*4]}, 82f918ec44SMatthew G. Knepley {J[i+Q*2], J[i+Q*5]}}; 83f918ec44SMatthew G. Knepley // Modulus of dxxdX column vectors 84f918ec44SMatthew G. Knepley const CeedScalar modg1 = PetscSqrtReal(dxxdX[0][0]*dxxdX[0][0] + dxxdX[1][0]*dxxdX[1][0] + dxxdX[2][0]*dxxdX[2][0]); 85f918ec44SMatthew G. Knepley const CeedScalar modg2 = PetscSqrtReal(dxxdX[0][1]*dxxdX[0][1] + dxxdX[1][1]*dxxdX[1][1] + dxxdX[2][1]*dxxdX[2][1]); 86f918ec44SMatthew G. Knepley // Use normalized column vectors of dxxdX as rows of dxdxx 87f918ec44SMatthew G. Knepley const CeedScalar dxdxx[2][3] = {{dxxdX[0][0] / modg1, dxxdX[1][0] / modg1, dxxdX[2][0] / modg1}, 88f918ec44SMatthew G. Knepley {dxxdX[0][1] / modg2, dxxdX[1][1] / modg2, dxxdX[2][1] / modg2}}; 89f918ec44SMatthew G. Knepley 90f918ec44SMatthew G. Knepley CeedScalar dxdX[2][2]; 91f918ec44SMatthew G. Knepley for (int j = 0; j < 2; ++j) 92f918ec44SMatthew G. Knepley for (int k = 0; k < 2; ++k) { 93f918ec44SMatthew G. Knepley dxdX[j][k] = 0; 94f918ec44SMatthew G. Knepley for (int l = 0; l < 3; ++l) 95f918ec44SMatthew G. Knepley dxdX[j][k] += dxdxx[j][l]*dxxdX[l][k]; 96f918ec44SMatthew G. Knepley } 97f918ec44SMatthew G. Knepley qdata[i+Q*0] = (dxdX[0][0]*dxdX[1][1] - dxdX[1][0]*dxdX[0][1]) * w[i]; /* det J * weight */ 98f918ec44SMatthew G. Knepley } 99f918ec44SMatthew G. Knepley return 0; 100f918ec44SMatthew G. Knepley } 101f918ec44SMatthew G. Knepley 102f918ec44SMatthew G. Knepley /* 103f918ec44SMatthew G. Knepley // Reference (parent) 2D coordinates: X \in [-1, 1]^2 104f918ec44SMatthew G. Knepley // 105f918ec44SMatthew G. Knepley // Global 3D physical coordinates given by the mesh: xx \in [-R, R]^3 106f918ec44SMatthew G. Knepley // with R radius of the sphere 107f918ec44SMatthew G. Knepley // 108f918ec44SMatthew G. Knepley // Local 3D physical coordinates on the 2D manifold: x \in [-l, l]^3 109f918ec44SMatthew G. Knepley // with l half edge of the cube inscribed in the sphere 110f918ec44SMatthew G. Knepley // 111f918ec44SMatthew G. Knepley // Change of coordinates matrix computed by the library: 112f918ec44SMatthew G. Knepley // (physical 3D coords relative to reference 2D coords) 113f918ec44SMatthew G. Knepley // dxx_j/dX_i (indicial notation) [3 * 2] 114f918ec44SMatthew G. Knepley // 115f918ec44SMatthew G. Knepley // Change of coordinates x (on the 2D manifold) relative to xx (phyisical 3D): 116f918ec44SMatthew G. Knepley // dx_i/dxx_j (indicial notation) [3 * 3] 117f918ec44SMatthew G. Knepley // 118f918ec44SMatthew G. Knepley // Change of coordinates x (on the 2D manifold) relative to X (reference 2D): 119f918ec44SMatthew G. Knepley // (by chain rule) 120f918ec44SMatthew G. Knepley // dx_i/dX_j = dx_i/dxx_k * dxx_k/dX_j [3 * 2] 121f918ec44SMatthew G. Knepley // 122f918ec44SMatthew G. Knepley // modJ is given by the magnitude of the cross product of the columns of dx_i/dX_j 123f918ec44SMatthew G. Knepley // 124f918ec44SMatthew G. Knepley // The quadrature data is stored in the array qdata. 125f918ec44SMatthew G. Knepley // 126f918ec44SMatthew G. Knepley // We require the determinant of the Jacobian to properly compute integrals of 127f918ec44SMatthew G. Knepley // the form: int(u v) 128f918ec44SMatthew G. Knepley // 129f918ec44SMatthew G. Knepley // Qdata: modJ * w 130f918ec44SMatthew G. Knepley */ 131f918ec44SMatthew G. Knepley CEED_QFUNCTION(SetupMassGeoSphere)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 132f918ec44SMatthew G. Knepley const CeedScalar *X = in[0], *J = in[1], *w = in[2]; 133f918ec44SMatthew G. Knepley CeedScalar *qdata = out[0]; 134f918ec44SMatthew G. Knepley 135f918ec44SMatthew G. Knepley CeedPragmaSIMD 136f918ec44SMatthew G. Knepley for (CeedInt i = 0; i < Q; ++i) { 137f918ec44SMatthew G. Knepley const CeedScalar xx[3][1] = {{X[i+0*Q]}, {X[i+1*Q]}, {X[i+2*Q]}}; 138f918ec44SMatthew G. Knepley // Read dxxdX Jacobian entries, stored as [[0 3], [1 4], [2 5]] 139f918ec44SMatthew G. Knepley const CeedScalar dxxdX[3][2] = {{J[i+Q*0], J[i+Q*3]}, 140f918ec44SMatthew G. Knepley {J[i+Q*1], J[i+Q*4]}, 141f918ec44SMatthew G. Knepley {J[i+Q*2], J[i+Q*5]}}; 142f918ec44SMatthew G. Knepley // Setup 143f918ec44SMatthew G. Knepley const CeedScalar modxxsq = xx[0][0]*xx[0][0]+xx[1][0]*xx[1][0]+xx[2][0]*xx[2][0]; 144f918ec44SMatthew G. Knepley CeedScalar xxsq[3][3]; 145f918ec44SMatthew G. Knepley for (int j = 0; j < 3; ++j) 146f918ec44SMatthew G. Knepley for (int k = 0; k < 3; ++k) { 147f918ec44SMatthew G. Knepley xxsq[j][k] = 0.; 148f918ec44SMatthew G. Knepley for (int l = 0; l < 1; ++l) 149f918ec44SMatthew G. Knepley xxsq[j][k] += xx[j][l]*xx[k][l] / (sqrt(modxxsq) * modxxsq); 150f918ec44SMatthew G. Knepley } 151f918ec44SMatthew G. Knepley 152f918ec44SMatthew G. Knepley const CeedScalar dxdxx[3][3] = {{1./sqrt(modxxsq) - xxsq[0][0], -xxsq[0][1], -xxsq[0][2]}, 153f918ec44SMatthew G. Knepley {-xxsq[1][0], 1./sqrt(modxxsq) - xxsq[1][1], -xxsq[1][2]}, 154f918ec44SMatthew G. Knepley {-xxsq[2][0], -xxsq[2][1], 1./sqrt(modxxsq) - xxsq[2][2]}}; 155f918ec44SMatthew G. Knepley 156f918ec44SMatthew G. Knepley CeedScalar dxdX[3][2]; 157f918ec44SMatthew G. Knepley for (int j = 0; j < 3; ++j) 158f918ec44SMatthew G. Knepley for (int k = 0; k < 2; ++k) { 159f918ec44SMatthew G. Knepley dxdX[j][k] = 0.; 160f918ec44SMatthew G. Knepley for (int l = 0; l < 3; ++l) 161f918ec44SMatthew G. Knepley dxdX[j][k] += dxdxx[j][l]*dxxdX[l][k]; 162f918ec44SMatthew G. Knepley } 163f918ec44SMatthew G. Knepley // J is given by the cross product of the columns of dxdX 164f918ec44SMatthew G. Knepley const CeedScalar J[3][1] = {{dxdX[1][0]*dxdX[2][1] - dxdX[2][0]*dxdX[1][1]}, 165f918ec44SMatthew G. Knepley {dxdX[2][0]*dxdX[0][1] - dxdX[0][0]*dxdX[2][1]}, 166f918ec44SMatthew G. Knepley {dxdX[0][0]*dxdX[1][1] - dxdX[1][0]*dxdX[0][1]}}; 167f918ec44SMatthew G. Knepley // Use the magnitude of J as our detJ (volume scaling factor) 168f918ec44SMatthew G. Knepley const CeedScalar modJ = sqrt(J[0][0]*J[0][0]+J[1][0]*J[1][0]+J[2][0]*J[2][0]); 169f918ec44SMatthew G. Knepley qdata[i+Q*0] = modJ * w[i]; 170f918ec44SMatthew G. Knepley } 171f918ec44SMatthew G. Knepley return 0; 172f918ec44SMatthew G. Knepley } 173f918ec44SMatthew G. Knepley 174f918ec44SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *ctx) 175f918ec44SMatthew G. Knepley { 176f918ec44SMatthew G. Knepley DMPlexShape shape = DM_SHAPE_UNKNOWN; 177f918ec44SMatthew G. Knepley 178f918ec44SMatthew G. Knepley PetscFunctionBeginUser; 179d0609cedSBarry Smith PetscOptionsBegin(comm, "", "libCEED Test Options", "DMPLEX"); 180d0609cedSBarry Smith PetscOptionsEnd(); 1819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-dm_plex_shape", DMPlexShapes, (PetscEnum *) &shape, NULL)); 18230602db0SMatthew G. Knepley ctx->setupgeo = NULL; 18330602db0SMatthew G. Knepley ctx->setupgeofname = NULL; 18430602db0SMatthew G. Knepley ctx->apply = Mass; 18530602db0SMatthew G. Knepley ctx->applyfname = Mass_loc; 18630602db0SMatthew G. Knepley ctx->areaExact = 0.0; 187f918ec44SMatthew G. Knepley switch (shape) { 188f918ec44SMatthew G. Knepley case DM_SHAPE_BOX_SURFACE: 189f918ec44SMatthew G. Knepley ctx->setupgeo = SetupMassGeoCube; 190f918ec44SMatthew G. Knepley ctx->setupgeofname = SetupMassGeoCube_loc; 191f918ec44SMatthew G. Knepley ctx->areaExact = 6.0; 192f918ec44SMatthew G. Knepley break; 193f918ec44SMatthew G. Knepley case DM_SHAPE_SPHERE: 194f918ec44SMatthew G. Knepley ctx->setupgeo = SetupMassGeoSphere; 195f918ec44SMatthew G. Knepley ctx->setupgeofname = SetupMassGeoSphere_loc; 196f918ec44SMatthew G. Knepley ctx->areaExact = 4.0*M_PI; 197f918ec44SMatthew G. Knepley break; 19830602db0SMatthew G. Knepley default: break; 199f918ec44SMatthew G. Knepley } 200f918ec44SMatthew G. Knepley PetscFunctionReturn(0); 201f918ec44SMatthew G. Knepley } 202f918ec44SMatthew G. Knepley 203f918ec44SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 204f918ec44SMatthew G. Knepley { 205f918ec44SMatthew G. Knepley PetscFunctionBegin; 2069566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2079566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2089566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2099566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 21030602db0SMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 21130602db0SMatthew G. Knepley { 21230602db0SMatthew G. Knepley Ceed ceed; 21330602db0SMatthew G. Knepley const char *usedresource; 21430602db0SMatthew G. Knepley 2159566063dSJacob Faibussowitsch PetscCall(DMGetCeed(*dm, &ceed)); 2169566063dSJacob Faibussowitsch PetscCall(CeedGetResource(ceed, &usedresource)); 2179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) *dm), "libCEED Backend: %s\n", usedresource)); 21830602db0SMatthew G. Knepley } 21930602db0SMatthew G. Knepley #endif 220f918ec44SMatthew G. Knepley PetscFunctionReturn(0); 221f918ec44SMatthew G. Knepley } 222f918ec44SMatthew G. Knepley 223f918ec44SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm) 224f918ec44SMatthew G. Knepley { 225f918ec44SMatthew G. Knepley DM cdm; 226a2c9b50fSJeremy L Thompson PetscFE fe, cfe; 227a2c9b50fSJeremy L Thompson PetscInt dim, cnc; 228f918ec44SMatthew G. Knepley PetscBool simplex; 229f918ec44SMatthew G. Knepley 230f918ec44SMatthew G. Knepley PetscFunctionBeginUser; 2319566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2329566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2339566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe)); 2349566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "indicator")); 2359566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, NULL, (PetscObject) fe)); 2369566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2379566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2389566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 2399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cnc)); 2409566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, cnc, simplex, NULL, PETSC_DETERMINE, &cfe)); 2419566063dSJacob Faibussowitsch PetscCall(DMProjectCoordinates(dm, cfe)); 2429566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&cfe)); 2439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 2449566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL)); 245f918ec44SMatthew G. Knepley PetscFunctionReturn(0); 246f918ec44SMatthew G. Knepley } 247f918ec44SMatthew G. Knepley 248f918ec44SMatthew G. Knepley static PetscErrorCode LibCeedSetupByDegree(DM dm, AppCtx *ctx, CeedData *data) 249f918ec44SMatthew G. Knepley { 250f918ec44SMatthew G. Knepley PetscDS ds; 251f918ec44SMatthew G. Knepley PetscFE fe, cfe; 252f918ec44SMatthew G. Knepley Ceed ceed; 253f918ec44SMatthew G. Knepley CeedElemRestriction Erestrictx, Erestrictu, Erestrictq; 254f918ec44SMatthew G. Knepley CeedQFunction qf_setupgeo; 255f918ec44SMatthew G. Knepley CeedOperator op_setupgeo; 256f918ec44SMatthew G. Knepley CeedVector xcoord; 257f918ec44SMatthew G. Knepley CeedBasis basisu, basisx; 258f918ec44SMatthew G. Knepley CeedInt Nqdata = 1; 259f918ec44SMatthew G. Knepley CeedInt nqpts, nqptsx; 260f918ec44SMatthew G. Knepley DM cdm; 261f918ec44SMatthew G. Knepley Vec coords; 262f918ec44SMatthew G. Knepley const PetscScalar *coordArray; 263f918ec44SMatthew G. Knepley PetscInt dim, cdim, cStart, cEnd, Ncell; 264f918ec44SMatthew G. Knepley 265f918ec44SMatthew G. Knepley PetscFunctionBeginUser; 2669566063dSJacob Faibussowitsch PetscCall(DMGetCeed(dm, &ceed)); 2679566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2689566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 2699566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 270f918ec44SMatthew G. Knepley Ncell = cEnd - cStart; 271f918ec44SMatthew G. Knepley // CEED bases 2729566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2739566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *) &fe)); 2749566063dSJacob Faibussowitsch PetscCall(PetscFEGetCeedBasis(fe, &basisu)); 2759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 2769566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &ds)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *) &cfe)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFEGetCeedBasis(cfe, &basisx)); 279f918ec44SMatthew G. Knepley 2809566063dSJacob Faibussowitsch PetscCall(DMPlexGetCeedRestriction(cdm, NULL, 0, 0, 0, &Erestrictx)); 2819566063dSJacob Faibussowitsch PetscCall(DMPlexGetCeedRestriction(dm, NULL, 0, 0, 0, &Erestrictu)); 2829566063dSJacob Faibussowitsch PetscCall(CeedBasisGetNumQuadraturePoints(basisu, &nqpts)); 2839566063dSJacob Faibussowitsch PetscCall(CeedBasisGetNumQuadraturePoints(basisx, &nqptsx)); 28463a3b9bcSJacob Faibussowitsch 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); 2859566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateStrided(ceed, Ncell, nqpts, Nqdata, Nqdata*Ncell*nqpts, CEED_STRIDES_BACKEND, &Erestrictq)); 286f918ec44SMatthew G. Knepley 2879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 2889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords, &coordArray)); 2899566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictx, &xcoord, NULL)); 2909566063dSJacob Faibussowitsch PetscCall(CeedVectorSetArray(xcoord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *) coordArray)); 2919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords, &coordArray)); 292f918ec44SMatthew G. Knepley 293f918ec44SMatthew G. Knepley // Create the vectors that will be needed in setup and apply 2949566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->uceed, NULL)); 2959566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictu, &data->vceed, NULL)); 2969566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionCreateVector(Erestrictq, &data->qdata, NULL)); 297f918ec44SMatthew G. Knepley 298f918ec44SMatthew G. Knepley // Create the Q-function that builds the operator (i.e. computes its quadrature data) and set its context data 2999566063dSJacob Faibussowitsch PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->setupgeo, ctx->setupgeofname, &qf_setupgeo)); 3009566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(qf_setupgeo, "x", cdim, CEED_EVAL_INTERP)); 3019566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(qf_setupgeo, "dx", cdim*dim, CEED_EVAL_GRAD)); 3029566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(qf_setupgeo, "weight", 1, CEED_EVAL_WEIGHT)); 3039566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddOutput(qf_setupgeo, "qdata", Nqdata, CEED_EVAL_NONE)); 304f918ec44SMatthew G. Knepley 305f918ec44SMatthew G. Knepley // Set up the mass operator 3069566063dSJacob Faibussowitsch PetscCall(CeedQFunctionCreateInterior(ceed, 1, ctx->apply, ctx->applyfname, &data->qf_apply)); 3079566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(data->qf_apply, "u", 1, CEED_EVAL_INTERP)); 3089566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddInput(data->qf_apply, "qdata", Nqdata, CEED_EVAL_NONE)); 3099566063dSJacob Faibussowitsch PetscCall(CeedQFunctionAddOutput(data->qf_apply, "v", 1, CEED_EVAL_INTERP)); 310f918ec44SMatthew G. Knepley 311f918ec44SMatthew G. Knepley // Create the operator that builds the quadrature data for the operator 3129566063dSJacob Faibussowitsch PetscCall(CeedOperatorCreate(ceed, qf_setupgeo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setupgeo)); 3139566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(op_setupgeo, "x", Erestrictx, basisx, CEED_VECTOR_ACTIVE)); 3149566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(op_setupgeo, "dx", Erestrictx, basisx, CEED_VECTOR_ACTIVE)); 3159566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(op_setupgeo, "weight", CEED_ELEMRESTRICTION_NONE, basisx, CEED_VECTOR_NONE)); 3169566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(op_setupgeo, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 317f918ec44SMatthew G. Knepley 318f918ec44SMatthew G. Knepley // Create the mass operator 3199566063dSJacob Faibussowitsch PetscCall(CeedOperatorCreate(ceed, data->qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &data->op_apply)); 3209566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(data->op_apply, "u", Erestrictu, basisu, CEED_VECTOR_ACTIVE)); 3219566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(data->op_apply, "qdata", Erestrictq, CEED_BASIS_COLLOCATED, data->qdata)); 3229566063dSJacob Faibussowitsch PetscCall(CeedOperatorSetField(data->op_apply, "v", Erestrictu, basisu, CEED_VECTOR_ACTIVE)); 323f918ec44SMatthew G. Knepley 324f918ec44SMatthew G. Knepley // Setup qdata 3259566063dSJacob Faibussowitsch PetscCall(CeedOperatorApply(op_setupgeo, xcoord, data->qdata, CEED_REQUEST_IMMEDIATE)); 326f918ec44SMatthew G. Knepley 3279566063dSJacob Faibussowitsch PetscCall(CeedElemRestrictionDestroy(&Erestrictq)); 3289566063dSJacob Faibussowitsch PetscCall(CeedQFunctionDestroy(&qf_setupgeo)); 3299566063dSJacob Faibussowitsch PetscCall(CeedOperatorDestroy(&op_setupgeo)); 3309566063dSJacob Faibussowitsch PetscCall(CeedVectorDestroy(&xcoord)); 331f918ec44SMatthew G. Knepley PetscFunctionReturn(0); 332f918ec44SMatthew G. Knepley } 333f918ec44SMatthew G. Knepley 334f918ec44SMatthew G. Knepley int main(int argc, char **argv) 335f918ec44SMatthew G. Knepley { 336f918ec44SMatthew G. Knepley MPI_Comm comm; 337f918ec44SMatthew G. Knepley DM dm; 338f918ec44SMatthew G. Knepley AppCtx ctx; 339f918ec44SMatthew G. Knepley Vec U, Uloc, V, Vloc; 340f918ec44SMatthew G. Knepley PetscScalar *v; 341f918ec44SMatthew G. Knepley PetscScalar area; 342f918ec44SMatthew G. Knepley CeedData ceeddata; 343f918ec44SMatthew G. Knepley 344*327415f7SBarry Smith PetscFunctionBeginUser; 3459566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 346f918ec44SMatthew G. Knepley comm = PETSC_COMM_WORLD; 3479566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &ctx)); 3489566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &ctx, &dm)); 3499566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm)); 350f918ec44SMatthew G. Knepley 3519566063dSJacob Faibussowitsch PetscCall(LibCeedSetupByDegree(dm, &ctx, &ceeddata)); 352f918ec44SMatthew G. Knepley 3539566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &U)); 3549566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, &Uloc)); 3559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &V)); 3569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Uloc, &Vloc)); 357f918ec44SMatthew G. Knepley 35830602db0SMatthew G. Knepley /**/ 3599566063dSJacob Faibussowitsch PetscCall(VecSet(Uloc, 1.)); 3609566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(V)); 3619566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Vloc)); 3629566063dSJacob Faibussowitsch PetscCall(VecGetArray(Vloc, &v)); 3639566063dSJacob Faibussowitsch PetscCall(CeedVectorSetArray(ceeddata.vceed, CEED_MEM_HOST, CEED_USE_POINTER, v)); 3649566063dSJacob Faibussowitsch PetscCall(CeedVectorSetValue(ceeddata.uceed, 1.0)); 3659566063dSJacob Faibussowitsch PetscCall(CeedOperatorApply(ceeddata.op_apply, ceeddata.uceed, ceeddata.vceed, CEED_REQUEST_IMMEDIATE)); 3669566063dSJacob Faibussowitsch PetscCall(CeedVectorTakeArray(ceeddata.vceed, CEED_MEM_HOST, NULL)); 3679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Vloc, &v)); 3689566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, Vloc, ADD_VALUES, V)); 3699566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, Vloc, ADD_VALUES, V)); 370f918ec44SMatthew G. Knepley 3719566063dSJacob Faibussowitsch PetscCall(VecSum(V, &area)); 372f918ec44SMatthew G. Knepley if (ctx.areaExact > 0.) { 373f918ec44SMatthew G. Knepley PetscReal error = PetscAbsReal(area - ctx.areaExact); 374f918ec44SMatthew G. Knepley PetscReal tol = PETSC_SMALL; 375f918ec44SMatthew G. Knepley 37663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Exact mesh surface area : % .*f\n", PetscAbsReal(ctx.areaExact - round(ctx.areaExact)) > 1E-15 ? 14 : 1, (double) ctx.areaExact)); 37763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Computed mesh surface area : % .*f\n", PetscAbsScalar(area - round(area)) > 1E-15 ? 14 : 1, (double)PetscRealPart(area))); 378f918ec44SMatthew G. Knepley if (error > tol) { 37963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Area error : % .14g\n", (double) error)); 380f918ec44SMatthew G. Knepley } else { 38163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Area verifies!\n")); 382f918ec44SMatthew G. Knepley } 383f918ec44SMatthew G. Knepley } 384f918ec44SMatthew G. Knepley 3859566063dSJacob Faibussowitsch PetscCall(CeedDataDestroy(&ceeddata)); 3869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Uloc)); 3889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V)); 3899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Vloc)); 3909566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 391f918ec44SMatthew G. Knepley return PetscFinalize(); 392f918ec44SMatthew G. Knepley } 393f918ec44SMatthew G. Knepley 394f918ec44SMatthew G. Knepley /*TEST 395f918ec44SMatthew G. Knepley 396f918ec44SMatthew G. Knepley build: 397f918ec44SMatthew G. Knepley requires: libceed 398f918ec44SMatthew G. Knepley 399f918ec44SMatthew G. Knepley testset: 400e600fa54SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -dm_view -dm_petscds_view \ 401f918ec44SMatthew G. Knepley -petscfe_default_quadrature_order 4 -coord_dm_default_quadrature_order 4 402f918ec44SMatthew G. Knepley 403f918ec44SMatthew G. Knepley test: 404f918ec44SMatthew G. Knepley suffix: cube_3 405f918ec44SMatthew G. Knepley args: -dm_plex_shape box_surface -dm_refine 2 406f918ec44SMatthew G. Knepley 407f918ec44SMatthew G. Knepley test: 408f918ec44SMatthew G. Knepley suffix: cube_3_p4 409f918ec44SMatthew G. Knepley nsize: 4 410fb796b39SJed Brown args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape box_surface -dm_refine 1 411f918ec44SMatthew G. Knepley 412f918ec44SMatthew G. Knepley test: 413f918ec44SMatthew G. Knepley suffix: sphere_3 414f918ec44SMatthew G. Knepley args: -dm_plex_shape sphere -dm_refine 3 415f918ec44SMatthew G. Knepley 416f918ec44SMatthew G. Knepley test: 417f918ec44SMatthew G. Knepley suffix: sphere_3_p4 418f918ec44SMatthew G. Knepley nsize: 4 419fb796b39SJed Brown args: -petscpartitioner_type simple -dm_refine_pre 1 -dm_plex_shape sphere -dm_refine 2 420f918ec44SMatthew G. Knepley 421f918ec44SMatthew G. Knepley TEST*/ 422