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 89371c9d4SSatish Balay static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem) { 92cf64d79SDave May const PetscInt Nc = 1; 102cf64d79SDave May PetscQuadrature q, fq; 112cf64d79SDave May DM K; 122cf64d79SDave May PetscSpace P; 132cf64d79SDave May PetscDualSpace Q; 142cf64d79SDave May PetscInt order, quadPointsPerEdge; 152cf64d79SDave May PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 162cf64d79SDave May 172cf64d79SDave May PetscFunctionBegin; 182cf64d79SDave May /* Create space */ 199566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)dm), &P)); 209566063dSJacob Faibussowitsch /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix)); */ 219566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialSetTensor(P, tensor)); 229566063dSJacob Faibussowitsch /* PetscCall(PetscSpaceSetFromOptions(P)); */ 239566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 249566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(P, 1, PETSC_DETERMINE)); 259566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(P, Nc)); 269566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(P, dim)); 279566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(P)); 289566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(P, &order, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialGetTensor(P, &tensor)); 302cf64d79SDave May /* Create dual space */ 319566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)dm), &Q)); 329566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 339566063dSJacob Faibussowitsch /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix)); */ 349566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K)); 359566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(Q, K)); 369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K)); 379566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 389566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(Q, order)); 399566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTensor(Q, tensor)); 409566063dSJacob Faibussowitsch /* PetscCall(PetscDualSpaceSetFromOptions(Q)); */ 419566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 429566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(Q)); 432cf64d79SDave May /* Create element */ 449566063dSJacob Faibussowitsch PetscCall(PetscFECreate(PetscObjectComm((PetscObject)dm), fem)); 459566063dSJacob Faibussowitsch /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix)); */ 469566063dSJacob Faibussowitsch /* PetscCall(PetscFESetFromOptions(*fem)); */ 479566063dSJacob Faibussowitsch PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 489566063dSJacob Faibussowitsch PetscCall(PetscFESetBasisSpace(*fem, P)); 499566063dSJacob Faibussowitsch PetscCall(PetscFESetDualSpace(*fem, Q)); 509566063dSJacob Faibussowitsch PetscCall(PetscFESetNumComponents(*fem, Nc)); 519566063dSJacob Faibussowitsch PetscCall(PetscFESetUp(*fem)); 529566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&P)); 539566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&Q)); 542cf64d79SDave May /* Create quadrature (with specified order if given) */ 552cf64d79SDave May qorder = qorder >= 0 ? qorder : order; 562cf64d79SDave May quadPointsPerEdge = PetscMax(qorder + 1, 1); 572cf64d79SDave May if (isSimplex) { 589566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); 599566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); 609371c9d4SSatish Balay } else { 619566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); 629566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); 632cf64d79SDave May } 649566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(*fem, q)); 659566063dSJacob Faibussowitsch PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 669566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 679566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&fq)); 682cf64d79SDave May PetscFunctionReturn(0); 692cf64d79SDave May } 702cf64d79SDave May 719371c9d4SSatish Balay PetscErrorCode subdivide_triangle(PetscReal v1[2], PetscReal v2[2], PetscReal v3[2], PetscInt depth, PetscInt max, PetscReal xi[], PetscInt *np) { 720e2ec84fSDave May PetscReal v12[2], v23[2], v31[2]; 730e2ec84fSDave May PetscInt i; 740e2ec84fSDave May 750e2ec84fSDave May PetscFunctionBegin; 760e2ec84fSDave May if (depth == max) { 770e2ec84fSDave May PetscReal cx[2]; 780e2ec84fSDave May 790e2ec84fSDave May cx[0] = (v1[0] + v2[0] + v3[0]) / 3.0; 800e2ec84fSDave May cx[1] = (v1[1] + v2[1] + v3[1]) / 3.0; 810e2ec84fSDave May 820e2ec84fSDave May xi[2 * (*np) + 0] = cx[0]; 830e2ec84fSDave May xi[2 * (*np) + 1] = cx[1]; 840e2ec84fSDave May (*np)++; 850e2ec84fSDave May PetscFunctionReturn(0); 860e2ec84fSDave May } 870e2ec84fSDave May 880e2ec84fSDave May /* calculate midpoints of each side */ 890e2ec84fSDave May for (i = 0; i < 2; i++) { 900e2ec84fSDave May v12[i] = (v1[i] + v2[i]) / 2.0; 910e2ec84fSDave May v23[i] = (v2[i] + v3[i]) / 2.0; 920e2ec84fSDave May v31[i] = (v3[i] + v1[i]) / 2.0; 930e2ec84fSDave May } 940e2ec84fSDave May 950e2ec84fSDave May /* recursively subdivide new triangles */ 969566063dSJacob Faibussowitsch PetscCall(subdivide_triangle(v1, v12, v31, depth + 1, max, xi, np)); 979566063dSJacob Faibussowitsch PetscCall(subdivide_triangle(v2, v23, v12, depth + 1, max, xi, np)); 989566063dSJacob Faibussowitsch PetscCall(subdivide_triangle(v3, v31, v23, depth + 1, max, xi, np)); 999566063dSJacob Faibussowitsch PetscCall(subdivide_triangle(v12, v23, v31, depth + 1, max, xi, np)); 1000e2ec84fSDave May PetscFunctionReturn(0); 1010e2ec84fSDave May } 1020e2ec84fSDave May 1039371c9d4SSatish Balay PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm, DM dmc, PetscInt nsub) { 1040e2ec84fSDave May const PetscInt dim = 2; 1050e2ec84fSDave May PetscInt q, npoints_q, e, nel, npe, pcnt, ps, pe, d, k, depth; 1060e2ec84fSDave May PetscReal *xi; 1070e2ec84fSDave May PetscReal **basis; 1080e2ec84fSDave May Vec coorlocal; 1090e2ec84fSDave May PetscSection coordSection; 1100e2ec84fSDave May PetscScalar *elcoor = NULL; 1110e2ec84fSDave May PetscReal *swarm_coor; 1120e2ec84fSDave May PetscInt *swarm_cellid; 1130e2ec84fSDave May PetscReal v1[2], v2[2], v3[2]; 1140e2ec84fSDave May 1150e2ec84fSDave May PetscFunctionBegin; 1160e2ec84fSDave May npoints_q = 1; 1170e2ec84fSDave May for (d = 0; d < nsub; d++) { npoints_q *= 4; } 1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints_q, &xi)); 1190e2ec84fSDave May 1209371c9d4SSatish Balay v1[0] = 0.0; 1219371c9d4SSatish Balay v1[1] = 0.0; 1229371c9d4SSatish Balay v2[0] = 1.0; 1239371c9d4SSatish Balay v2[1] = 0.0; 1249371c9d4SSatish Balay v3[0] = 0.0; 1259371c9d4SSatish Balay v3[1] = 1.0; 1260e2ec84fSDave May depth = 0; 1270e2ec84fSDave May pcnt = 0; 1289566063dSJacob Faibussowitsch PetscCall(subdivide_triangle(v1, v2, v3, depth, nsub, xi, &pcnt)); 1290e2ec84fSDave May 1300e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 1319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints_q, &basis)); 1320e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 1339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npe, &basis[q])); 1340e2ec84fSDave May 1350e2ec84fSDave May basis[q][0] = 1.0 - xi[dim * q + 0] - xi[dim * q + 1]; 1360e2ec84fSDave May basis[q][1] = xi[dim * q + 0]; 1370e2ec84fSDave May basis[q][2] = xi[dim * q + 1]; 1380e2ec84fSDave May } 1390e2ec84fSDave May 1400e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 1419566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 1420e2ec84fSDave May nel = pe - ps; 1430e2ec84fSDave May 1449566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 1459566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1469566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1470e2ec84fSDave May 1489566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 1499566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 1500e2ec84fSDave May 1510e2ec84fSDave May pcnt = 0; 1520e2ec84fSDave May for (e = 0; e < nel; e++) { 1539566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 1540e2ec84fSDave May 1550e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 1560e2ec84fSDave May for (d = 0; d < dim; d++) { 1570e2ec84fSDave May swarm_coor[dim * pcnt + d] = 0.0; 1589371c9d4SSatish Balay for (k = 0; k < npe; k++) { swarm_coor[dim * pcnt + d] += basis[q][k] * PetscRealPart(elcoor[dim * k + d]); } 1590e2ec84fSDave May } 1600e2ec84fSDave May swarm_cellid[pcnt] = e; 1610e2ec84fSDave May pcnt++; 1620e2ec84fSDave May } 1639566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 1640e2ec84fSDave May } 1659566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1669566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1670e2ec84fSDave May 1689566063dSJacob Faibussowitsch PetscCall(PetscFree(xi)); 169*48a46eb9SPierre Jolivet for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q])); 1709566063dSJacob Faibussowitsch PetscCall(PetscFree(basis)); 1710e2ec84fSDave May PetscFunctionReturn(0); 1720e2ec84fSDave May } 1730e2ec84fSDave May 1749371c9d4SSatish Balay PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm, DM dmc, PetscInt nsub) { 17578185223SDave May PetscInt dim, nfaces, nbasis; 17678185223SDave May PetscInt q, npoints_q, e, nel, pcnt, ps, pe, d, k, r; 177ef0bb6c7SMatthew G. Knepley PetscTabulation T; 17878185223SDave May Vec coorlocal; 17978185223SDave May PetscSection coordSection; 18078185223SDave May PetscScalar *elcoor = NULL; 18178185223SDave May PetscReal *swarm_coor; 18278185223SDave May PetscInt *swarm_cellid; 18378185223SDave May const PetscReal *xiq; 18478185223SDave May PetscQuadrature quadrature; 18578185223SDave May PetscFE fe, feRef; 18678185223SDave May PetscBool is_simplex; 18778185223SDave May 18878185223SDave May PetscFunctionBegin; 1899566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim)); 19078185223SDave May is_simplex = PETSC_FALSE; 1919566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 1929566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 19378185223SDave May if (nfaces == (dim + 1)) { is_simplex = PETSC_TRUE; } 19478185223SDave May 1959566063dSJacob Faibussowitsch PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); 19678185223SDave May 19778185223SDave May for (r = 0; r < nsub; r++) { 1989566063dSJacob Faibussowitsch PetscCall(PetscFERefine(fe, &feRef)); 1999566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(feRef, fe)); 2009566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feRef)); 20178185223SDave May } 20278185223SDave May 2039566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 2049566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL)); 2059566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fe, &nbasis)); 2069566063dSJacob Faibussowitsch PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); 20778185223SDave May 20878185223SDave May /* 0->cell, 1->edge, 2->vert */ 2099566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 21078185223SDave May nel = pe - ps; 21178185223SDave May 2129566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 2139566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 2149566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 21578185223SDave May 2169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 2179566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 21878185223SDave May 21978185223SDave May pcnt = 0; 22078185223SDave May for (e = 0; e < nel; e++) { 2219566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 22278185223SDave May 22378185223SDave May for (q = 0; q < npoints_q; q++) { 22478185223SDave May for (d = 0; d < dim; d++) { 22578185223SDave May swarm_coor[dim * pcnt + d] = 0.0; 2269371c9d4SSatish Balay for (k = 0; k < nbasis; k++) { swarm_coor[dim * pcnt + d] += T->T[0][q * nbasis + k] * PetscRealPart(elcoor[dim * k + d]); } 22778185223SDave May } 22878185223SDave May swarm_cellid[pcnt] = e; 22978185223SDave May pcnt++; 23078185223SDave May } 2319566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 23278185223SDave May } 2339566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 2349566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 23578185223SDave May 2369566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 23778185223SDave May PetscFunctionReturn(0); 23878185223SDave May } 23978185223SDave May 2409371c9d4SSatish Balay PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm, DM dmc, PetscInt npoints) { 241bc53056eSDave May PetscInt dim; 242bc53056eSDave May PetscInt ii, jj, q, npoints_q, e, nel, npe, pcnt, ps, pe, d, k, nfaces; 2430e2ec84fSDave May PetscReal *xi, ds, ds2; 2440e2ec84fSDave May PetscReal **basis; 2450e2ec84fSDave May Vec coorlocal; 2460e2ec84fSDave May PetscSection coordSection; 2470e2ec84fSDave May PetscScalar *elcoor = NULL; 2480e2ec84fSDave May PetscReal *swarm_coor; 2490e2ec84fSDave May PetscInt *swarm_cellid; 250bc53056eSDave May PetscBool is_simplex; 2510e2ec84fSDave May 2520e2ec84fSDave May PetscFunctionBegin; 2539566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim)); 25408401ef6SPierre Jolivet PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only 2D is supported"); 255bc53056eSDave May is_simplex = PETSC_FALSE; 2569566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 2579566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 258bc53056eSDave May if (nfaces == (dim + 1)) { is_simplex = PETSC_TRUE; } 25928b400f6SJacob Faibussowitsch PetscCheck(is_simplex, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only the simplex is supported"); 260bc53056eSDave May 2619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints * npoints, &xi)); 2620e2ec84fSDave May pcnt = 0; 2630e2ec84fSDave May ds = 1.0 / ((PetscReal)(npoints - 1)); 2640e2ec84fSDave May ds2 = 1.0 / ((PetscReal)(npoints)); 2650e2ec84fSDave May for (jj = 0; jj < npoints; jj++) { 2660e2ec84fSDave May for (ii = 0; ii < npoints - jj; ii++) { 2670e2ec84fSDave May xi[dim * pcnt + 0] = ii * ds; 2680e2ec84fSDave May xi[dim * pcnt + 1] = jj * ds; 2690e2ec84fSDave May 2700e2ec84fSDave May xi[dim * pcnt + 0] *= (1.0 - 1.2 * ds2); 2710e2ec84fSDave May xi[dim * pcnt + 1] *= (1.0 - 1.2 * ds2); 2720e2ec84fSDave May 2730e2ec84fSDave May xi[dim * pcnt + 0] += 0.35 * ds2; 2740e2ec84fSDave May xi[dim * pcnt + 1] += 0.35 * ds2; 2750e2ec84fSDave May pcnt++; 2760e2ec84fSDave May } 2770e2ec84fSDave May } 2780e2ec84fSDave May npoints_q = pcnt; 2790e2ec84fSDave May 2800e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 2819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints_q, &basis)); 2820e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npe, &basis[q])); 2840e2ec84fSDave May 2850e2ec84fSDave May basis[q][0] = 1.0 - xi[dim * q + 0] - xi[dim * q + 1]; 2860e2ec84fSDave May basis[q][1] = xi[dim * q + 0]; 2870e2ec84fSDave May basis[q][2] = xi[dim * q + 1]; 2880e2ec84fSDave May } 2890e2ec84fSDave May 2900e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 2919566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 2920e2ec84fSDave May nel = pe - ps; 2930e2ec84fSDave May 2949566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 2959566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 2969566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 2970e2ec84fSDave May 2989566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 2999566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 3000e2ec84fSDave May 3010e2ec84fSDave May pcnt = 0; 3020e2ec84fSDave May for (e = 0; e < nel; e++) { 3039566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 3040e2ec84fSDave May 3050e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 3060e2ec84fSDave May for (d = 0; d < dim; d++) { 3070e2ec84fSDave May swarm_coor[dim * pcnt + d] = 0.0; 3089371c9d4SSatish Balay for (k = 0; k < npe; k++) { swarm_coor[dim * pcnt + d] += basis[q][k] * PetscRealPart(elcoor[dim * k + d]); } 3090e2ec84fSDave May } 3100e2ec84fSDave May swarm_cellid[pcnt] = e; 3110e2ec84fSDave May pcnt++; 3120e2ec84fSDave May } 3139566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor)); 3140e2ec84fSDave May } 3159566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 3169566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 3170e2ec84fSDave May 3189566063dSJacob Faibussowitsch PetscCall(PetscFree(xi)); 319*48a46eb9SPierre Jolivet for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q])); 3209566063dSJacob Faibussowitsch PetscCall(PetscFree(basis)); 3210e2ec84fSDave May PetscFunctionReturn(0); 3220e2ec84fSDave May } 3230e2ec84fSDave May 3249371c9d4SSatish Balay PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param) { 3250e2ec84fSDave May PetscInt dim; 3260e2ec84fSDave May 3270e2ec84fSDave May PetscFunctionBegin; 3289566063dSJacob Faibussowitsch PetscCall(DMGetDimension(celldm, &dim)); 3290e2ec84fSDave May switch (layout) { 3300e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 33108401ef6SPierre Jolivet PetscCheck(dim != 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No 3D support for REGULAR+PLEX"); 3329566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm, celldm, layout_param)); 3330e2ec84fSDave May break; 3349371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_GAUSS: { 335b12d2206SDave May PetscInt npoints, npoints1, ps, pe, nfaces; 336b12d2206SDave May const PetscReal *xi; 337b12d2206SDave May PetscBool is_simplex; 338b12d2206SDave May PetscQuadrature quadrature; 339b12d2206SDave May 340b12d2206SDave May is_simplex = PETSC_FALSE; 3419566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe)); 3429566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(celldm, ps, &nfaces)); 343b12d2206SDave May if (nfaces == (dim + 1)) { is_simplex = PETSC_TRUE; } 344b12d2206SDave May 345b12d2206SDave May npoints1 = layout_param; 346b12d2206SDave May if (is_simplex) { 3479566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, npoints1, -1.0, 1.0, &quadrature)); 348b12d2206SDave May } else { 3499566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, npoints1, -1.0, 1.0, &quadrature)); 350b12d2206SDave May } 3519566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints, &xi, NULL)); 3529566063dSJacob Faibussowitsch PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, (PetscReal *)xi)); 3539566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quadrature)); 3549371c9d4SSatish Balay } break; 3559371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_SUBDIVISION: PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm, celldm, layout_param)); break; 3560e2ec84fSDave May } 3570e2ec84fSDave May PetscFunctionReturn(0); 3580e2ec84fSDave May } 359e3cd5995SDave May 360e3cd5995SDave May /* 361e3cd5995SDave May typedef struct { 362e3cd5995SDave May PetscReal x,y; 363e3cd5995SDave May } Point2d; 364e3cd5995SDave May 365e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s) 366e3cd5995SDave May { 367e3cd5995SDave May PetscFunctionBegin; 368e3cd5995SDave May *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y); 369e3cd5995SDave May PetscFunctionReturn(0); 370e3cd5995SDave May } 371e3cd5995SDave May */ 372e3cd5995SDave May /* 373e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v) 374e3cd5995SDave May { 375e3cd5995SDave May PetscReal s1,s2,s3; 376e3cd5995SDave May PetscBool b1, b2, b3; 377e3cd5995SDave May 378e3cd5995SDave May PetscFunctionBegin; 379e3cd5995SDave May signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f; 380e3cd5995SDave May signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f; 381e3cd5995SDave May signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f; 382e3cd5995SDave May 383e3cd5995SDave May *v = ((b1 == b2) && (b2 == b3)); 384e3cd5995SDave May PetscFunctionReturn(0); 385e3cd5995SDave May } 386e3cd5995SDave May */ 387e3cd5995SDave May /* 388e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ) 389e3cd5995SDave May { 390e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 391e3cd5995SDave May PetscReal c,b[2],A[2][2],inv[2][2],detJ,od; 392e3cd5995SDave May 393e3cd5995SDave May PetscFunctionBegin; 394e3cd5995SDave May x1 = coords[2*0+0]; 395e3cd5995SDave May x2 = coords[2*1+0]; 396e3cd5995SDave May x3 = coords[2*2+0]; 397e3cd5995SDave May 398e3cd5995SDave May y1 = coords[2*0+1]; 399e3cd5995SDave May y2 = coords[2*1+1]; 400e3cd5995SDave May y3 = coords[2*2+1]; 401e3cd5995SDave May 402e3cd5995SDave May c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3; 403e3cd5995SDave May b[0] = xp[0] - c; 404e3cd5995SDave May c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3; 405e3cd5995SDave May b[1] = xp[1] - c; 406e3cd5995SDave May 407e3cd5995SDave May A[0][0] = -0.5*x1 + 0.5*x2; A[0][1] = -0.5*x1 + 0.5*x3; 408e3cd5995SDave May A[1][0] = -0.5*y1 + 0.5*y2; A[1][1] = -0.5*y1 + 0.5*y3; 409e3cd5995SDave May 410e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 411e3cd5995SDave May *dJ = PetscAbsReal(detJ); 412e3cd5995SDave May od = 1.0/detJ; 413e3cd5995SDave May 414e3cd5995SDave May inv[0][0] = A[1][1] * od; 415e3cd5995SDave May inv[0][1] = -A[0][1] * od; 416e3cd5995SDave May inv[1][0] = -A[1][0] * od; 417e3cd5995SDave May inv[1][1] = A[0][0] * od; 418e3cd5995SDave May 419e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 420e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 421e3cd5995SDave May PetscFunctionReturn(0); 422e3cd5995SDave May } 423e3cd5995SDave May */ 424e3cd5995SDave May 4259371c9d4SSatish Balay static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[], PetscScalar coords[], PetscReal xip[], PetscReal *dJ) { 426e3cd5995SDave May PetscReal x1, y1, x2, y2, x3, y3; 427e3cd5995SDave May PetscReal b[2], A[2][2], inv[2][2], detJ, od; 428e3cd5995SDave May 429e3cd5995SDave May PetscFunctionBegin; 430a5f152d1SDave May x1 = PetscRealPart(coords[2 * 0 + 0]); 431a5f152d1SDave May x2 = PetscRealPart(coords[2 * 1 + 0]); 432a5f152d1SDave May x3 = PetscRealPart(coords[2 * 2 + 0]); 433e3cd5995SDave May 434a5f152d1SDave May y1 = PetscRealPart(coords[2 * 0 + 1]); 435a5f152d1SDave May y2 = PetscRealPart(coords[2 * 1 + 1]); 436a5f152d1SDave May y3 = PetscRealPart(coords[2 * 2 + 1]); 437e3cd5995SDave May 438e3cd5995SDave May b[0] = xp[0] - x1; 439e3cd5995SDave May b[1] = xp[1] - y1; 440e3cd5995SDave May 4419371c9d4SSatish Balay A[0][0] = x2 - x1; 4429371c9d4SSatish Balay A[0][1] = x3 - x1; 4439371c9d4SSatish Balay A[1][0] = y2 - y1; 4449371c9d4SSatish Balay A[1][1] = y3 - y1; 445e3cd5995SDave May 446e3cd5995SDave May detJ = A[0][0] * A[1][1] - A[0][1] * A[1][0]; 447e3cd5995SDave May *dJ = PetscAbsReal(detJ); 448e3cd5995SDave May od = 1.0 / detJ; 449e3cd5995SDave May 450e3cd5995SDave May inv[0][0] = A[1][1] * od; 451e3cd5995SDave May inv[0][1] = -A[0][1] * od; 452e3cd5995SDave May inv[1][0] = -A[1][0] * od; 453e3cd5995SDave May inv[1][1] = A[0][0] * od; 454e3cd5995SDave May 455e3cd5995SDave May xip[0] = inv[0][0] * b[0] + inv[0][1] * b[1]; 456e3cd5995SDave May xip[1] = inv[1][0] * b[0] + inv[1][1] * b[1]; 457e3cd5995SDave May PetscFunctionReturn(0); 458e3cd5995SDave May } 459e3cd5995SDave May 4609371c9d4SSatish Balay PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm, PetscReal *swarm_field, DM dm, Vec v_field) { 461e3cd5995SDave May const PetscReal PLEX_C_EPS = 1.0e-8; 462e3cd5995SDave May Vec v_field_l, denom_l, coor_l, denom; 463e3cd5995SDave May PetscInt k, p, e, npoints; 464e3cd5995SDave May PetscInt *mpfield_cell; 465e3cd5995SDave May PetscReal *mpfield_coor; 466a5f152d1SDave May PetscReal xi_p[2]; 467a5f152d1SDave May PetscScalar Ni[3]; 468e3cd5995SDave May PetscSection coordSection; 469e3cd5995SDave May PetscScalar *elcoor = NULL; 470e3cd5995SDave May 471e3cd5995SDave May PetscFunctionBegin; 4729566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field)); 473e3cd5995SDave May 4749566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v_field_l)); 4759566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &denom)); 4769566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &denom_l)); 4779566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field_l)); 4789566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom)); 4799566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom_l)); 480e3cd5995SDave May 4819566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coor_l)); 4829566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 483e3cd5995SDave May 4849566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(swarm, &npoints)); 4859566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 4869566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 487e3cd5995SDave May 488e3cd5995SDave May for (p = 0; p < npoints; p++) { 489a5f152d1SDave May PetscReal *coor_p, dJ; 490a5f152d1SDave May PetscScalar elfield[3]; 491e3cd5995SDave May PetscBool point_located; 492e3cd5995SDave May 493e3cd5995SDave May e = mpfield_cell[p]; 494e3cd5995SDave May coor_p = &mpfield_coor[2 * p]; 495e3cd5995SDave May 4969566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 497e3cd5995SDave May 498e3cd5995SDave May /* 499e3cd5995SDave May while (!point_located && (failed_counter < 25)) { 5009566063dSJacob Faibussowitsch PetscCall(PointInTriangle(point, coords[0], coords[1], coords[2], &point_located)); 501e3cd5995SDave May point.x = coor_p[0]; 502e3cd5995SDave May point.y = coor_p[1]; 503e3cd5995SDave May point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 504e3cd5995SDave May point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 505e3cd5995SDave May failed_counter++; 506e3cd5995SDave May } 507e3cd5995SDave May 508e3cd5995SDave May if (!point_located) { 50963a3b9bcSJacob 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); 510e3cd5995SDave May } 511e3cd5995SDave May 51263a3b9bcSJacob 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); 513e3cd5995SDave May else { 5149566063dSJacob Faibussowitsch PetscCall(_ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ)); 515e3cd5995SDave May xi_p[0] = 0.5*(xi_p[0] + 1.0); 516e3cd5995SDave May xi_p[1] = 0.5*(xi_p[1] + 1.0); 517e3cd5995SDave May 51863a3b9bcSJacob 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]); 519e3cd5995SDave May 520e3cd5995SDave May } 521e3cd5995SDave May */ 522e3cd5995SDave May 5239566063dSJacob Faibussowitsch PetscCall(ComputeLocalCoordinateAffine2d(coor_p, elcoor, xi_p, &dJ)); 524e3cd5995SDave May /* 52563a3b9bcSJacob 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]); 526e3cd5995SDave May */ 527e3cd5995SDave May /* 528e3cd5995SDave May point_located = PETSC_TRUE; 529e3cd5995SDave May if (xi_p[0] < 0.0) { 530e3cd5995SDave May if (xi_p[0] > -PLEX_C_EPS) { 531e3cd5995SDave May xi_p[0] = 0.0; 532e3cd5995SDave May } else { 533e3cd5995SDave May point_located = PETSC_FALSE; 534e3cd5995SDave May } 535e3cd5995SDave May } 536e3cd5995SDave May if (xi_p[1] < 0.0) { 537e3cd5995SDave May if (xi_p[1] > -PLEX_C_EPS) { 538e3cd5995SDave May xi_p[1] = 0.0; 539e3cd5995SDave May } else { 540e3cd5995SDave May point_located = PETSC_FALSE; 541e3cd5995SDave May } 542e3cd5995SDave May } 543e3cd5995SDave May if (xi_p[1] > (1.0-xi_p[0])) { 544e3cd5995SDave May if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { 545e3cd5995SDave May xi_p[1] = 1.0 - xi_p[0]; 546e3cd5995SDave May } else { 547e3cd5995SDave May point_located = PETSC_FALSE; 548e3cd5995SDave May } 549e3cd5995SDave May } 550e3cd5995SDave May if (!point_located) { 551e3cd5995SDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); 55263a3b9bcSJacob 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]); 553e3cd5995SDave May } 55463a3b9bcSJacob 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); 555e3cd5995SDave May */ 556e3cd5995SDave May 557e3cd5995SDave May Ni[0] = 1.0 - xi_p[0] - xi_p[1]; 558e3cd5995SDave May Ni[1] = xi_p[0]; 559e3cd5995SDave May Ni[2] = xi_p[1]; 560e3cd5995SDave May 561e3cd5995SDave May point_located = PETSC_TRUE; 562e3cd5995SDave May for (k = 0; k < 3; k++) { 563a5f152d1SDave May if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; 564a5f152d1SDave May if (PetscRealPart(Ni[k]) > (1.0 + PLEX_C_EPS)) point_located = PETSC_FALSE; 565e3cd5995SDave May } 566e3cd5995SDave May if (!point_located) { 5679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[Error] xi,eta = %+1.8e, %+1.8e\n", (double)xi_p[0], (double)xi_p[1])); 56863a3b9bcSJacob 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]))); 569e3cd5995SDave May } 57063a3b9bcSJacob 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); 571e3cd5995SDave May 572e3cd5995SDave May for (k = 0; k < 3; k++) { 573e3cd5995SDave May Ni[k] = Ni[k] * dJ; 574e3cd5995SDave May elfield[k] = Ni[k] * swarm_field[p]; 575e3cd5995SDave May } 5769566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coor_l, e, NULL, &elcoor)); 577e3cd5995SDave May 5789566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, v_field_l, e, elfield, ADD_VALUES)); 5799566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, denom_l, e, Ni, ADD_VALUES)); 580e3cd5995SDave May } 581e3cd5995SDave May 5829566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 5839566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 584e3cd5995SDave May 5859566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, v_field_l, ADD_VALUES, v_field)); 5869566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, v_field_l, ADD_VALUES, v_field)); 5879566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, denom_l, ADD_VALUES, denom)); 5889566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, denom_l, ADD_VALUES, denom)); 589e3cd5995SDave May 5909566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(v_field, v_field, denom)); 591e3cd5995SDave May 5929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v_field_l)); 5939566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &denom_l)); 5949566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &denom)); 595e3cd5995SDave May PetscFunctionReturn(0); 596e3cd5995SDave May } 597e3cd5995SDave May 5989371c9d4SSatish Balay PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]) { 599e3cd5995SDave May PetscInt f, dim; 600e3cd5995SDave May 601e3cd5995SDave May PetscFunctionBegin; 6029566063dSJacob Faibussowitsch PetscCall(DMGetDimension(swarm, &dim)); 603e3cd5995SDave May switch (dim) { 604e3cd5995SDave May case 2: 605e3cd5995SDave May for (f = 0; f < nfields; f++) { 606e3cd5995SDave May PetscReal *swarm_field; 607e3cd5995SDave May 6089566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(dfield[f], (void **)&swarm_field)); 6099566063dSJacob Faibussowitsch PetscCall(DMSwarmProjectField_ApproxP1_PLEX_2D(swarm, swarm_field, celldm, vecs[f])); 610e3cd5995SDave May } 611e3cd5995SDave May break; 6129371c9d4SSatish Balay case 3: SETERRQ(PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "No support for 3D"); 6139371c9d4SSatish Balay default: break; 614e3cd5995SDave May } 615e3cd5995SDave May PetscFunctionReturn(0); 616e3cd5995SDave May } 617431382f2SDave May 6189371c9d4SSatish Balay PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm, DM dmc, PetscInt npoints, PetscReal xi[]) { 619431382f2SDave May PetscBool is_simplex, is_tensorcell; 62092e40656SDave May PetscInt dim, nfaces, ps, pe, p, d, nbasis, pcnt, e, k, nel; 62192e40656SDave May PetscFE fe; 62292e40656SDave May PetscQuadrature quadrature; 623ef0bb6c7SMatthew G. Knepley PetscTabulation T; 624ef0bb6c7SMatthew G. Knepley PetscReal *xiq; 62592e40656SDave May Vec coorlocal; 62692e40656SDave May PetscSection coordSection; 62792e40656SDave May PetscScalar *elcoor = NULL; 62892e40656SDave May PetscReal *swarm_coor; 62992e40656SDave May PetscInt *swarm_cellid; 630431382f2SDave May 63192e40656SDave May PetscFunctionBegin; 6329566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmc, &dim)); 633431382f2SDave May 634431382f2SDave May is_simplex = PETSC_FALSE; 635431382f2SDave May is_tensorcell = PETSC_FALSE; 6369566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 6379566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); 638431382f2SDave May 639431382f2SDave May if (nfaces == (dim + 1)) { is_simplex = PETSC_TRUE; } 640431382f2SDave May 641431382f2SDave May switch (dim) { 642431382f2SDave May case 2: 643431382f2SDave May if (nfaces == 4) { is_tensorcell = PETSC_TRUE; } 644431382f2SDave May break; 645431382f2SDave May case 3: 646431382f2SDave May if (nfaces == 6) { is_tensorcell = PETSC_TRUE; } 647431382f2SDave May break; 6489371c9d4SSatish Balay default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for 2D, 3D"); 649431382f2SDave May } 650431382f2SDave May 65192e40656SDave May /* check points provided fail inside the reference cell */ 652431382f2SDave May if (is_simplex) { 65392e40656SDave May for (p = 0; p < npoints; p++) { 65492e40656SDave May PetscReal sum; 6559371c9d4SSatish Balay 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"); } 65692e40656SDave May sum = 0.0; 6579371c9d4SSatish Balay for (d = 0; d < dim; d++) { sum += xi[dim * p + d]; } 65808401ef6SPierre Jolivet PetscCheck(sum <= 0.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain"); 65992e40656SDave May } 660431382f2SDave May } else if (is_tensorcell) { 66192e40656SDave May for (p = 0; p < npoints; p++) { 6629371c9d4SSatish Balay 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"); } 66392e40656SDave May } 664431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for d-simplex and d-tensorcell"); 665431382f2SDave May 6669566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm), &quadrature)); 6679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints * dim, &xiq)); 6689566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(xiq, xi, npoints * dim)); 6699566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(quadrature, dim, 1, npoints, (const PetscReal *)xiq, NULL)); 6709566063dSJacob Faibussowitsch PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); 6719566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(fe, quadrature)); 6729566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fe, &nbasis)); 6739566063dSJacob Faibussowitsch PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); 67492e40656SDave May 67592e40656SDave May /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ 67692e40656SDave May /* 0->cell, 1->edge, 2->vert */ 6779566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe)); 67892e40656SDave May nel = pe - ps; 67992e40656SDave May 6809566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints * nel, -1)); 6819566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 6829566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 68392e40656SDave May 6849566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal)); 6859566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmc, &coordSection)); 68692e40656SDave May 68792e40656SDave May pcnt = 0; 68892e40656SDave May for (e = 0; e < nel; e++) { 6899566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 69092e40656SDave May 69192e40656SDave May for (p = 0; p < npoints; p++) { 69292e40656SDave May for (d = 0; d < dim; d++) { 69392e40656SDave May swarm_coor[dim * pcnt + d] = 0.0; 6949371c9d4SSatish Balay for (k = 0; k < nbasis; k++) { swarm_coor[dim * pcnt + d] += T->T[0][p * nbasis + k] * PetscRealPart(elcoor[dim * k + d]); } 69592e40656SDave May } 69692e40656SDave May swarm_cellid[pcnt] = e; 69792e40656SDave May pcnt++; 69892e40656SDave May } 6999566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor)); 70092e40656SDave May } 7019566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 7029566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 70392e40656SDave May 7049566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quadrature)); 7059566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 706431382f2SDave May PetscFunctionReturn(0); 707431382f2SDave May } 708