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
private_PetscFECreateDefault_scalar_pk1(DM dm,PetscInt dim,PetscBool isSimplex,PetscInt qorder,PetscFE * fem)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
private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm,DM dmc,PetscInt nsub)83ba38deedSJacob Faibussowitsch static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm, DM dmc, PetscInt nsub)
84d71ae5a4SJacob Faibussowitsch {
8578185223SDave May PetscInt dim, nfaces, nbasis;
8619307e5cSMatthew G. Knepley PetscInt q, npoints_q, e, nel, pcnt, ps, pe, d, k, r, Nfc;
8719307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
88ef0bb6c7SMatthew G. Knepley PetscTabulation T;
8978185223SDave May Vec coorlocal;
9078185223SDave May PetscSection coordSection;
9178185223SDave May PetscScalar *elcoor = NULL;
9278185223SDave May PetscReal *swarm_coor;
9378185223SDave May PetscInt *swarm_cellid;
9478185223SDave May const PetscReal *xiq;
9578185223SDave May PetscQuadrature quadrature;
9678185223SDave May PetscFE fe, feRef;
9778185223SDave May PetscBool is_simplex;
9819307e5cSMatthew G. Knepley const char **coordFields, *cellid;
9978185223SDave May
10078185223SDave May PetscFunctionBegin;
1019566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim));
10278185223SDave May is_simplex = PETSC_FALSE;
1039566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
1049566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
105ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
10678185223SDave May
1079566063dSJacob Faibussowitsch PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe));
10878185223SDave May
10978185223SDave May for (r = 0; r < nsub; r++) {
1109566063dSJacob Faibussowitsch PetscCall(PetscFERefine(fe, &feRef));
1119566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(feRef, fe));
1129566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feRef));
11378185223SDave May }
11478185223SDave May
1159566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quadrature));
1169566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL));
1179566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fe, &nbasis));
1189566063dSJacob Faibussowitsch PetscCall(PetscFEGetCellTabulation(fe, 1, &T));
11978185223SDave May
12078185223SDave May /* 0->cell, 1->edge, 2->vert */
1219566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
12278185223SDave May nel = pe - ps;
12378185223SDave May
12419307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(dmc, &celldm));
12519307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
12619307e5cSMatthew G. Knepley PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
12719307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
12819307e5cSMatthew G. Knepley
1299566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
13019307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
13119307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
13278185223SDave May
1339566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
1349566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection));
13578185223SDave May
13678185223SDave May pcnt = 0;
13778185223SDave May for (e = 0; e < nel; e++) {
1389566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
13978185223SDave May
14078185223SDave May for (q = 0; q < npoints_q; q++) {
14178185223SDave May for (d = 0; d < dim; d++) {
14278185223SDave May swarm_coor[dim * pcnt + d] = 0.0;
143ad540459SPierre Jolivet for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][q * nbasis + k] * PetscRealPart(elcoor[dim * k + d]);
14478185223SDave May }
14578185223SDave May swarm_cellid[pcnt] = e;
14678185223SDave May pcnt++;
14778185223SDave May }
1489566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
14978185223SDave May }
15019307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
15119307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
15278185223SDave May
1539566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
15578185223SDave May }
15678185223SDave May
private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints)157ba38deedSJacob Faibussowitsch static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm, DM dmc, PetscInt npoints)
158d71ae5a4SJacob Faibussowitsch {
159bc53056eSDave May PetscInt dim;
16019307e5cSMatthew G. Knepley PetscInt ii, jj, q, npoints_q, e, nel, npe, pcnt, ps, pe, d, k, nfaces, Nfc;
1610e2ec84fSDave May PetscReal *xi, ds, ds2;
1620e2ec84fSDave May PetscReal **basis;
16319307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
1640e2ec84fSDave May Vec coorlocal;
1650e2ec84fSDave May PetscSection coordSection;
1660e2ec84fSDave May PetscScalar *elcoor = NULL;
1670e2ec84fSDave May PetscReal *swarm_coor;
1680e2ec84fSDave May PetscInt *swarm_cellid;
169bc53056eSDave May PetscBool is_simplex;
17019307e5cSMatthew G. Knepley const char **coordFields, *cellid;
1710e2ec84fSDave May
1720e2ec84fSDave May PetscFunctionBegin;
1739566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim));
17408401ef6SPierre Jolivet PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only 2D is supported");
175bc53056eSDave May is_simplex = PETSC_FALSE;
1769566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
1779566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
178ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
17928b400f6SJacob Faibussowitsch PetscCheck(is_simplex, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only the simplex is supported");
180bc53056eSDave May
1819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints * npoints, &xi));
1820e2ec84fSDave May pcnt = 0;
18357508eceSPierre Jolivet ds = 1.0 / (PetscReal)(npoints - 1);
18457508eceSPierre Jolivet ds2 = 1.0 / (PetscReal)npoints;
1850e2ec84fSDave May for (jj = 0; jj < npoints; jj++) {
1860e2ec84fSDave May for (ii = 0; ii < npoints - jj; ii++) {
1870e2ec84fSDave May xi[dim * pcnt + 0] = ii * ds;
1880e2ec84fSDave May xi[dim * pcnt + 1] = jj * ds;
1890e2ec84fSDave May
1900e2ec84fSDave May xi[dim * pcnt + 0] *= (1.0 - 1.2 * ds2);
1910e2ec84fSDave May xi[dim * pcnt + 1] *= (1.0 - 1.2 * ds2);
1920e2ec84fSDave May
1930e2ec84fSDave May xi[dim * pcnt + 0] += 0.35 * ds2;
1940e2ec84fSDave May xi[dim * pcnt + 1] += 0.35 * ds2;
1950e2ec84fSDave May pcnt++;
1960e2ec84fSDave May }
1970e2ec84fSDave May }
1980e2ec84fSDave May npoints_q = pcnt;
1990e2ec84fSDave May
2000e2ec84fSDave May npe = 3; /* nodes per element (triangle) */
2019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints_q, &basis));
2020e2ec84fSDave May for (q = 0; q < npoints_q; q++) {
2039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npe, &basis[q]));
2040e2ec84fSDave May
2050e2ec84fSDave May basis[q][0] = 1.0 - xi[dim * q + 0] - xi[dim * q + 1];
2060e2ec84fSDave May basis[q][1] = xi[dim * q + 0];
2070e2ec84fSDave May basis[q][2] = xi[dim * q + 1];
2080e2ec84fSDave May }
2090e2ec84fSDave May
2100e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */
2119566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
2120e2ec84fSDave May nel = pe - ps;
2130e2ec84fSDave May
21419307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(dmc, &celldm));
21519307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
21619307e5cSMatthew G. Knepley PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
21719307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
21819307e5cSMatthew G. Knepley
2199566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
22019307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
22119307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
2220e2ec84fSDave May
2239566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
2249566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection));
2250e2ec84fSDave May
2260e2ec84fSDave May pcnt = 0;
2270e2ec84fSDave May for (e = 0; e < nel; e++) {
2289566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor));
2290e2ec84fSDave May
2300e2ec84fSDave May for (q = 0; q < npoints_q; q++) {
2310e2ec84fSDave May for (d = 0; d < dim; d++) {
2320e2ec84fSDave May swarm_coor[dim * pcnt + d] = 0.0;
233ad540459SPierre Jolivet for (k = 0; k < npe; k++) swarm_coor[dim * pcnt + d] += basis[q][k] * PetscRealPart(elcoor[dim * k + d]);
2340e2ec84fSDave May }
2350e2ec84fSDave May swarm_cellid[pcnt] = e;
2360e2ec84fSDave May pcnt++;
2370e2ec84fSDave May }
2389566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor));
2390e2ec84fSDave May }
24019307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
24119307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
2420e2ec84fSDave May
2439566063dSJacob Faibussowitsch PetscCall(PetscFree(xi));
24448a46eb9SPierre Jolivet for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q]));
2459566063dSJacob Faibussowitsch PetscCall(PetscFree(basis));
2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2470e2ec84fSDave May }
2480e2ec84fSDave May
private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param)249d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param)
250d71ae5a4SJacob Faibussowitsch {
2510e2ec84fSDave May PetscInt dim;
2520e2ec84fSDave May
2530e2ec84fSDave May PetscFunctionBegin;
2549566063dSJacob Faibussowitsch PetscCall(DMGetDimension(celldm, &dim));
2550e2ec84fSDave May switch (layout) {
2560e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR:
25708401ef6SPierre Jolivet PetscCheck(dim != 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No 3D support for REGULAR+PLEX");
2589566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm, celldm, layout_param));
2590e2ec84fSDave May break;
2609371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_GAUSS: {
261013f9f9eSMatthew G. Knepley PetscQuadrature quad, facequad;
262b12d2206SDave May const PetscReal *xi;
263013f9f9eSMatthew G. Knepley DMPolytopeType ct;
264013f9f9eSMatthew G. Knepley PetscInt cStart, Nq;
265b12d2206SDave May
266013f9f9eSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(celldm, 0, &cStart, NULL));
267013f9f9eSMatthew G. Knepley PetscCall(DMPlexGetCellType(celldm, cStart, &ct));
268013f9f9eSMatthew G. Knepley PetscCall(PetscDTCreateDefaultQuadrature(ct, layout_param, &quad, &facequad));
269013f9f9eSMatthew G. Knepley PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xi, NULL));
270013f9f9eSMatthew G. Knepley PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, Nq, (PetscReal *)xi));
271013f9f9eSMatthew G. Knepley PetscCall(PetscQuadratureDestroy(&quad));
272013f9f9eSMatthew G. Knepley PetscCall(PetscQuadratureDestroy(&facequad));
2739371c9d4SSatish Balay } break;
274d71ae5a4SJacob Faibussowitsch case DMSWARMPIC_LAYOUT_SUBDIVISION:
275d71ae5a4SJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm, celldm, layout_param));
276d71ae5a4SJacob Faibussowitsch break;
2770e2ec84fSDave May }
2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2790e2ec84fSDave May }
280e3cd5995SDave May
private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm,DM dmc,PetscInt npoints,PetscReal xi[])281d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm, DM dmc, PetscInt npoints, PetscReal xi[])
282d71ae5a4SJacob Faibussowitsch {
283431382f2SDave May PetscBool is_simplex, is_tensorcell;
284*5d7528f4SZach Atkins PetscInt dim, ps, pe, nel, nfaces, Nfc;
28519307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
28692e40656SDave May PetscReal *swarm_coor;
28792e40656SDave May PetscInt *swarm_cellid;
28819307e5cSMatthew G. Knepley const char **coordFields, *cellid;
289431382f2SDave May
29092e40656SDave May PetscFunctionBegin;
2919566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim));
292431382f2SDave May
293431382f2SDave May is_simplex = PETSC_FALSE;
294431382f2SDave May is_tensorcell = PETSC_FALSE;
2959566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
2969566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
297431382f2SDave May
298ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
299431382f2SDave May
300431382f2SDave May switch (dim) {
301431382f2SDave May case 2:
302ad540459SPierre Jolivet if (nfaces == 4) is_tensorcell = PETSC_TRUE;
303431382f2SDave May break;
304431382f2SDave May case 3:
305ad540459SPierre Jolivet if (nfaces == 6) is_tensorcell = PETSC_TRUE;
306431382f2SDave May break;
307d71ae5a4SJacob Faibussowitsch default:
308d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for 2D, 3D");
309431382f2SDave May }
310431382f2SDave May
31192e40656SDave May /* check points provided fail inside the reference cell */
312431382f2SDave May if (is_simplex) {
313*5d7528f4SZach Atkins for (PetscInt p = 0; p < npoints; p++) {
31492e40656SDave May PetscReal sum;
315*5d7528f4SZach Atkins for (PetscInt 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");
31692e40656SDave May sum = 0.0;
317*5d7528f4SZach Atkins for (PetscInt d = 0; d < dim; d++) sum += xi[dim * p + d];
31808401ef6SPierre Jolivet PetscCheck(sum <= 0.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain");
31992e40656SDave May }
320431382f2SDave May } else if (is_tensorcell) {
321*5d7528f4SZach Atkins for (PetscInt p = 0; p < npoints; p++) {
322*5d7528f4SZach Atkins for (PetscInt 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");
32392e40656SDave May }
324431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for d-simplex and d-tensorcell");
325431382f2SDave May
3269566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
32792e40656SDave May nel = pe - ps;
32892e40656SDave May
32919307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
33019307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
33119307e5cSMatthew G. Knepley PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
33219307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
33319307e5cSMatthew G. Knepley
334*5d7528f4SZach Atkins PetscCall(DMSwarmSetLocalSizes(dm, npoints * nel, PETSC_DECIDE));
33519307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
33619307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
33792e40656SDave May
338*5d7528f4SZach Atkins // Use DMPlexReferenceToCoordinates so that arbitrary discretizations work
339*5d7528f4SZach Atkins for (PetscInt e = 0; e < nel; e++) {
340*5d7528f4SZach Atkins PetscCall(DMPlexReferenceToCoordinates(dmc, e + ps, npoints, xi, &swarm_coor[npoints * dim * e]));
341*5d7528f4SZach Atkins for (PetscInt p = 0; p < npoints; p++) swarm_cellid[e * npoints + p] = e + ps;
34292e40656SDave May }
34319307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
34419307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
346431382f2SDave May }
347