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; 1740e2ec84fSDave May ds = 1.0 / ((PetscReal)(npoints - 1)); 1750e2ec84fSDave May 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: { 247*013f9f9eSMatthew G. Knepley PetscQuadrature quad, facequad; 248b12d2206SDave May const PetscReal *xi; 249*013f9f9eSMatthew G. Knepley DMPolytopeType ct; 250*013f9f9eSMatthew G. Knepley PetscInt cStart, Nq; 251b12d2206SDave May 252*013f9f9eSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(celldm, 0, &cStart, NULL)); 253*013f9f9eSMatthew G. Knepley PetscCall(DMPlexGetCellType(celldm, cStart, &ct)); 254*013f9f9eSMatthew G. Knepley PetscCall(PetscDTCreateDefaultQuadrature(ct, layout_param, &quad, &facequad)); 255*013f9f9eSMatthew G. Knepley PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xi, NULL)); 256*013f9f9eSMatthew G. Knepley PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, Nq, (PetscReal *)xi)); 257*013f9f9eSMatthew G. Knepley PetscCall(PetscQuadratureDestroy(&quad)); 258*013f9f9eSMatthew 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 267e3cd5995SDave May /* 268e3cd5995SDave May typedef struct { 269e3cd5995SDave May PetscReal x,y; 270e3cd5995SDave May } Point2d; 271e3cd5995SDave May 272e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s) 273e3cd5995SDave May { 274e3cd5995SDave May PetscFunctionBegin; 275e3cd5995SDave May *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y); 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277e3cd5995SDave May } 278e3cd5995SDave May */ 279e3cd5995SDave May /* 280e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v) 281e3cd5995SDave May { 282e3cd5995SDave May PetscReal s1,s2,s3; 283e3cd5995SDave May PetscBool b1, b2, b3; 284e3cd5995SDave May 285e3cd5995SDave May PetscFunctionBegin; 286e3cd5995SDave May signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f; 287e3cd5995SDave May signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f; 288e3cd5995SDave May signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f; 289e3cd5995SDave May 290e3cd5995SDave May *v = ((b1 == b2) && (b2 == b3)); 2913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 292e3cd5995SDave May } 293e3cd5995SDave May */ 294e3cd5995SDave May /* 295e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ) 296e3cd5995SDave May { 297e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 298e3cd5995SDave May PetscReal c,b[2],A[2][2],inv[2][2],detJ,od; 299e3cd5995SDave May 300e3cd5995SDave May PetscFunctionBegin; 301e3cd5995SDave May x1 = coords[2*0+0]; 302e3cd5995SDave May x2 = coords[2*1+0]; 303e3cd5995SDave May x3 = coords[2*2+0]; 304e3cd5995SDave May 305e3cd5995SDave May y1 = coords[2*0+1]; 306e3cd5995SDave May y2 = coords[2*1+1]; 307e3cd5995SDave May y3 = coords[2*2+1]; 308e3cd5995SDave May 309e3cd5995SDave May c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3; 310e3cd5995SDave May b[0] = xp[0] - c; 311e3cd5995SDave May c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3; 312e3cd5995SDave May b[1] = xp[1] - c; 313e3cd5995SDave May 314e3cd5995SDave May A[0][0] = -0.5*x1 + 0.5*x2; A[0][1] = -0.5*x1 + 0.5*x3; 315e3cd5995SDave May A[1][0] = -0.5*y1 + 0.5*y2; A[1][1] = -0.5*y1 + 0.5*y3; 316e3cd5995SDave May 317e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 318e3cd5995SDave May *dJ = PetscAbsReal(detJ); 319e3cd5995SDave May od = 1.0/detJ; 320e3cd5995SDave May 321e3cd5995SDave May inv[0][0] = A[1][1] * od; 322e3cd5995SDave May inv[0][1] = -A[0][1] * od; 323e3cd5995SDave May inv[1][0] = -A[1][0] * od; 324e3cd5995SDave May inv[1][1] = A[0][0] * od; 325e3cd5995SDave May 326e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 327e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 329e3cd5995SDave May } 330e3cd5995SDave May */ 331e3cd5995SDave May 332d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[], PetscScalar coords[], PetscReal xip[], PetscReal *dJ) 333d71ae5a4SJacob Faibussowitsch { 334e3cd5995SDave May PetscReal x1, y1, x2, y2, x3, y3; 335e3cd5995SDave May PetscReal b[2], A[2][2], inv[2][2], detJ, od; 336e3cd5995SDave May 337e3cd5995SDave May PetscFunctionBegin; 338a5f152d1SDave May x1 = PetscRealPart(coords[2 * 0 + 0]); 339a5f152d1SDave May x2 = PetscRealPart(coords[2 * 1 + 0]); 340a5f152d1SDave May x3 = PetscRealPart(coords[2 * 2 + 0]); 341e3cd5995SDave May 342a5f152d1SDave May y1 = PetscRealPart(coords[2 * 0 + 1]); 343a5f152d1SDave May y2 = PetscRealPart(coords[2 * 1 + 1]); 344a5f152d1SDave May y3 = PetscRealPart(coords[2 * 2 + 1]); 345e3cd5995SDave May 346e3cd5995SDave May b[0] = xp[0] - x1; 347e3cd5995SDave May b[1] = xp[1] - y1; 348e3cd5995SDave May 3499371c9d4SSatish Balay A[0][0] = x2 - x1; 3509371c9d4SSatish Balay A[0][1] = x3 - x1; 3519371c9d4SSatish Balay A[1][0] = y2 - y1; 3529371c9d4SSatish Balay A[1][1] = y3 - y1; 353e3cd5995SDave May 354e3cd5995SDave May detJ = A[0][0] * A[1][1] - A[0][1] * A[1][0]; 355e3cd5995SDave May *dJ = PetscAbsReal(detJ); 356e3cd5995SDave May od = 1.0 / detJ; 357e3cd5995SDave May 358e3cd5995SDave May inv[0][0] = A[1][1] * od; 359e3cd5995SDave May inv[0][1] = -A[0][1] * od; 360e3cd5995SDave May inv[1][0] = -A[1][0] * od; 361e3cd5995SDave May inv[1][1] = A[0][0] * od; 362e3cd5995SDave May 363e3cd5995SDave May xip[0] = inv[0][0] * b[0] + inv[0][1] * b[1]; 364e3cd5995SDave May xip[1] = inv[1][0] * b[0] + inv[1][1] * b[1]; 3653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 366e3cd5995SDave May } 367e3cd5995SDave May 368ba38deedSJacob Faibussowitsch static PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm, PetscReal *swarm_field, DM dm, Vec v_field) 369d71ae5a4SJacob Faibussowitsch { 370e3cd5995SDave May const PetscReal PLEX_C_EPS = 1.0e-8; 371e3cd5995SDave May Vec v_field_l, denom_l, coor_l, denom; 372e3cd5995SDave May PetscInt k, p, e, npoints; 373e3cd5995SDave May PetscInt *mpfield_cell; 374e3cd5995SDave May PetscReal *mpfield_coor; 375a5f152d1SDave May PetscReal xi_p[2]; 376a5f152d1SDave May PetscScalar Ni[3]; 377e3cd5995SDave May PetscSection coordSection; 378e3cd5995SDave May PetscScalar *elcoor = NULL; 379e3cd5995SDave May 380e3cd5995SDave May PetscFunctionBegin; 3819566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field)); 382e3cd5995SDave May 3839566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v_field_l)); 3849566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &denom)); 3859566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &denom_l)); 3869566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field_l)); 3879566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom)); 3889566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom_l)); 389e3cd5995SDave May 3909566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coor_l)); 3919566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 392e3cd5995SDave May 3939566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(swarm, &npoints)); 3949566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 3959566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 396e3cd5995SDave May 397e3cd5995SDave May for (p = 0; p < npoints; p++) { 398a5f152d1SDave May PetscReal *coor_p, dJ; 399a5f152d1SDave May PetscScalar elfield[3]; 400e3cd5995SDave May PetscBool point_located; 401e3cd5995SDave May 402e3cd5995SDave May e = mpfield_cell[p]; 403e3cd5995SDave May coor_p = &mpfield_coor[2 * p]; 404e3cd5995SDave May 4059566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 406e3cd5995SDave May 407e3cd5995SDave May /* 408e3cd5995SDave May while (!point_located && (failed_counter < 25)) { 4099566063dSJacob Faibussowitsch PetscCall(PointInTriangle(point, coords[0], coords[1], coords[2], &point_located)); 410e3cd5995SDave May point.x = coor_p[0]; 411e3cd5995SDave May point.y = coor_p[1]; 412e3cd5995SDave May point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 413e3cd5995SDave May point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 414e3cd5995SDave May failed_counter++; 415e3cd5995SDave May } 416e3cd5995SDave May 417e3cd5995SDave May if (!point_located) { 41863a3b9bcSJacob Faibussowitsch PetscPrintf(PETSC_COMM_SELF,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ") with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e) in %" PetscInt_FMT " iterations\n",point.x,point.y,e,coords[0].x,coords[0].y,coords[1].x,coords[1].y,coords[2].x,coords[2].y,failed_counter); 419e3cd5995SDave May } 420e3cd5995SDave May 42163a3b9bcSJacob Faibussowitsch PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ")",point.x,point.y,e); 422e3cd5995SDave May else { 4239566063dSJacob Faibussowitsch PetscCall(_ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ)); 424e3cd5995SDave May xi_p[0] = 0.5*(xi_p[0] + 1.0); 425e3cd5995SDave May xi_p[1] = 0.5*(xi_p[1] + 1.0); 426e3cd5995SDave May 42763a3b9bcSJacob Faibussowitsch PetscPrintf(PETSC_COMM_SELF,"[p=%" PetscInt_FMT "] x(%+1.4e,%+1.4e) -> mapped to element %" PetscInt_FMT " xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]); 428e3cd5995SDave May 429e3cd5995SDave May } 430e3cd5995SDave May */ 431e3cd5995SDave May 4329566063dSJacob Faibussowitsch PetscCall(ComputeLocalCoordinateAffine2d(coor_p, elcoor, xi_p, &dJ)); 433e3cd5995SDave May /* 43463a3b9bcSJacob Faibussowitsch PetscPrintf(PETSC_COMM_SELF,"[p=%" PetscInt_FMT "] x(%+1.4e,%+1.4e) -> mapped to element %" PetscInt_FMT " xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]); 435e3cd5995SDave May */ 436e3cd5995SDave May /* 437e3cd5995SDave May point_located = PETSC_TRUE; 438e3cd5995SDave May if (xi_p[0] < 0.0) { 439e3cd5995SDave May if (xi_p[0] > -PLEX_C_EPS) { 440e3cd5995SDave May xi_p[0] = 0.0; 441e3cd5995SDave May } else { 442e3cd5995SDave May point_located = PETSC_FALSE; 443e3cd5995SDave May } 444e3cd5995SDave May } 445e3cd5995SDave May if (xi_p[1] < 0.0) { 446e3cd5995SDave May if (xi_p[1] > -PLEX_C_EPS) { 447e3cd5995SDave May xi_p[1] = 0.0; 448e3cd5995SDave May } else { 449e3cd5995SDave May point_located = PETSC_FALSE; 450e3cd5995SDave May } 451e3cd5995SDave May } 452e3cd5995SDave May if (xi_p[1] > (1.0-xi_p[0])) { 453e3cd5995SDave May if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { 454e3cd5995SDave May xi_p[1] = 1.0 - xi_p[0]; 455e3cd5995SDave May } else { 456e3cd5995SDave May point_located = PETSC_FALSE; 457e3cd5995SDave May } 458e3cd5995SDave May } 459e3cd5995SDave May if (!point_located) { 460e3cd5995SDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); 46163a3b9bcSJacob Faibussowitsch PetscPrintf(PETSC_COMM_SELF,"[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ") with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e)\n",coor_p[0],coor_p[1],e,elcoor[0],elcoor[1],elcoor[2],elcoor[3],elcoor[4],elcoor[5]); 462e3cd5995SDave May } 46363a3b9bcSJacob Faibussowitsch PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ")",coor_p[0],coor_p[1],e); 464e3cd5995SDave May */ 465e3cd5995SDave May 466e3cd5995SDave May Ni[0] = 1.0 - xi_p[0] - xi_p[1]; 467e3cd5995SDave May Ni[1] = xi_p[0]; 468e3cd5995SDave May Ni[2] = xi_p[1]; 469e3cd5995SDave May 470e3cd5995SDave May point_located = PETSC_TRUE; 471e3cd5995SDave May for (k = 0; k < 3; k++) { 472a5f152d1SDave May if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; 473a5f152d1SDave May if (PetscRealPart(Ni[k]) > (1.0 + PLEX_C_EPS)) point_located = PETSC_FALSE; 474e3cd5995SDave May } 475e3cd5995SDave May if (!point_located) { 4769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[Error] xi,eta = %+1.8e, %+1.8e\n", (double)xi_p[0], (double)xi_p[1])); 47763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ") with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e)\n", (double)coor_p[0], (double)coor_p[1], e, (double)PetscRealPart(elcoor[0]), (double)PetscRealPart(elcoor[1]), (double)PetscRealPart(elcoor[2]), (double)PetscRealPart(elcoor[3]), (double)PetscRealPart(elcoor[4]), (double)PetscRealPart(elcoor[5]))); 478e3cd5995SDave May } 47963a3b9bcSJacob Faibussowitsch PetscCheck(point_located, PETSC_COMM_SELF, PETSC_ERR_SUP, "Failed to locate point (%1.8e,%1.8e) in local mesh (cell %" PetscInt_FMT ")", (double)coor_p[0], (double)coor_p[1], e); 480e3cd5995SDave May 481e3cd5995SDave May for (k = 0; k < 3; k++) { 482e3cd5995SDave May Ni[k] = Ni[k] * dJ; 483e3cd5995SDave May elfield[k] = Ni[k] * swarm_field[p]; 484e3cd5995SDave May } 4859566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 486e3cd5995SDave May 4879566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, v_field_l, e, elfield, ADD_VALUES)); 4889566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, denom_l, e, Ni, ADD_VALUES)); 489e3cd5995SDave May } 490e3cd5995SDave May 4919566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 4929566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 493e3cd5995SDave May 4949566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, v_field_l, ADD_VALUES, v_field)); 4959566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, v_field_l, ADD_VALUES, v_field)); 4969566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, denom_l, ADD_VALUES, denom)); 4979566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, denom_l, ADD_VALUES, denom)); 498e3cd5995SDave May 4999566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(v_field, v_field, denom)); 500e3cd5995SDave May 5019566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v_field_l)); 5029566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &denom_l)); 5039566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &denom)); 5043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 505e3cd5995SDave May } 506e3cd5995SDave May 507d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]) 508d71ae5a4SJacob Faibussowitsch { 509e3cd5995SDave May PetscInt f, dim; 510e3cd5995SDave May 511e3cd5995SDave May PetscFunctionBegin; 5129566063dSJacob Faibussowitsch PetscCall(DMGetDimension(swarm, &dim)); 513e3cd5995SDave May switch (dim) { 514e3cd5995SDave May case 2: 515e3cd5995SDave May for (f = 0; f < nfields; f++) { 516e3cd5995SDave May PetscReal *swarm_field; 517e3cd5995SDave May 5189566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(dfield[f], (void **)&swarm_field)); 5199566063dSJacob Faibussowitsch PetscCall(DMSwarmProjectField_ApproxP1_PLEX_2D(swarm, swarm_field, celldm, vecs[f])); 520e3cd5995SDave May } 521e3cd5995SDave May break; 522d71ae5a4SJacob Faibussowitsch case 3: 523d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "No support for 3D"); 524d71ae5a4SJacob Faibussowitsch default: 525d71ae5a4SJacob Faibussowitsch break; 526e3cd5995SDave May } 5273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 528e3cd5995SDave May } 529431382f2SDave May 530d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm, DM dmc, PetscInt npoints, PetscReal xi[]) 531d71ae5a4SJacob Faibussowitsch { 532431382f2SDave May PetscBool is_simplex, is_tensorcell; 53392e40656SDave May PetscInt dim, nfaces, ps, pe, p, d, nbasis, pcnt, e, k, nel; 53492e40656SDave May PetscFE fe; 53592e40656SDave May PetscQuadrature quadrature; 536ef0bb6c7SMatthew G. Knepley PetscTabulation T; 537ef0bb6c7SMatthew G. Knepley PetscReal *xiq; 53892e40656SDave May Vec coorlocal; 53992e40656SDave May PetscSection coordSection; 54092e40656SDave May PetscScalar *elcoor = NULL; 54192e40656SDave May PetscReal *swarm_coor; 54292e40656SDave May PetscInt *swarm_cellid; 543431382f2SDave May 54492e40656SDave May PetscFunctionBegin; 5459566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim)); 546431382f2SDave May 547431382f2SDave May is_simplex = PETSC_FALSE; 548431382f2SDave May is_tensorcell = PETSC_FALSE; 5499566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 5509566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 551431382f2SDave May 552ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 553431382f2SDave May 554431382f2SDave May switch (dim) { 555431382f2SDave May case 2: 556ad540459SPierre Jolivet if (nfaces == 4) is_tensorcell = PETSC_TRUE; 557431382f2SDave May break; 558431382f2SDave May case 3: 559ad540459SPierre Jolivet if (nfaces == 6) is_tensorcell = PETSC_TRUE; 560431382f2SDave May break; 561d71ae5a4SJacob Faibussowitsch default: 562d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for 2D, 3D"); 563431382f2SDave May } 564431382f2SDave May 56592e40656SDave May /* check points provided fail inside the reference cell */ 566431382f2SDave May if (is_simplex) { 56792e40656SDave May for (p = 0; p < npoints; p++) { 56892e40656SDave May PetscReal sum; 569ad540459SPierre 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"); 57092e40656SDave May sum = 0.0; 571ad540459SPierre Jolivet for (d = 0; d < dim; d++) sum += xi[dim * p + d]; 57208401ef6SPierre Jolivet PetscCheck(sum <= 0.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain"); 57392e40656SDave May } 574431382f2SDave May } else if (is_tensorcell) { 57592e40656SDave May for (p = 0; p < npoints; p++) { 576ad540459SPierre 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"); 57792e40656SDave May } 578431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for d-simplex and d-tensorcell"); 579431382f2SDave May 5809566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm), &quadrature)); 5819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints * dim, &xiq)); 5829566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(xiq, xi, npoints * dim)); 5839566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(quadrature, dim, 1, npoints, (const PetscReal *)xiq, NULL)); 5849566063dSJacob Faibussowitsch PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); 5859566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(fe, quadrature)); 5869566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fe, &nbasis)); 5879566063dSJacob Faibussowitsch PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); 58892e40656SDave May 58992e40656SDave May /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ 59092e40656SDave May /* 0->cell, 1->edge, 2->vert */ 5919566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 59292e40656SDave May nel = pe - ps; 59392e40656SDave May 5949566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints * nel, -1)); 5959566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 5969566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 59792e40656SDave May 5989566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 5999566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 60092e40656SDave May 60192e40656SDave May pcnt = 0; 60292e40656SDave May for (e = 0; e < nel; e++) { 6039566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 60492e40656SDave May 60592e40656SDave May for (p = 0; p < npoints; p++) { 60692e40656SDave May for (d = 0; d < dim; d++) { 60792e40656SDave May swarm_coor[dim * pcnt + d] = 0.0; 608ad540459SPierre Jolivet for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][p * nbasis + k] * PetscRealPart(elcoor[dim * k + d]); 60992e40656SDave May } 61092e40656SDave May swarm_cellid[pcnt] = e; 61192e40656SDave May pcnt++; 61292e40656SDave May } 6139566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 61492e40656SDave May } 6159566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 6169566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 61792e40656SDave May 6189566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quadrature)); 6199566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 6203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 621431382f2SDave May } 622