xref: /petsc/src/dm/impls/swarm/swarmpic_plex.c (revision 57508ece14a6b1339c0bbf016ecd72f673a062b0)
10e2ec84fSDave May #include <petscdm.h>
20e2ec84fSDave May #include <petscdmplex.h>
30e2ec84fSDave May #include <petscdmswarm.h>
4279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
50e2ec84fSDave May 
6f39ec787SMatthew G. Knepley PetscBool  SwarmProjcite       = PETSC_FALSE;
7f39ec787SMatthew G. Knepley const char SwarmProjCitation[] = "@article{PusztayKnepleyAdams2022,\n"
8f39ec787SMatthew G. Knepley                                  "title   = {Conservative Projection Between FEM and Particle Bases},\n"
9f39ec787SMatthew G. Knepley                                  "author  = {Joseph V. Pusztay and Matthew G. Knepley and Mark F. Adams},\n"
10f39ec787SMatthew G. Knepley                                  "journal = {SIAM Journal on Scientific Computing},\n"
11f39ec787SMatthew G. Knepley                                  "volume  = {44},\n"
12f39ec787SMatthew G. Knepley                                  "number  = {4},\n"
13f39ec787SMatthew G. Knepley                                  "pages   = {C310--C319},\n"
14f39ec787SMatthew G. Knepley                                  "doi     = {10.1137/21M145407},\n"
15f39ec787SMatthew G. Knepley                                  "year    = {2022}\n}\n";
16f39ec787SMatthew G. Knepley 
17b12d2206SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *xi);
18b12d2206SDave May 
19d71ae5a4SJacob Faibussowitsch static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem)
20d71ae5a4SJacob Faibussowitsch {
212cf64d79SDave May   const PetscInt  Nc = 1;
222cf64d79SDave May   PetscQuadrature q, fq;
232cf64d79SDave May   DM              K;
242cf64d79SDave May   PetscSpace      P;
252cf64d79SDave May   PetscDualSpace  Q;
262cf64d79SDave May   PetscInt        order, quadPointsPerEdge;
272cf64d79SDave May   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
282cf64d79SDave May 
292cf64d79SDave May   PetscFunctionBegin;
302cf64d79SDave May   /* Create space */
319566063dSJacob Faibussowitsch   PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)dm), &P));
329566063dSJacob Faibussowitsch   /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix)); */
339566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialSetTensor(P, tensor));
349566063dSJacob Faibussowitsch   /* PetscCall(PetscSpaceSetFromOptions(P)); */
359566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
369566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetDegree(P, 1, PETSC_DETERMINE));
379566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumComponents(P, Nc));
389566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumVariables(P, dim));
399566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetUp(P));
409566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(P, &order, NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialGetTensor(P, &tensor));
422cf64d79SDave May   /* Create dual space */
439566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)dm), &Q));
449566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
459566063dSJacob Faibussowitsch   /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix)); */
469566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
479566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(Q, K));
489566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&K));
499566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetNumComponents(Q, Nc));
509566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetOrder(Q, order));
519566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetTensor(Q, tensor));
529566063dSJacob Faibussowitsch   /* PetscCall(PetscDualSpaceSetFromOptions(Q)); */
539566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
549566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(Q));
552cf64d79SDave May   /* Create element */
569566063dSJacob Faibussowitsch   PetscCall(PetscFECreate(PetscObjectComm((PetscObject)dm), fem));
579566063dSJacob Faibussowitsch   /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix)); */
589566063dSJacob Faibussowitsch   /* PetscCall(PetscFESetFromOptions(*fem)); */
599566063dSJacob Faibussowitsch   PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
609566063dSJacob Faibussowitsch   PetscCall(PetscFESetBasisSpace(*fem, P));
619566063dSJacob Faibussowitsch   PetscCall(PetscFESetDualSpace(*fem, Q));
629566063dSJacob Faibussowitsch   PetscCall(PetscFESetNumComponents(*fem, Nc));
639566063dSJacob Faibussowitsch   PetscCall(PetscFESetUp(*fem));
649566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&P));
659566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&Q));
662cf64d79SDave May   /* Create quadrature (with specified order if given) */
672cf64d79SDave May   qorder            = qorder >= 0 ? qorder : order;
682cf64d79SDave May   quadPointsPerEdge = PetscMax(qorder + 1, 1);
692cf64d79SDave May   if (isSimplex) {
709566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q));
719566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq));
729371c9d4SSatish Balay   } else {
739566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q));
749566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq));
752cf64d79SDave May   }
769566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(*fem, q));
779566063dSJacob Faibussowitsch   PetscCall(PetscFESetFaceQuadrature(*fem, fq));
789566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&q));
799566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&fq));
803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
812cf64d79SDave May }
822cf64d79SDave May 
83ba38deedSJacob Faibussowitsch static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm, DM dmc, PetscInt nsub)
84d71ae5a4SJacob Faibussowitsch {
8578185223SDave May   PetscInt         dim, nfaces, nbasis;
8678185223SDave May   PetscInt         q, npoints_q, e, nel, pcnt, ps, pe, d, k, r;
87ef0bb6c7SMatthew G. Knepley   PetscTabulation  T;
8878185223SDave May   Vec              coorlocal;
8978185223SDave May   PetscSection     coordSection;
9078185223SDave May   PetscScalar     *elcoor = NULL;
9178185223SDave May   PetscReal       *swarm_coor;
9278185223SDave May   PetscInt        *swarm_cellid;
9378185223SDave May   const PetscReal *xiq;
9478185223SDave May   PetscQuadrature  quadrature;
9578185223SDave May   PetscFE          fe, feRef;
9678185223SDave May   PetscBool        is_simplex;
9778185223SDave May 
9878185223SDave May   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmc, &dim));
10078185223SDave May   is_simplex = PETSC_FALSE;
1019566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
1029566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
103ad540459SPierre Jolivet   if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
10478185223SDave May 
1059566063dSJacob Faibussowitsch   PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe));
10678185223SDave May 
10778185223SDave May   for (r = 0; r < nsub; r++) {
1089566063dSJacob Faibussowitsch     PetscCall(PetscFERefine(fe, &feRef));
1099566063dSJacob Faibussowitsch     PetscCall(PetscFECopyQuadrature(feRef, fe));
1109566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&feRef));
11178185223SDave May   }
11278185223SDave May 
1139566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
1149566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL));
1159566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDimension(fe, &nbasis));
1169566063dSJacob Faibussowitsch   PetscCall(PetscFEGetCellTabulation(fe, 1, &T));
11778185223SDave May 
11878185223SDave May   /* 0->cell, 1->edge, 2->vert */
1199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
12078185223SDave May   nel = pe - ps;
12178185223SDave May 
1229566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
1239566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
1249566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
12578185223SDave May 
1269566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
1279566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dmc, &coordSection));
12878185223SDave May 
12978185223SDave May   pcnt = 0;
13078185223SDave May   for (e = 0; e < nel; e++) {
1319566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
13278185223SDave May 
13378185223SDave May     for (q = 0; q < npoints_q; q++) {
13478185223SDave May       for (d = 0; d < dim; d++) {
13578185223SDave May         swarm_coor[dim * pcnt + d] = 0.0;
136ad540459SPierre Jolivet         for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][q * nbasis + k] * PetscRealPart(elcoor[dim * k + d]);
13778185223SDave May       }
13878185223SDave May       swarm_cellid[pcnt] = e;
13978185223SDave May       pcnt++;
14078185223SDave May     }
1419566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
14278185223SDave May   }
1439566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
1449566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
14578185223SDave May 
1469566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14878185223SDave May }
14978185223SDave May 
150ba38deedSJacob Faibussowitsch static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm, DM dmc, PetscInt npoints)
151d71ae5a4SJacob Faibussowitsch {
152bc53056eSDave May   PetscInt     dim;
153bc53056eSDave May   PetscInt     ii, jj, q, npoints_q, e, nel, npe, pcnt, ps, pe, d, k, nfaces;
1540e2ec84fSDave May   PetscReal   *xi, ds, ds2;
1550e2ec84fSDave May   PetscReal  **basis;
1560e2ec84fSDave May   Vec          coorlocal;
1570e2ec84fSDave May   PetscSection coordSection;
1580e2ec84fSDave May   PetscScalar *elcoor = NULL;
1590e2ec84fSDave May   PetscReal   *swarm_coor;
1600e2ec84fSDave May   PetscInt    *swarm_cellid;
161bc53056eSDave May   PetscBool    is_simplex;
1620e2ec84fSDave May 
1630e2ec84fSDave May   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmc, &dim));
16508401ef6SPierre Jolivet   PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only 2D is supported");
166bc53056eSDave May   is_simplex = PETSC_FALSE;
1679566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
1689566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
169ad540459SPierre Jolivet   if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
17028b400f6SJacob Faibussowitsch   PetscCheck(is_simplex, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only the simplex is supported");
171bc53056eSDave May 
1729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * npoints * npoints, &xi));
1730e2ec84fSDave May   pcnt = 0;
174*57508eceSPierre Jolivet   ds   = 1.0 / (PetscReal)(npoints - 1);
175*57508eceSPierre Jolivet   ds2  = 1.0 / (PetscReal)npoints;
1760e2ec84fSDave May   for (jj = 0; jj < npoints; jj++) {
1770e2ec84fSDave May     for (ii = 0; ii < npoints - jj; ii++) {
1780e2ec84fSDave May       xi[dim * pcnt + 0] = ii * ds;
1790e2ec84fSDave May       xi[dim * pcnt + 1] = jj * ds;
1800e2ec84fSDave May 
1810e2ec84fSDave May       xi[dim * pcnt + 0] *= (1.0 - 1.2 * ds2);
1820e2ec84fSDave May       xi[dim * pcnt + 1] *= (1.0 - 1.2 * ds2);
1830e2ec84fSDave May 
1840e2ec84fSDave May       xi[dim * pcnt + 0] += 0.35 * ds2;
1850e2ec84fSDave May       xi[dim * pcnt + 1] += 0.35 * ds2;
1860e2ec84fSDave May       pcnt++;
1870e2ec84fSDave May     }
1880e2ec84fSDave May   }
1890e2ec84fSDave May   npoints_q = pcnt;
1900e2ec84fSDave May 
1910e2ec84fSDave May   npe = 3; /* nodes per element (triangle) */
1929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints_q, &basis));
1930e2ec84fSDave May   for (q = 0; q < npoints_q; q++) {
1949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npe, &basis[q]));
1950e2ec84fSDave May 
1960e2ec84fSDave May     basis[q][0] = 1.0 - xi[dim * q + 0] - xi[dim * q + 1];
1970e2ec84fSDave May     basis[q][1] = xi[dim * q + 0];
1980e2ec84fSDave May     basis[q][2] = xi[dim * q + 1];
1990e2ec84fSDave May   }
2000e2ec84fSDave May 
2010e2ec84fSDave May   /* 0->cell, 1->edge, 2->vert */
2029566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
2030e2ec84fSDave May   nel = pe - ps;
2040e2ec84fSDave May 
2059566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
2069566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
2079566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
2080e2ec84fSDave May 
2099566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
2109566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dmc, &coordSection));
2110e2ec84fSDave May 
2120e2ec84fSDave May   pcnt = 0;
2130e2ec84fSDave May   for (e = 0; e < nel; e++) {
2149566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor));
2150e2ec84fSDave May 
2160e2ec84fSDave May     for (q = 0; q < npoints_q; q++) {
2170e2ec84fSDave May       for (d = 0; d < dim; d++) {
2180e2ec84fSDave May         swarm_coor[dim * pcnt + d] = 0.0;
219ad540459SPierre Jolivet         for (k = 0; k < npe; k++) swarm_coor[dim * pcnt + d] += basis[q][k] * PetscRealPart(elcoor[dim * k + d]);
2200e2ec84fSDave May       }
2210e2ec84fSDave May       swarm_cellid[pcnt] = e;
2220e2ec84fSDave May       pcnt++;
2230e2ec84fSDave May     }
2249566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor));
2250e2ec84fSDave May   }
2269566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
2279566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
2280e2ec84fSDave May 
2299566063dSJacob Faibussowitsch   PetscCall(PetscFree(xi));
23048a46eb9SPierre Jolivet   for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q]));
2319566063dSJacob Faibussowitsch   PetscCall(PetscFree(basis));
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2330e2ec84fSDave May }
2340e2ec84fSDave May 
235d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param)
236d71ae5a4SJacob Faibussowitsch {
2370e2ec84fSDave May   PetscInt dim;
2380e2ec84fSDave May 
2390e2ec84fSDave May   PetscFunctionBegin;
2409566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(celldm, &dim));
2410e2ec84fSDave May   switch (layout) {
2420e2ec84fSDave May   case DMSWARMPIC_LAYOUT_REGULAR:
24308401ef6SPierre Jolivet     PetscCheck(dim != 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No 3D support for REGULAR+PLEX");
2449566063dSJacob Faibussowitsch     PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm, celldm, layout_param));
2450e2ec84fSDave May     break;
2469371c9d4SSatish Balay   case DMSWARMPIC_LAYOUT_GAUSS: {
247013f9f9eSMatthew G. Knepley     PetscQuadrature  quad, facequad;
248b12d2206SDave May     const PetscReal *xi;
249013f9f9eSMatthew G. Knepley     DMPolytopeType   ct;
250013f9f9eSMatthew G. Knepley     PetscInt         cStart, Nq;
251b12d2206SDave May 
252013f9f9eSMatthew G. Knepley     PetscCall(DMPlexGetHeightStratum(celldm, 0, &cStart, NULL));
253013f9f9eSMatthew G. Knepley     PetscCall(DMPlexGetCellType(celldm, cStart, &ct));
254013f9f9eSMatthew G. Knepley     PetscCall(PetscDTCreateDefaultQuadrature(ct, layout_param, &quad, &facequad));
255013f9f9eSMatthew G. Knepley     PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xi, NULL));
256013f9f9eSMatthew G. Knepley     PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, Nq, (PetscReal *)xi));
257013f9f9eSMatthew G. Knepley     PetscCall(PetscQuadratureDestroy(&quad));
258013f9f9eSMatthew G. Knepley     PetscCall(PetscQuadratureDestroy(&facequad));
2599371c9d4SSatish Balay   } break;
260d71ae5a4SJacob Faibussowitsch   case DMSWARMPIC_LAYOUT_SUBDIVISION:
261d71ae5a4SJacob Faibussowitsch     PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm, celldm, layout_param));
262d71ae5a4SJacob Faibussowitsch     break;
2630e2ec84fSDave May   }
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2650e2ec84fSDave May }
266e3cd5995SDave May 
267d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm, DM dmc, PetscInt npoints, PetscReal xi[])
268d71ae5a4SJacob Faibussowitsch {
269431382f2SDave May   PetscBool       is_simplex, is_tensorcell;
27092e40656SDave May   PetscInt        dim, nfaces, ps, pe, p, d, nbasis, pcnt, e, k, nel;
27192e40656SDave May   PetscFE         fe;
27292e40656SDave May   PetscQuadrature quadrature;
273ef0bb6c7SMatthew G. Knepley   PetscTabulation T;
274ef0bb6c7SMatthew G. Knepley   PetscReal      *xiq;
27592e40656SDave May   Vec             coorlocal;
27692e40656SDave May   PetscSection    coordSection;
27792e40656SDave May   PetscScalar    *elcoor = NULL;
27892e40656SDave May   PetscReal      *swarm_coor;
27992e40656SDave May   PetscInt       *swarm_cellid;
280431382f2SDave May 
28192e40656SDave May   PetscFunctionBegin;
2829566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmc, &dim));
283431382f2SDave May 
284431382f2SDave May   is_simplex    = PETSC_FALSE;
285431382f2SDave May   is_tensorcell = PETSC_FALSE;
2869566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
2879566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
288431382f2SDave May 
289ad540459SPierre Jolivet   if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
290431382f2SDave May 
291431382f2SDave May   switch (dim) {
292431382f2SDave May   case 2:
293ad540459SPierre Jolivet     if (nfaces == 4) is_tensorcell = PETSC_TRUE;
294431382f2SDave May     break;
295431382f2SDave May   case 3:
296ad540459SPierre Jolivet     if (nfaces == 6) is_tensorcell = PETSC_TRUE;
297431382f2SDave May     break;
298d71ae5a4SJacob Faibussowitsch   default:
299d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for 2D, 3D");
300431382f2SDave May   }
301431382f2SDave May 
30292e40656SDave May   /* check points provided fail inside the reference cell */
303431382f2SDave May   if (is_simplex) {
30492e40656SDave May     for (p = 0; p < npoints; p++) {
30592e40656SDave May       PetscReal sum;
306ad540459SPierre Jolivet       for (d = 0; d < dim; d++) PetscCheck(xi[dim * p + d] >= -1.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain");
30792e40656SDave May       sum = 0.0;
308ad540459SPierre Jolivet       for (d = 0; d < dim; d++) sum += xi[dim * p + d];
30908401ef6SPierre Jolivet       PetscCheck(sum <= 0.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain");
31092e40656SDave May     }
311431382f2SDave May   } else if (is_tensorcell) {
31292e40656SDave May     for (p = 0; p < npoints; p++) {
313ad540459SPierre Jolivet       for (d = 0; d < dim; d++) PetscCheck(PetscAbsReal(xi[dim * p + d]) <= 1.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the tensor domain [-1,1]^d");
31492e40656SDave May     }
315431382f2SDave May   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for d-simplex and d-tensorcell");
316431382f2SDave May 
3179566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm), &quadrature));
3189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints * dim, &xiq));
3199566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(xiq, xi, npoints * dim));
3209566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(quadrature, dim, 1, npoints, (const PetscReal *)xiq, NULL));
3219566063dSJacob Faibussowitsch   PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe));
3229566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(fe, quadrature));
3239566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDimension(fe, &nbasis));
3249566063dSJacob Faibussowitsch   PetscCall(PetscFEGetCellTabulation(fe, 1, &T));
32592e40656SDave May 
32692e40656SDave May   /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */
32792e40656SDave May   /* 0->cell, 1->edge, 2->vert */
3289566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
32992e40656SDave May   nel = pe - ps;
33092e40656SDave May 
3319566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, npoints * nel, -1));
3329566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
3339566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
33492e40656SDave May 
3359566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
3369566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dmc, &coordSection));
33792e40656SDave May 
33892e40656SDave May   pcnt = 0;
33992e40656SDave May   for (e = 0; e < nel; e++) {
3409566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
34192e40656SDave May 
34292e40656SDave May     for (p = 0; p < npoints; p++) {
34392e40656SDave May       for (d = 0; d < dim; d++) {
34492e40656SDave May         swarm_coor[dim * pcnt + d] = 0.0;
345ad540459SPierre Jolivet         for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][p * nbasis + k] * PetscRealPart(elcoor[dim * k + d]);
34692e40656SDave May       }
34792e40656SDave May       swarm_cellid[pcnt] = e;
34892e40656SDave May       pcnt++;
34992e40656SDave May     }
3509566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
35192e40656SDave May   }
3529566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
3539566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
35492e40656SDave May 
3559566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quadrature));
3569566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
358431382f2SDave May }
359