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 6b12d2206SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *xi); 7b12d2206SDave May 8d71ae5a4SJacob Faibussowitsch static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem) 9d71ae5a4SJacob Faibussowitsch { 102cf64d79SDave May const PetscInt Nc = 1; 112cf64d79SDave May PetscQuadrature q, fq; 122cf64d79SDave May DM K; 132cf64d79SDave May PetscSpace P; 142cf64d79SDave May PetscDualSpace Q; 152cf64d79SDave May PetscInt order, quadPointsPerEdge; 162cf64d79SDave May PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 172cf64d79SDave May 182cf64d79SDave May PetscFunctionBegin; 192cf64d79SDave May /* Create space */ 209566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)dm), &P)); 219566063dSJacob Faibussowitsch /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix)); */ 229566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialSetTensor(P, tensor)); 239566063dSJacob Faibussowitsch /* PetscCall(PetscSpaceSetFromOptions(P)); */ 249566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 259566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(P, 1, PETSC_DETERMINE)); 269566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(P, Nc)); 279566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(P, dim)); 289566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(P)); 299566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(P, &order, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialGetTensor(P, &tensor)); 312cf64d79SDave May /* Create dual space */ 329566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)dm), &Q)); 339566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 349566063dSJacob Faibussowitsch /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix)); */ 359566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K)); 369566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(Q, K)); 379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K)); 389566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 399566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(Q, order)); 409566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTensor(Q, tensor)); 419566063dSJacob Faibussowitsch /* PetscCall(PetscDualSpaceSetFromOptions(Q)); */ 429566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 439566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(Q)); 442cf64d79SDave May /* Create element */ 459566063dSJacob Faibussowitsch PetscCall(PetscFECreate(PetscObjectComm((PetscObject)dm), fem)); 469566063dSJacob Faibussowitsch /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix)); */ 479566063dSJacob Faibussowitsch /* PetscCall(PetscFESetFromOptions(*fem)); */ 489566063dSJacob Faibussowitsch PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 499566063dSJacob Faibussowitsch PetscCall(PetscFESetBasisSpace(*fem, P)); 509566063dSJacob Faibussowitsch PetscCall(PetscFESetDualSpace(*fem, Q)); 519566063dSJacob Faibussowitsch PetscCall(PetscFESetNumComponents(*fem, Nc)); 529566063dSJacob Faibussowitsch PetscCall(PetscFESetUp(*fem)); 539566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&P)); 549566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&Q)); 552cf64d79SDave May /* Create quadrature (with specified order if given) */ 562cf64d79SDave May qorder = qorder >= 0 ? qorder : order; 572cf64d79SDave May quadPointsPerEdge = PetscMax(qorder + 1, 1); 582cf64d79SDave May if (isSimplex) { 599566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); 609566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); 619371c9d4SSatish Balay } else { 629566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); 639566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); 642cf64d79SDave May } 659566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(*fem, q)); 669566063dSJacob Faibussowitsch PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 679566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 689566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&fq)); 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 702cf64d79SDave May } 712cf64d79SDave May 72*ba38deedSJacob Faibussowitsch static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm, DM dmc, PetscInt nsub) 73d71ae5a4SJacob Faibussowitsch { 7478185223SDave May PetscInt dim, nfaces, nbasis; 7578185223SDave May PetscInt q, npoints_q, e, nel, pcnt, ps, pe, d, k, r; 76ef0bb6c7SMatthew G. Knepley PetscTabulation T; 7778185223SDave May Vec coorlocal; 7878185223SDave May PetscSection coordSection; 7978185223SDave May PetscScalar *elcoor = NULL; 8078185223SDave May PetscReal *swarm_coor; 8178185223SDave May PetscInt *swarm_cellid; 8278185223SDave May const PetscReal *xiq; 8378185223SDave May PetscQuadrature quadrature; 8478185223SDave May PetscFE fe, feRef; 8578185223SDave May PetscBool is_simplex; 8678185223SDave May 8778185223SDave May PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim)); 8978185223SDave May is_simplex = PETSC_FALSE; 909566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 919566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 92ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 9378185223SDave May 949566063dSJacob Faibussowitsch PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); 9578185223SDave May 9678185223SDave May for (r = 0; r < nsub; r++) { 979566063dSJacob Faibussowitsch PetscCall(PetscFERefine(fe, &feRef)); 989566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(feRef, fe)); 999566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feRef)); 10078185223SDave May } 10178185223SDave May 1029566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 1039566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fe, &nbasis)); 1059566063dSJacob Faibussowitsch PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); 10678185223SDave May 10778185223SDave May /* 0->cell, 1->edge, 2->vert */ 1089566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 10978185223SDave May nel = pe - ps; 11078185223SDave May 1119566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 1129566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1139566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 11478185223SDave May 1159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 1169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 11778185223SDave May 11878185223SDave May pcnt = 0; 11978185223SDave May for (e = 0; e < nel; e++) { 1209566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 12178185223SDave May 12278185223SDave May for (q = 0; q < npoints_q; q++) { 12378185223SDave May for (d = 0; d < dim; d++) { 12478185223SDave May swarm_coor[dim * pcnt + d] = 0.0; 125ad540459SPierre Jolivet for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][q * nbasis + k] * PetscRealPart(elcoor[dim * k + d]); 12678185223SDave May } 12778185223SDave May swarm_cellid[pcnt] = e; 12878185223SDave May pcnt++; 12978185223SDave May } 1309566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 13178185223SDave May } 1329566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1339566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 13478185223SDave May 1359566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13778185223SDave May } 13878185223SDave May 139*ba38deedSJacob Faibussowitsch static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm, DM dmc, PetscInt npoints) 140d71ae5a4SJacob Faibussowitsch { 141bc53056eSDave May PetscInt dim; 142bc53056eSDave May PetscInt ii, jj, q, npoints_q, e, nel, npe, pcnt, ps, pe, d, k, nfaces; 1430e2ec84fSDave May PetscReal *xi, ds, ds2; 1440e2ec84fSDave May PetscReal **basis; 1450e2ec84fSDave May Vec coorlocal; 1460e2ec84fSDave May PetscSection coordSection; 1470e2ec84fSDave May PetscScalar *elcoor = NULL; 1480e2ec84fSDave May PetscReal *swarm_coor; 1490e2ec84fSDave May PetscInt *swarm_cellid; 150bc53056eSDave May PetscBool is_simplex; 1510e2ec84fSDave May 1520e2ec84fSDave May PetscFunctionBegin; 1539566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim)); 15408401ef6SPierre Jolivet PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only 2D is supported"); 155bc53056eSDave May is_simplex = PETSC_FALSE; 1569566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 1579566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 158ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 15928b400f6SJacob Faibussowitsch PetscCheck(is_simplex, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only the simplex is supported"); 160bc53056eSDave May 1619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints * npoints, &xi)); 1620e2ec84fSDave May pcnt = 0; 1630e2ec84fSDave May ds = 1.0 / ((PetscReal)(npoints - 1)); 1640e2ec84fSDave May ds2 = 1.0 / ((PetscReal)(npoints)); 1650e2ec84fSDave May for (jj = 0; jj < npoints; jj++) { 1660e2ec84fSDave May for (ii = 0; ii < npoints - jj; ii++) { 1670e2ec84fSDave May xi[dim * pcnt + 0] = ii * ds; 1680e2ec84fSDave May xi[dim * pcnt + 1] = jj * ds; 1690e2ec84fSDave May 1700e2ec84fSDave May xi[dim * pcnt + 0] *= (1.0 - 1.2 * ds2); 1710e2ec84fSDave May xi[dim * pcnt + 1] *= (1.0 - 1.2 * ds2); 1720e2ec84fSDave May 1730e2ec84fSDave May xi[dim * pcnt + 0] += 0.35 * ds2; 1740e2ec84fSDave May xi[dim * pcnt + 1] += 0.35 * ds2; 1750e2ec84fSDave May pcnt++; 1760e2ec84fSDave May } 1770e2ec84fSDave May } 1780e2ec84fSDave May npoints_q = pcnt; 1790e2ec84fSDave May 1800e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 1819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints_q, &basis)); 1820e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 1839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npe, &basis[q])); 1840e2ec84fSDave May 1850e2ec84fSDave May basis[q][0] = 1.0 - xi[dim * q + 0] - xi[dim * q + 1]; 1860e2ec84fSDave May basis[q][1] = xi[dim * q + 0]; 1870e2ec84fSDave May basis[q][2] = xi[dim * q + 1]; 1880e2ec84fSDave May } 1890e2ec84fSDave May 1900e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 1919566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 1920e2ec84fSDave May nel = pe - ps; 1930e2ec84fSDave May 1949566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 1959566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1969566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1970e2ec84fSDave May 1989566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 1999566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 2000e2ec84fSDave May 2010e2ec84fSDave May pcnt = 0; 2020e2ec84fSDave May for (e = 0; e < nel; e++) { 2039566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 2040e2ec84fSDave May 2050e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 2060e2ec84fSDave May for (d = 0; d < dim; d++) { 2070e2ec84fSDave May swarm_coor[dim * pcnt + d] = 0.0; 208ad540459SPierre Jolivet for (k = 0; k < npe; k++) swarm_coor[dim * pcnt + d] += basis[q][k] * PetscRealPart(elcoor[dim * k + d]); 2090e2ec84fSDave May } 2100e2ec84fSDave May swarm_cellid[pcnt] = e; 2110e2ec84fSDave May pcnt++; 2120e2ec84fSDave May } 2139566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 2140e2ec84fSDave May } 2159566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 2169566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 2170e2ec84fSDave May 2189566063dSJacob Faibussowitsch PetscCall(PetscFree(xi)); 21948a46eb9SPierre Jolivet for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q])); 2209566063dSJacob Faibussowitsch PetscCall(PetscFree(basis)); 2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2220e2ec84fSDave May } 2230e2ec84fSDave May 224d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param) 225d71ae5a4SJacob Faibussowitsch { 2260e2ec84fSDave May PetscInt dim; 2270e2ec84fSDave May 2280e2ec84fSDave May PetscFunctionBegin; 2299566063dSJacob Faibussowitsch PetscCall(DMGetDimension(celldm, &dim)); 2300e2ec84fSDave May switch (layout) { 2310e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 23208401ef6SPierre Jolivet PetscCheck(dim != 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No 3D support for REGULAR+PLEX"); 2339566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm, celldm, layout_param)); 2340e2ec84fSDave May break; 2359371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_GAUSS: { 236b12d2206SDave May PetscInt npoints, npoints1, ps, pe, nfaces; 237b12d2206SDave May const PetscReal *xi; 238b12d2206SDave May PetscBool is_simplex; 239b12d2206SDave May PetscQuadrature quadrature; 240b12d2206SDave May 241b12d2206SDave May is_simplex = PETSC_FALSE; 2429566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe)); 2439566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(celldm, ps, &nfaces)); 244ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 245b12d2206SDave May 246b12d2206SDave May npoints1 = layout_param; 247b12d2206SDave May if (is_simplex) { 2489566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, npoints1, -1.0, 1.0, &quadrature)); 249b12d2206SDave May } else { 2509566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, npoints1, -1.0, 1.0, &quadrature)); 251b12d2206SDave May } 2529566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints, &xi, NULL)); 2539566063dSJacob Faibussowitsch PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, (PetscReal *)xi)); 2549566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quadrature)); 2559371c9d4SSatish Balay } break; 256d71ae5a4SJacob Faibussowitsch case DMSWARMPIC_LAYOUT_SUBDIVISION: 257d71ae5a4SJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm, celldm, layout_param)); 258d71ae5a4SJacob Faibussowitsch break; 2590e2ec84fSDave May } 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2610e2ec84fSDave May } 262e3cd5995SDave May 263e3cd5995SDave May /* 264e3cd5995SDave May typedef struct { 265e3cd5995SDave May PetscReal x,y; 266e3cd5995SDave May } Point2d; 267e3cd5995SDave May 268e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s) 269e3cd5995SDave May { 270e3cd5995SDave May PetscFunctionBegin; 271e3cd5995SDave May *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y); 2723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 273e3cd5995SDave May } 274e3cd5995SDave May */ 275e3cd5995SDave May /* 276e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v) 277e3cd5995SDave May { 278e3cd5995SDave May PetscReal s1,s2,s3; 279e3cd5995SDave May PetscBool b1, b2, b3; 280e3cd5995SDave May 281e3cd5995SDave May PetscFunctionBegin; 282e3cd5995SDave May signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f; 283e3cd5995SDave May signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f; 284e3cd5995SDave May signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f; 285e3cd5995SDave May 286e3cd5995SDave May *v = ((b1 == b2) && (b2 == b3)); 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 288e3cd5995SDave May } 289e3cd5995SDave May */ 290e3cd5995SDave May /* 291e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ) 292e3cd5995SDave May { 293e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 294e3cd5995SDave May PetscReal c,b[2],A[2][2],inv[2][2],detJ,od; 295e3cd5995SDave May 296e3cd5995SDave May PetscFunctionBegin; 297e3cd5995SDave May x1 = coords[2*0+0]; 298e3cd5995SDave May x2 = coords[2*1+0]; 299e3cd5995SDave May x3 = coords[2*2+0]; 300e3cd5995SDave May 301e3cd5995SDave May y1 = coords[2*0+1]; 302e3cd5995SDave May y2 = coords[2*1+1]; 303e3cd5995SDave May y3 = coords[2*2+1]; 304e3cd5995SDave May 305e3cd5995SDave May c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3; 306e3cd5995SDave May b[0] = xp[0] - c; 307e3cd5995SDave May c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3; 308e3cd5995SDave May b[1] = xp[1] - c; 309e3cd5995SDave May 310e3cd5995SDave May A[0][0] = -0.5*x1 + 0.5*x2; A[0][1] = -0.5*x1 + 0.5*x3; 311e3cd5995SDave May A[1][0] = -0.5*y1 + 0.5*y2; A[1][1] = -0.5*y1 + 0.5*y3; 312e3cd5995SDave May 313e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 314e3cd5995SDave May *dJ = PetscAbsReal(detJ); 315e3cd5995SDave May od = 1.0/detJ; 316e3cd5995SDave May 317e3cd5995SDave May inv[0][0] = A[1][1] * od; 318e3cd5995SDave May inv[0][1] = -A[0][1] * od; 319e3cd5995SDave May inv[1][0] = -A[1][0] * od; 320e3cd5995SDave May inv[1][1] = A[0][0] * od; 321e3cd5995SDave May 322e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 323e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 3243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 325e3cd5995SDave May } 326e3cd5995SDave May */ 327e3cd5995SDave May 328d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[], PetscScalar coords[], PetscReal xip[], PetscReal *dJ) 329d71ae5a4SJacob Faibussowitsch { 330e3cd5995SDave May PetscReal x1, y1, x2, y2, x3, y3; 331e3cd5995SDave May PetscReal b[2], A[2][2], inv[2][2], detJ, od; 332e3cd5995SDave May 333e3cd5995SDave May PetscFunctionBegin; 334a5f152d1SDave May x1 = PetscRealPart(coords[2 * 0 + 0]); 335a5f152d1SDave May x2 = PetscRealPart(coords[2 * 1 + 0]); 336a5f152d1SDave May x3 = PetscRealPart(coords[2 * 2 + 0]); 337e3cd5995SDave May 338a5f152d1SDave May y1 = PetscRealPart(coords[2 * 0 + 1]); 339a5f152d1SDave May y2 = PetscRealPart(coords[2 * 1 + 1]); 340a5f152d1SDave May y3 = PetscRealPart(coords[2 * 2 + 1]); 341e3cd5995SDave May 342e3cd5995SDave May b[0] = xp[0] - x1; 343e3cd5995SDave May b[1] = xp[1] - y1; 344e3cd5995SDave May 3459371c9d4SSatish Balay A[0][0] = x2 - x1; 3469371c9d4SSatish Balay A[0][1] = x3 - x1; 3479371c9d4SSatish Balay A[1][0] = y2 - y1; 3489371c9d4SSatish Balay A[1][1] = y3 - y1; 349e3cd5995SDave May 350e3cd5995SDave May detJ = A[0][0] * A[1][1] - A[0][1] * A[1][0]; 351e3cd5995SDave May *dJ = PetscAbsReal(detJ); 352e3cd5995SDave May od = 1.0 / detJ; 353e3cd5995SDave May 354e3cd5995SDave May inv[0][0] = A[1][1] * od; 355e3cd5995SDave May inv[0][1] = -A[0][1] * od; 356e3cd5995SDave May inv[1][0] = -A[1][0] * od; 357e3cd5995SDave May inv[1][1] = A[0][0] * od; 358e3cd5995SDave May 359e3cd5995SDave May xip[0] = inv[0][0] * b[0] + inv[0][1] * b[1]; 360e3cd5995SDave May xip[1] = inv[1][0] * b[0] + inv[1][1] * b[1]; 3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 362e3cd5995SDave May } 363e3cd5995SDave May 364*ba38deedSJacob Faibussowitsch static PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm, PetscReal *swarm_field, DM dm, Vec v_field) 365d71ae5a4SJacob Faibussowitsch { 366e3cd5995SDave May const PetscReal PLEX_C_EPS = 1.0e-8; 367e3cd5995SDave May Vec v_field_l, denom_l, coor_l, denom; 368e3cd5995SDave May PetscInt k, p, e, npoints; 369e3cd5995SDave May PetscInt *mpfield_cell; 370e3cd5995SDave May PetscReal *mpfield_coor; 371a5f152d1SDave May PetscReal xi_p[2]; 372a5f152d1SDave May PetscScalar Ni[3]; 373e3cd5995SDave May PetscSection coordSection; 374e3cd5995SDave May PetscScalar *elcoor = NULL; 375e3cd5995SDave May 376e3cd5995SDave May PetscFunctionBegin; 3779566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field)); 378e3cd5995SDave May 3799566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v_field_l)); 3809566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &denom)); 3819566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &denom_l)); 3829566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field_l)); 3839566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom)); 3849566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom_l)); 385e3cd5995SDave May 3869566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coor_l)); 3879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 388e3cd5995SDave May 3899566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(swarm, &npoints)); 3909566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 3919566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 392e3cd5995SDave May 393e3cd5995SDave May for (p = 0; p < npoints; p++) { 394a5f152d1SDave May PetscReal *coor_p, dJ; 395a5f152d1SDave May PetscScalar elfield[3]; 396e3cd5995SDave May PetscBool point_located; 397e3cd5995SDave May 398e3cd5995SDave May e = mpfield_cell[p]; 399e3cd5995SDave May coor_p = &mpfield_coor[2 * p]; 400e3cd5995SDave May 4019566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 402e3cd5995SDave May 403e3cd5995SDave May /* 404e3cd5995SDave May while (!point_located && (failed_counter < 25)) { 4059566063dSJacob Faibussowitsch PetscCall(PointInTriangle(point, coords[0], coords[1], coords[2], &point_located)); 406e3cd5995SDave May point.x = coor_p[0]; 407e3cd5995SDave May point.y = coor_p[1]; 408e3cd5995SDave May point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 409e3cd5995SDave May point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 410e3cd5995SDave May failed_counter++; 411e3cd5995SDave May } 412e3cd5995SDave May 413e3cd5995SDave May if (!point_located) { 41463a3b9bcSJacob 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); 415e3cd5995SDave May } 416e3cd5995SDave May 41763a3b9bcSJacob 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); 418e3cd5995SDave May else { 4199566063dSJacob Faibussowitsch PetscCall(_ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ)); 420e3cd5995SDave May xi_p[0] = 0.5*(xi_p[0] + 1.0); 421e3cd5995SDave May xi_p[1] = 0.5*(xi_p[1] + 1.0); 422e3cd5995SDave May 42363a3b9bcSJacob 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]); 424e3cd5995SDave May 425e3cd5995SDave May } 426e3cd5995SDave May */ 427e3cd5995SDave May 4289566063dSJacob Faibussowitsch PetscCall(ComputeLocalCoordinateAffine2d(coor_p, elcoor, xi_p, &dJ)); 429e3cd5995SDave May /* 43063a3b9bcSJacob 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]); 431e3cd5995SDave May */ 432e3cd5995SDave May /* 433e3cd5995SDave May point_located = PETSC_TRUE; 434e3cd5995SDave May if (xi_p[0] < 0.0) { 435e3cd5995SDave May if (xi_p[0] > -PLEX_C_EPS) { 436e3cd5995SDave May xi_p[0] = 0.0; 437e3cd5995SDave May } else { 438e3cd5995SDave May point_located = PETSC_FALSE; 439e3cd5995SDave May } 440e3cd5995SDave May } 441e3cd5995SDave May if (xi_p[1] < 0.0) { 442e3cd5995SDave May if (xi_p[1] > -PLEX_C_EPS) { 443e3cd5995SDave May xi_p[1] = 0.0; 444e3cd5995SDave May } else { 445e3cd5995SDave May point_located = PETSC_FALSE; 446e3cd5995SDave May } 447e3cd5995SDave May } 448e3cd5995SDave May if (xi_p[1] > (1.0-xi_p[0])) { 449e3cd5995SDave May if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { 450e3cd5995SDave May xi_p[1] = 1.0 - xi_p[0]; 451e3cd5995SDave May } else { 452e3cd5995SDave May point_located = PETSC_FALSE; 453e3cd5995SDave May } 454e3cd5995SDave May } 455e3cd5995SDave May if (!point_located) { 456e3cd5995SDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); 45763a3b9bcSJacob 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]); 458e3cd5995SDave May } 45963a3b9bcSJacob 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); 460e3cd5995SDave May */ 461e3cd5995SDave May 462e3cd5995SDave May Ni[0] = 1.0 - xi_p[0] - xi_p[1]; 463e3cd5995SDave May Ni[1] = xi_p[0]; 464e3cd5995SDave May Ni[2] = xi_p[1]; 465e3cd5995SDave May 466e3cd5995SDave May point_located = PETSC_TRUE; 467e3cd5995SDave May for (k = 0; k < 3; k++) { 468a5f152d1SDave May if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; 469a5f152d1SDave May if (PetscRealPart(Ni[k]) > (1.0 + PLEX_C_EPS)) point_located = PETSC_FALSE; 470e3cd5995SDave May } 471e3cd5995SDave May if (!point_located) { 4729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[Error] xi,eta = %+1.8e, %+1.8e\n", (double)xi_p[0], (double)xi_p[1])); 47363a3b9bcSJacob 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]))); 474e3cd5995SDave May } 47563a3b9bcSJacob 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); 476e3cd5995SDave May 477e3cd5995SDave May for (k = 0; k < 3; k++) { 478e3cd5995SDave May Ni[k] = Ni[k] * dJ; 479e3cd5995SDave May elfield[k] = Ni[k] * swarm_field[p]; 480e3cd5995SDave May } 4819566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 482e3cd5995SDave May 4839566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, v_field_l, e, elfield, ADD_VALUES)); 4849566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, denom_l, e, Ni, ADD_VALUES)); 485e3cd5995SDave May } 486e3cd5995SDave May 4879566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 4889566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 489e3cd5995SDave May 4909566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, v_field_l, ADD_VALUES, v_field)); 4919566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, v_field_l, ADD_VALUES, v_field)); 4929566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, denom_l, ADD_VALUES, denom)); 4939566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, denom_l, ADD_VALUES, denom)); 494e3cd5995SDave May 4959566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(v_field, v_field, denom)); 496e3cd5995SDave May 4979566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v_field_l)); 4989566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &denom_l)); 4999566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &denom)); 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 501e3cd5995SDave May } 502e3cd5995SDave May 503d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]) 504d71ae5a4SJacob Faibussowitsch { 505e3cd5995SDave May PetscInt f, dim; 506e3cd5995SDave May 507e3cd5995SDave May PetscFunctionBegin; 5089566063dSJacob Faibussowitsch PetscCall(DMGetDimension(swarm, &dim)); 509e3cd5995SDave May switch (dim) { 510e3cd5995SDave May case 2: 511e3cd5995SDave May for (f = 0; f < nfields; f++) { 512e3cd5995SDave May PetscReal *swarm_field; 513e3cd5995SDave May 5149566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(dfield[f], (void **)&swarm_field)); 5159566063dSJacob Faibussowitsch PetscCall(DMSwarmProjectField_ApproxP1_PLEX_2D(swarm, swarm_field, celldm, vecs[f])); 516e3cd5995SDave May } 517e3cd5995SDave May break; 518d71ae5a4SJacob Faibussowitsch case 3: 519d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "No support for 3D"); 520d71ae5a4SJacob Faibussowitsch default: 521d71ae5a4SJacob Faibussowitsch break; 522e3cd5995SDave May } 5233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 524e3cd5995SDave May } 525431382f2SDave May 526d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm, DM dmc, PetscInt npoints, PetscReal xi[]) 527d71ae5a4SJacob Faibussowitsch { 528431382f2SDave May PetscBool is_simplex, is_tensorcell; 52992e40656SDave May PetscInt dim, nfaces, ps, pe, p, d, nbasis, pcnt, e, k, nel; 53092e40656SDave May PetscFE fe; 53192e40656SDave May PetscQuadrature quadrature; 532ef0bb6c7SMatthew G. Knepley PetscTabulation T; 533ef0bb6c7SMatthew G. Knepley PetscReal *xiq; 53492e40656SDave May Vec coorlocal; 53592e40656SDave May PetscSection coordSection; 53692e40656SDave May PetscScalar *elcoor = NULL; 53792e40656SDave May PetscReal *swarm_coor; 53892e40656SDave May PetscInt *swarm_cellid; 539431382f2SDave May 54092e40656SDave May PetscFunctionBegin; 5419566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim)); 542431382f2SDave May 543431382f2SDave May is_simplex = PETSC_FALSE; 544431382f2SDave May is_tensorcell = PETSC_FALSE; 5459566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 5469566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 547431382f2SDave May 548ad540459SPierre Jolivet if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE; 549431382f2SDave May 550431382f2SDave May switch (dim) { 551431382f2SDave May case 2: 552ad540459SPierre Jolivet if (nfaces == 4) is_tensorcell = PETSC_TRUE; 553431382f2SDave May break; 554431382f2SDave May case 3: 555ad540459SPierre Jolivet if (nfaces == 6) is_tensorcell = PETSC_TRUE; 556431382f2SDave May break; 557d71ae5a4SJacob Faibussowitsch default: 558d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for 2D, 3D"); 559431382f2SDave May } 560431382f2SDave May 56192e40656SDave May /* check points provided fail inside the reference cell */ 562431382f2SDave May if (is_simplex) { 56392e40656SDave May for (p = 0; p < npoints; p++) { 56492e40656SDave May PetscReal sum; 565ad540459SPierre 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"); 56692e40656SDave May sum = 0.0; 567ad540459SPierre Jolivet for (d = 0; d < dim; d++) sum += xi[dim * p + d]; 56808401ef6SPierre Jolivet PetscCheck(sum <= 0.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain"); 56992e40656SDave May } 570431382f2SDave May } else if (is_tensorcell) { 57192e40656SDave May for (p = 0; p < npoints; p++) { 572ad540459SPierre 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"); 57392e40656SDave May } 574431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for d-simplex and d-tensorcell"); 575431382f2SDave May 5769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm), &quadrature)); 5779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints * dim, &xiq)); 5789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(xiq, xi, npoints * dim)); 5799566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(quadrature, dim, 1, npoints, (const PetscReal *)xiq, NULL)); 5809566063dSJacob Faibussowitsch PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); 5819566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(fe, quadrature)); 5829566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fe, &nbasis)); 5839566063dSJacob Faibussowitsch PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); 58492e40656SDave May 58592e40656SDave May /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ 58692e40656SDave May /* 0->cell, 1->edge, 2->vert */ 5879566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 58892e40656SDave May nel = pe - ps; 58992e40656SDave May 5909566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints * nel, -1)); 5919566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 5929566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 59392e40656SDave May 5949566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 5959566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 59692e40656SDave May 59792e40656SDave May pcnt = 0; 59892e40656SDave May for (e = 0; e < nel; e++) { 5999566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 60092e40656SDave May 60192e40656SDave May for (p = 0; p < npoints; p++) { 60292e40656SDave May for (d = 0; d < dim; d++) { 60392e40656SDave May swarm_coor[dim * pcnt + d] = 0.0; 604ad540459SPierre Jolivet for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][p * nbasis + k] * PetscRealPart(elcoor[dim * k + d]); 60592e40656SDave May } 60692e40656SDave May swarm_cellid[pcnt] = e; 60792e40656SDave May pcnt++; 60892e40656SDave May } 6099566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 61092e40656SDave May } 6119566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 6129566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 61392e40656SDave May 6149566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quadrature)); 6159566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 6163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 617431382f2SDave May } 618