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