xref: /petsc/src/dm/impls/swarm/swarmpic_plex.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 
82cf64d79SDave May static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem)
92cf64d79SDave May {
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 */
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P));
215f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) P, prefix));*/
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpacePolynomialSetTensor(P, tensor));
235f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscSpaceSetFromOptions(P));*/
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetDegree(P,1,PETSC_DETERMINE));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetNumComponents(P, Nc));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetNumVariables(P, dim));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetUp(P));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceGetDegree(P, &order, NULL));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpacePolynomialGetTensor(P, &tensor));
312cf64d79SDave May   /* Create dual space */
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE));
345f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix));*/
355f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetDM(Q, K));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&K));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetNumComponents(Q, Nc));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetOrder(Q, order));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceLagrangeSetTensor(Q, tensor));
415f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscDualSpaceSetFromOptions(Q));*/
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetUp(Q));
442cf64d79SDave May   /* Create element */
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreate(PetscObjectComm((PetscObject) dm), fem));
465f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix));*/
475f80ce2aSJacob Faibussowitsch   /*CHKERRQ(PetscFESetFromOptions(*fem));*/
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetType(*fem,PETSCFEBASIC));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetBasisSpace(*fem, P));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetDualSpace(*fem, Q));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetNumComponents(*fem, Nc));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetUp(*fem));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceDestroy(&P));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(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) {
595f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTStroudConicalQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q));
605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq));
612cf64d79SDave May   }
622cf64d79SDave May   else {
635f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q));
645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq));
652cf64d79SDave May   }
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetQuadrature(*fem, q));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetFaceQuadrature(*fem, fq));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&q));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&fq));
702cf64d79SDave May   PetscFunctionReturn(0);
712cf64d79SDave May }
722cf64d79SDave May 
736335e310SSatish Balay PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[2],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np)
740e2ec84fSDave May {
750e2ec84fSDave May   PetscReal      v12[2],v23[2],v31[2];
760e2ec84fSDave May   PetscInt       i;
770e2ec84fSDave May 
780e2ec84fSDave May   PetscFunctionBegin;
790e2ec84fSDave May   if (depth == max) {
800e2ec84fSDave May     PetscReal cx[2];
810e2ec84fSDave May 
820e2ec84fSDave May     cx[0] = (v1[0] + v2[0] + v3[0])/3.0;
830e2ec84fSDave May     cx[1] = (v1[1] + v2[1] + v3[1])/3.0;
840e2ec84fSDave May 
850e2ec84fSDave May     xi[2*(*np)+0] = cx[0];
860e2ec84fSDave May     xi[2*(*np)+1] = cx[1];
870e2ec84fSDave May     (*np)++;
880e2ec84fSDave May     PetscFunctionReturn(0);
890e2ec84fSDave May   }
900e2ec84fSDave May 
910e2ec84fSDave May   /* calculate midpoints of each side */
920e2ec84fSDave May   for (i = 0; i < 2; i++) {
930e2ec84fSDave May     v12[i] = (v1[i]+v2[i])/2.0;
940e2ec84fSDave May     v23[i] = (v2[i]+v3[i])/2.0;
950e2ec84fSDave May     v31[i] = (v3[i]+v1[i])/2.0;
960e2ec84fSDave May   }
970e2ec84fSDave May 
980e2ec84fSDave May   /* recursively subdivide new triangles */
995f80ce2aSJacob Faibussowitsch   CHKERRQ(subdivide_triangle(v1,v12,v31,depth+1,max,xi,np));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(subdivide_triangle(v2,v23,v12,depth+1,max,xi,np));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(subdivide_triangle(v3,v31,v23,depth+1,max,xi,np));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(subdivide_triangle(v12,v23,v31,depth+1,max,xi,np));
1030e2ec84fSDave May   PetscFunctionReturn(0);
1040e2ec84fSDave May }
1050e2ec84fSDave May 
1060e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub)
1070e2ec84fSDave May {
1080e2ec84fSDave May   const PetscInt dim = 2;
1090e2ec84fSDave May   PetscInt       q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth;
1100e2ec84fSDave May   PetscReal      *xi;
1110e2ec84fSDave May   PetscReal      **basis;
1120e2ec84fSDave May   Vec            coorlocal;
1130e2ec84fSDave May   PetscSection   coordSection;
1140e2ec84fSDave May   PetscScalar    *elcoor = NULL;
1150e2ec84fSDave May   PetscReal      *swarm_coor;
1160e2ec84fSDave May   PetscInt       *swarm_cellid;
1170e2ec84fSDave May   PetscReal      v1[2],v2[2],v3[2];
1180e2ec84fSDave May 
1190e2ec84fSDave May   PetscFunctionBegin;
1200e2ec84fSDave May   npoints_q = 1;
1210e2ec84fSDave May   for (d=0; d<nsub; d++) { npoints_q *= 4; }
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dim*npoints_q,&xi));
1230e2ec84fSDave May 
1240e2ec84fSDave May   v1[0] = 0.0;  v1[1] = 0.0;
1250e2ec84fSDave May   v2[0] = 1.0;  v2[1] = 0.0;
1260e2ec84fSDave May   v3[0] = 0.0;  v3[1] = 1.0;
1270e2ec84fSDave May   depth = 0;
1280e2ec84fSDave May   pcnt = 0;
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(subdivide_triangle(v1,v2,v3,depth,nsub,xi,&pcnt));
1300e2ec84fSDave May 
1310e2ec84fSDave May   npe = 3; /* nodes per element (triangle) */
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(npoints_q,&basis));
1330e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(npe,&basis[q]));
1350e2ec84fSDave May 
1360e2ec84fSDave May     basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1];
1370e2ec84fSDave May     basis[q][1] = xi[dim*q+0];
1380e2ec84fSDave May     basis[q][2] = xi[dim*q+1];
1390e2ec84fSDave May   }
1400e2ec84fSDave May 
1410e2ec84fSDave May   /* 0->cell, 1->edge, 2->vert */
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmc,0,&ps,&pe));
1430e2ec84fSDave May   nel = pe - ps;
1440e2ec84fSDave May 
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(dm,npoints_q*nel,-1));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
1480e2ec84fSDave May 
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dmc,&coorlocal));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dmc,&coordSection));
1510e2ec84fSDave May 
1520e2ec84fSDave May   pcnt = 0;
1530e2ec84fSDave May   for (e=0; e<nel; e++) {
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor));
1550e2ec84fSDave May 
1560e2ec84fSDave May     for (q=0; q<npoints_q; q++) {
1570e2ec84fSDave May       for (d=0; d<dim; d++) {
1580e2ec84fSDave May         swarm_coor[dim*pcnt+d] = 0.0;
1590e2ec84fSDave May         for (k=0; k<npe; k++) {
1600ca99263SDave May           swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]);
1610e2ec84fSDave May         }
1620e2ec84fSDave May       }
1630e2ec84fSDave May       swarm_cellid[pcnt] = e;
1640e2ec84fSDave May       pcnt++;
1650e2ec84fSDave May     }
1665f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor));
1670e2ec84fSDave May   }
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
1700e2ec84fSDave May 
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(xi));
1720e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(basis[q]));
1740e2ec84fSDave May   }
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(basis));
1760e2ec84fSDave May   PetscFunctionReturn(0);
1770e2ec84fSDave May }
1780e2ec84fSDave May 
17978185223SDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm,DM dmc,PetscInt nsub)
18078185223SDave May {
18178185223SDave May   PetscInt        dim,nfaces,nbasis;
18278185223SDave May   PetscInt        q,npoints_q,e,nel,pcnt,ps,pe,d,k,r;
183ef0bb6c7SMatthew G. Knepley   PetscTabulation T;
18478185223SDave May   Vec             coorlocal;
18578185223SDave May   PetscSection    coordSection;
18678185223SDave May   PetscScalar     *elcoor = NULL;
18778185223SDave May   PetscReal       *swarm_coor;
18878185223SDave May   PetscInt        *swarm_cellid;
18978185223SDave May   const PetscReal *xiq;
19078185223SDave May   PetscQuadrature quadrature;
19178185223SDave May   PetscFE         fe,feRef;
19278185223SDave May   PetscBool       is_simplex;
19378185223SDave May 
19478185223SDave May   PetscFunctionBegin;
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dmc,&dim));
19678185223SDave May   is_simplex = PETSC_FALSE;
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmc,0,&ps,&pe));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeSize(dmc, ps, &nfaces));
19978185223SDave May   if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
20078185223SDave May 
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe));
20278185223SDave May 
20378185223SDave May   for (r=0; r<nsub; r++) {
2045f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFERefine(fe,&feRef));
2055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFECopyQuadrature(feRef,fe));
2065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEDestroy(&feRef));
20778185223SDave May   }
20878185223SDave May 
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetQuadrature(fe,&quadrature));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetDimension(fe,&nbasis));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetCellTabulation(fe, 1, &T));
21378185223SDave May 
21478185223SDave May   /* 0->cell, 1->edge, 2->vert */
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmc,0,&ps,&pe));
21678185223SDave May   nel = pe - ps;
21778185223SDave May 
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(dm,npoints_q*nel,-1));
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
22178185223SDave May 
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dmc,&coorlocal));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dmc,&coordSection));
22478185223SDave May 
22578185223SDave May   pcnt = 0;
22678185223SDave May   for (e=0; e<nel; e++) {
2275f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor));
22878185223SDave May 
22978185223SDave May     for (q=0; q<npoints_q; q++) {
23078185223SDave May       for (d=0; d<dim; d++) {
23178185223SDave May         swarm_coor[dim*pcnt+d] = 0.0;
23278185223SDave May         for (k=0; k<nbasis; k++) {
233ef0bb6c7SMatthew G. Knepley           swarm_coor[dim*pcnt+d] += T->T[0][q*nbasis + k] * PetscRealPart(elcoor[dim*k+d]);
23478185223SDave May         }
23578185223SDave May       }
23678185223SDave May       swarm_cellid[pcnt] = e;
23778185223SDave May       pcnt++;
23878185223SDave May     }
2395f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor));
24078185223SDave May   }
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
2425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
24378185223SDave May 
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
24578185223SDave May   PetscFunctionReturn(0);
24678185223SDave May }
24778185223SDave May 
2480e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints)
2490e2ec84fSDave May {
250bc53056eSDave May   PetscInt       dim;
251bc53056eSDave May   PetscInt       ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,nfaces;
2520e2ec84fSDave May   PetscReal      *xi,ds,ds2;
2530e2ec84fSDave May   PetscReal      **basis;
2540e2ec84fSDave May   Vec            coorlocal;
2550e2ec84fSDave May   PetscSection   coordSection;
2560e2ec84fSDave May   PetscScalar    *elcoor = NULL;
2570e2ec84fSDave May   PetscReal      *swarm_coor;
2580e2ec84fSDave May   PetscInt       *swarm_cellid;
259bc53056eSDave May   PetscBool      is_simplex;
2600e2ec84fSDave May 
2610e2ec84fSDave May   PetscFunctionBegin;
2625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dmc,&dim));
2632c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dim != 2,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only 2D is supported");
264bc53056eSDave May   is_simplex = PETSC_FALSE;
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmc,0,&ps,&pe));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeSize(dmc, ps, &nfaces));
267bc53056eSDave May   if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
268*28b400f6SJacob Faibussowitsch   PetscCheck(is_simplex,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only the simplex is supported");
269bc53056eSDave May 
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dim*npoints*npoints,&xi));
2710e2ec84fSDave May   pcnt = 0;
2720e2ec84fSDave May   ds = 1.0/((PetscReal)(npoints-1));
2730e2ec84fSDave May   ds2 = 1.0/((PetscReal)(npoints));
2740e2ec84fSDave May   for (jj = 0; jj<npoints; jj++) {
2750e2ec84fSDave May     for (ii=0; ii<npoints-jj; ii++) {
2760e2ec84fSDave May       xi[dim*pcnt+0] = ii * ds;
2770e2ec84fSDave May       xi[dim*pcnt+1] = jj * ds;
2780e2ec84fSDave May 
2790e2ec84fSDave May       xi[dim*pcnt+0] *= (1.0 - 1.2*ds2);
2800e2ec84fSDave May       xi[dim*pcnt+1] *= (1.0 - 1.2*ds2);
2810e2ec84fSDave May 
2820e2ec84fSDave May       xi[dim*pcnt+0] += 0.35*ds2;
2830e2ec84fSDave May       xi[dim*pcnt+1] += 0.35*ds2;
2840e2ec84fSDave May       pcnt++;
2850e2ec84fSDave May     }
2860e2ec84fSDave May   }
2870e2ec84fSDave May   npoints_q = pcnt;
2880e2ec84fSDave May 
2890e2ec84fSDave May   npe = 3; /* nodes per element (triangle) */
2905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(npoints_q,&basis));
2910e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
2925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(npe,&basis[q]));
2930e2ec84fSDave May 
2940e2ec84fSDave May     basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1];
2950e2ec84fSDave May     basis[q][1] = xi[dim*q+0];
2960e2ec84fSDave May     basis[q][2] = xi[dim*q+1];
2970e2ec84fSDave May   }
2980e2ec84fSDave May 
2990e2ec84fSDave May   /* 0->cell, 1->edge, 2->vert */
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmc,0,&ps,&pe));
3010e2ec84fSDave May   nel = pe - ps;
3020e2ec84fSDave May 
3035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(dm,npoints_q*nel,-1));
3045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
3055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
3060e2ec84fSDave May 
3075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dmc,&coorlocal));
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dmc,&coordSection));
3090e2ec84fSDave May 
3100e2ec84fSDave May   pcnt = 0;
3110e2ec84fSDave May   for (e=0; e<nel; e++) {
3125f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor));
3130e2ec84fSDave May 
3140e2ec84fSDave May     for (q=0; q<npoints_q; q++) {
3150e2ec84fSDave May       for (d=0; d<dim; d++) {
3160e2ec84fSDave May         swarm_coor[dim*pcnt+d] = 0.0;
3170e2ec84fSDave May         for (k=0; k<npe; k++) {
3180ca99263SDave May           swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]);
3190e2ec84fSDave May         }
3200e2ec84fSDave May       }
3210e2ec84fSDave May       swarm_cellid[pcnt] = e;
3220e2ec84fSDave May       pcnt++;
3230e2ec84fSDave May     }
3245f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor));
3250e2ec84fSDave May   }
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
3275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
3280e2ec84fSDave May 
3295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(xi));
3300e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
3315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(basis[q]));
3320e2ec84fSDave May   }
3335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(basis));
3340e2ec84fSDave May   PetscFunctionReturn(0);
3350e2ec84fSDave May }
3360e2ec84fSDave May 
3370e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param)
3380e2ec84fSDave May {
3390e2ec84fSDave May   PetscInt       dim;
3400e2ec84fSDave May 
3410e2ec84fSDave May   PetscFunctionBegin;
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(celldm,&dim));
3430e2ec84fSDave May   switch (layout) {
3440e2ec84fSDave May     case DMSWARMPIC_LAYOUT_REGULAR:
3452c71b3e2SJacob Faibussowitsch       PetscCheckFalse(dim == 3,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No 3D support for REGULAR+PLEX");
3465f80ce2aSJacob Faibussowitsch       CHKERRQ(private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm,celldm,layout_param));
3470e2ec84fSDave May       break;
3480e2ec84fSDave May     case DMSWARMPIC_LAYOUT_GAUSS:
349b12d2206SDave May     {
350b12d2206SDave May       PetscInt npoints,npoints1,ps,pe,nfaces;
351b12d2206SDave May       const PetscReal *xi;
352b12d2206SDave May       PetscBool is_simplex;
353b12d2206SDave May       PetscQuadrature quadrature;
354b12d2206SDave May 
355b12d2206SDave May       is_simplex = PETSC_FALSE;
3565f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetHeightStratum(celldm,0,&ps,&pe));
3575f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetConeSize(celldm, ps, &nfaces));
358b12d2206SDave May       if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
359b12d2206SDave May 
360b12d2206SDave May       npoints1 = layout_param;
361b12d2206SDave May       if (is_simplex) {
3625f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscDTStroudConicalQuadrature(dim,1,npoints1,-1.0,1.0,&quadrature));
363b12d2206SDave May       } else {
3645f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscDTGaussTensorQuadrature(dim,1,npoints1,-1.0,1.0,&quadrature));
365b12d2206SDave May       }
3665f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscQuadratureGetData(quadrature,NULL,NULL,&npoints,&xi,NULL));
3675f80ce2aSJacob Faibussowitsch       CHKERRQ(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm,celldm,npoints,(PetscReal*)xi));
3685f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscQuadratureDestroy(&quadrature));
369b12d2206SDave May     }
3700e2ec84fSDave May       break;
3710e2ec84fSDave May     case DMSWARMPIC_LAYOUT_SUBDIVISION:
3725f80ce2aSJacob Faibussowitsch       CHKERRQ(private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm,celldm,layout_param));
3730e2ec84fSDave May       break;
3740e2ec84fSDave May   }
3750e2ec84fSDave May   PetscFunctionReturn(0);
3760e2ec84fSDave May }
377e3cd5995SDave May 
378e3cd5995SDave May /*
379e3cd5995SDave May typedef struct {
380e3cd5995SDave May   PetscReal x,y;
381e3cd5995SDave May } Point2d;
382e3cd5995SDave May 
383e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s)
384e3cd5995SDave May {
385e3cd5995SDave May   PetscFunctionBegin;
386e3cd5995SDave May   *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y);
387e3cd5995SDave May   PetscFunctionReturn(0);
388e3cd5995SDave May }
389e3cd5995SDave May */
390e3cd5995SDave May /*
391e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v)
392e3cd5995SDave May {
393e3cd5995SDave May   PetscReal s1,s2,s3;
394e3cd5995SDave May   PetscBool b1, b2, b3;
395e3cd5995SDave May 
396e3cd5995SDave May   PetscFunctionBegin;
397e3cd5995SDave May   signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f;
398e3cd5995SDave May   signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f;
399e3cd5995SDave May   signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f;
400e3cd5995SDave May 
401e3cd5995SDave May   *v = ((b1 == b2) && (b2 == b3));
402e3cd5995SDave May   PetscFunctionReturn(0);
403e3cd5995SDave May }
404e3cd5995SDave May */
405e3cd5995SDave May /*
406e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ)
407e3cd5995SDave May {
408e3cd5995SDave May   PetscReal x1,y1,x2,y2,x3,y3;
409e3cd5995SDave May   PetscReal c,b[2],A[2][2],inv[2][2],detJ,od;
410e3cd5995SDave May 
411e3cd5995SDave May   PetscFunctionBegin;
412e3cd5995SDave May   x1 = coords[2*0+0];
413e3cd5995SDave May   x2 = coords[2*1+0];
414e3cd5995SDave May   x3 = coords[2*2+0];
415e3cd5995SDave May 
416e3cd5995SDave May   y1 = coords[2*0+1];
417e3cd5995SDave May   y2 = coords[2*1+1];
418e3cd5995SDave May   y3 = coords[2*2+1];
419e3cd5995SDave May 
420e3cd5995SDave May   c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3;
421e3cd5995SDave May   b[0] = xp[0] - c;
422e3cd5995SDave May   c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3;
423e3cd5995SDave May   b[1] = xp[1] - c;
424e3cd5995SDave May 
425e3cd5995SDave May   A[0][0] = -0.5*x1 + 0.5*x2;   A[0][1] = -0.5*x1 + 0.5*x3;
426e3cd5995SDave May   A[1][0] = -0.5*y1 + 0.5*y2;   A[1][1] = -0.5*y1 + 0.5*y3;
427e3cd5995SDave May 
428e3cd5995SDave May   detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0];
429e3cd5995SDave May   *dJ = PetscAbsReal(detJ);
430e3cd5995SDave May   od = 1.0/detJ;
431e3cd5995SDave May 
432e3cd5995SDave May   inv[0][0] =  A[1][1] * od;
433e3cd5995SDave May   inv[0][1] = -A[0][1] * od;
434e3cd5995SDave May   inv[1][0] = -A[1][0] * od;
435e3cd5995SDave May   inv[1][1] =  A[0][0] * od;
436e3cd5995SDave May 
437e3cd5995SDave May   xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1];
438e3cd5995SDave May   xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1];
439e3cd5995SDave May   PetscFunctionReturn(0);
440e3cd5995SDave May }
441e3cd5995SDave May */
442e3cd5995SDave May 
443a5f152d1SDave May static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscScalar coords[],PetscReal xip[],PetscReal *dJ)
444e3cd5995SDave May {
445e3cd5995SDave May   PetscReal x1,y1,x2,y2,x3,y3;
446e3cd5995SDave May   PetscReal b[2],A[2][2],inv[2][2],detJ,od;
447e3cd5995SDave May 
448e3cd5995SDave May   PetscFunctionBegin;
449a5f152d1SDave May   x1 = PetscRealPart(coords[2*0+0]);
450a5f152d1SDave May   x2 = PetscRealPart(coords[2*1+0]);
451a5f152d1SDave May   x3 = PetscRealPart(coords[2*2+0]);
452e3cd5995SDave May 
453a5f152d1SDave May   y1 = PetscRealPart(coords[2*0+1]);
454a5f152d1SDave May   y2 = PetscRealPart(coords[2*1+1]);
455a5f152d1SDave May   y3 = PetscRealPart(coords[2*2+1]);
456e3cd5995SDave May 
457e3cd5995SDave May   b[0] = xp[0] - x1;
458e3cd5995SDave May   b[1] = xp[1] - y1;
459e3cd5995SDave May 
460e3cd5995SDave May   A[0][0] = x2-x1;   A[0][1] = x3-x1;
461e3cd5995SDave May   A[1][0] = y2-y1;   A[1][1] = y3-y1;
462e3cd5995SDave May 
463e3cd5995SDave May   detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0];
464e3cd5995SDave May   *dJ = PetscAbsReal(detJ);
465e3cd5995SDave May   od = 1.0/detJ;
466e3cd5995SDave May 
467e3cd5995SDave May   inv[0][0] =  A[1][1] * od;
468e3cd5995SDave May   inv[0][1] = -A[0][1] * od;
469e3cd5995SDave May   inv[1][0] = -A[1][0] * od;
470e3cd5995SDave May   inv[1][1] =  A[0][0] * od;
471e3cd5995SDave May 
472e3cd5995SDave May   xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1];
473e3cd5995SDave May   xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1];
474e3cd5995SDave May   PetscFunctionReturn(0);
475e3cd5995SDave May }
476e3cd5995SDave May 
477e3cd5995SDave May PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field)
478e3cd5995SDave May {
479e3cd5995SDave May   const PetscReal PLEX_C_EPS = 1.0e-8;
480e3cd5995SDave May   Vec             v_field_l,denom_l,coor_l,denom;
481e3cd5995SDave May   PetscInt        k,p,e,npoints;
482e3cd5995SDave May   PetscInt        *mpfield_cell;
483e3cd5995SDave May   PetscReal       *mpfield_coor;
484a5f152d1SDave May   PetscReal       xi_p[2];
485a5f152d1SDave May   PetscScalar     Ni[3];
486e3cd5995SDave May   PetscSection    coordSection;
487e3cd5995SDave May   PetscScalar     *elcoor = NULL;
488e3cd5995SDave May 
489e3cd5995SDave May   PetscFunctionBegin;
4905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(v_field));
491e3cd5995SDave May 
4925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm,&v_field_l));
4935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm,&denom));
4945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm,&denom_l));
4955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(v_field_l));
4965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(denom));
4975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(denom_l));
498e3cd5995SDave May 
4995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm,&coor_l));
5005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dm,&coordSection));
501e3cd5995SDave May 
5025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetLocalSize(swarm,&npoints));
5035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor));
5045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell));
505e3cd5995SDave May 
506e3cd5995SDave May   for (p=0; p<npoints; p++) {
507a5f152d1SDave May     PetscReal   *coor_p,dJ;
508a5f152d1SDave May     PetscScalar elfield[3];
509e3cd5995SDave May     PetscBool   point_located;
510e3cd5995SDave May 
511e3cd5995SDave May     e       = mpfield_cell[p];
512e3cd5995SDave May     coor_p  = &mpfield_coor[2*p];
513e3cd5995SDave May 
5145f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecGetClosure(dm,coordSection,coor_l,e,NULL,&elcoor));
515e3cd5995SDave May 
516e3cd5995SDave May /*
517e3cd5995SDave May     while (!point_located && (failed_counter < 25)) {
5185f80ce2aSJacob Faibussowitsch       CHKERRQ(PointInTriangle(point, coords[0], coords[1], coords[2], &point_located));
519e3cd5995SDave May       point.x = coor_p[0];
520e3cd5995SDave May       point.y = coor_p[1];
521e3cd5995SDave May       point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0);
522e3cd5995SDave May       point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0);
523e3cd5995SDave May       failed_counter++;
524e3cd5995SDave May     }
525e3cd5995SDave May 
526e3cd5995SDave May     if (!point_located) {
527e3cd5995SDave May         PetscPrintf(PETSC_COMM_SELF,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D) with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e) in %D 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);
528e3cd5995SDave May     }
529e3cd5995SDave May 
530*28b400f6SJacob Faibussowitsch     PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)",point.x,point.y,e);
531e3cd5995SDave May     else {
5325f80ce2aSJacob Faibussowitsch       CHKERRQ(_ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ));
533e3cd5995SDave May       xi_p[0] = 0.5*(xi_p[0] + 1.0);
534e3cd5995SDave May       xi_p[1] = 0.5*(xi_p[1] + 1.0);
535e3cd5995SDave May 
536e3cd5995SDave May       PetscPrintf(PETSC_COMM_SELF,"[p=%D] x(%+1.4e,%+1.4e) -> mapped to element %D xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]);
537e3cd5995SDave May 
538e3cd5995SDave May     }
539e3cd5995SDave May */
540e3cd5995SDave May 
5415f80ce2aSJacob Faibussowitsch     CHKERRQ(ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ));
542e3cd5995SDave May     /*
543e3cd5995SDave May     PetscPrintf(PETSC_COMM_SELF,"[p=%D] x(%+1.4e,%+1.4e) -> mapped to element %D xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]);
544e3cd5995SDave May     */
545e3cd5995SDave May     /*
546e3cd5995SDave May      point_located = PETSC_TRUE;
547e3cd5995SDave May     if (xi_p[0] < 0.0) {
548e3cd5995SDave May       if (xi_p[0] > -PLEX_C_EPS) {
549e3cd5995SDave May         xi_p[0] = 0.0;
550e3cd5995SDave May       } else {
551e3cd5995SDave May         point_located = PETSC_FALSE;
552e3cd5995SDave May       }
553e3cd5995SDave May     }
554e3cd5995SDave May     if (xi_p[1] < 0.0) {
555e3cd5995SDave May       if (xi_p[1] > -PLEX_C_EPS) {
556e3cd5995SDave May         xi_p[1] = 0.0;
557e3cd5995SDave May       } else {
558e3cd5995SDave May         point_located = PETSC_FALSE;
559e3cd5995SDave May       }
560e3cd5995SDave May     }
561e3cd5995SDave May     if (xi_p[1] > (1.0-xi_p[0])) {
562e3cd5995SDave May       if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) {
563e3cd5995SDave May         xi_p[1] = 1.0 - xi_p[0];
564e3cd5995SDave May       } else {
565e3cd5995SDave May         point_located = PETSC_FALSE;
566e3cd5995SDave May       }
567e3cd5995SDave May     }
568e3cd5995SDave May     if (!point_located) {
569e3cd5995SDave May       PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]);
570e3cd5995SDave May       PetscPrintf(PETSC_COMM_SELF,"[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D) 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]);
571e3cd5995SDave May     }
572*28b400f6SJacob Faibussowitsch     PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)",coor_p[0],coor_p[1],e);
573e3cd5995SDave May     */
574e3cd5995SDave May 
575e3cd5995SDave May     Ni[0] = 1.0 - xi_p[0] - xi_p[1];
576e3cd5995SDave May     Ni[1] = xi_p[0];
577e3cd5995SDave May     Ni[2] = xi_p[1];
578e3cd5995SDave May 
579e3cd5995SDave May     point_located = PETSC_TRUE;
580e3cd5995SDave May     for (k=0; k<3; k++) {
581a5f152d1SDave May       if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE;
582a5f152d1SDave May       if (PetscRealPart(Ni[k]) > (1.0+PLEX_C_EPS)) point_located = PETSC_FALSE;
583e3cd5995SDave May     }
584e3cd5995SDave May     if (!point_located) {
5855f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",(double)xi_p[0],(double)xi_p[1]));
5865f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D) 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])));
587e3cd5995SDave May     }
588*28b400f6SJacob Faibussowitsch     PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)",(double)coor_p[0],(double)coor_p[1],e);
589e3cd5995SDave May 
590e3cd5995SDave May     for (k=0; k<3; k++) {
591e3cd5995SDave May       Ni[k] = Ni[k] * dJ;
592e3cd5995SDave May       elfield[k] = Ni[k] * swarm_field[p];
593e3cd5995SDave May     }
5945f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecRestoreClosure(dm,coordSection,coor_l,e,NULL,&elcoor));
595e3cd5995SDave May 
5965f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecSetClosure(dm, NULL,v_field_l, e, elfield, ADD_VALUES));
5975f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecSetClosure(dm, NULL,denom_l, e, Ni, ADD_VALUES));
598e3cd5995SDave May   }
599e3cd5995SDave May 
6005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell));
6015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor));
602e3cd5995SDave May 
6035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field));
6045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field));
6055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom));
6065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom));
607e3cd5995SDave May 
6085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseDivide(v_field,v_field,denom));
609e3cd5995SDave May 
6105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm,&v_field_l));
6115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm,&denom_l));
6125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm,&denom));
613e3cd5995SDave May   PetscFunctionReturn(0);
614e3cd5995SDave May }
615e3cd5995SDave May 
61677048351SPatrick Sanan PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[])
617e3cd5995SDave May {
618e3cd5995SDave May   PetscInt       f,dim;
619e3cd5995SDave May 
620e3cd5995SDave May   PetscFunctionBegin;
6215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(swarm,&dim));
622e3cd5995SDave May   switch (dim) {
623e3cd5995SDave May     case 2:
624e3cd5995SDave May       for (f=0; f<nfields; f++) {
625e3cd5995SDave May         PetscReal *swarm_field;
626e3cd5995SDave May 
6275f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmDataFieldGetEntries(dfield[f],(void**)&swarm_field));
6285f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSwarmProjectField_ApproxP1_PLEX_2D(swarm,swarm_field,celldm,vecs[f]));
629e3cd5995SDave May       }
630e3cd5995SDave May       break;
631e3cd5995SDave May     case 3:
632e3cd5995SDave May       SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D");
633e3cd5995SDave May     default:
634e3cd5995SDave May       break;
635e3cd5995SDave May   }
636e3cd5995SDave May   PetscFunctionReturn(0);
637e3cd5995SDave May }
638431382f2SDave May 
639431382f2SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm,DM dmc,PetscInt npoints,PetscReal xi[])
640431382f2SDave May {
641431382f2SDave May   PetscBool       is_simplex,is_tensorcell;
64292e40656SDave May   PetscInt        dim,nfaces,ps,pe,p,d,nbasis,pcnt,e,k,nel;
64392e40656SDave May   PetscFE         fe;
64492e40656SDave May   PetscQuadrature quadrature;
645ef0bb6c7SMatthew G. Knepley   PetscTabulation T;
646ef0bb6c7SMatthew G. Knepley   PetscReal       *xiq;
64792e40656SDave May   Vec             coorlocal;
64892e40656SDave May   PetscSection    coordSection;
64992e40656SDave May   PetscScalar     *elcoor = NULL;
65092e40656SDave May   PetscReal       *swarm_coor;
65192e40656SDave May   PetscInt        *swarm_cellid;
652431382f2SDave May 
65392e40656SDave May   PetscFunctionBegin;
6545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dmc,&dim));
655431382f2SDave May 
656431382f2SDave May   is_simplex = PETSC_FALSE;
657431382f2SDave May   is_tensorcell = PETSC_FALSE;
6585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmc,0,&ps,&pe));
6595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeSize(dmc, ps, &nfaces));
660431382f2SDave May 
661431382f2SDave May   if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
662431382f2SDave May 
663431382f2SDave May   switch (dim) {
664431382f2SDave May     case 2:
665431382f2SDave May       if (nfaces == 4) { is_tensorcell = PETSC_TRUE; }
666431382f2SDave May       break;
667431382f2SDave May     case 3:
668431382f2SDave May       if (nfaces == 6) { is_tensorcell = PETSC_TRUE; }
669431382f2SDave May       break;
670431382f2SDave May     default:
671431382f2SDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for 2D, 3D");
672431382f2SDave May   }
673431382f2SDave May 
67492e40656SDave May   /* check points provided fail inside the reference cell */
675431382f2SDave May   if (is_simplex) {
67692e40656SDave May     for (p=0; p<npoints; p++) {
67792e40656SDave May       PetscReal sum;
67892e40656SDave May       for (d=0; d<dim; d++) {
6792c71b3e2SJacob Faibussowitsch         PetscCheckFalse(xi[dim*p+d] < -1.0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain");
68092e40656SDave May       }
68192e40656SDave May       sum = 0.0;
68292e40656SDave May       for (d=0; d<dim; d++) {
68392e40656SDave May         sum += xi[dim*p+d];
68492e40656SDave May       }
6852c71b3e2SJacob Faibussowitsch       PetscCheckFalse(sum > 0.0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain");
68692e40656SDave May     }
687431382f2SDave May   } else if (is_tensorcell) {
68892e40656SDave May     for (p=0; p<npoints; p++) {
68992e40656SDave May       for (d=0; d<dim; d++) {
6902c71b3e2SJacob Faibussowitsch         PetscCheckFalse(PetscAbsReal(xi[dim*p+d]) > 1.0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the tensor domain [-1,1]^d");
69192e40656SDave May       }
69292e40656SDave May     }
693431382f2SDave May   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for d-simplex and d-tensorcell");
694431382f2SDave May 
6955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm),&quadrature));
6965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(npoints*dim,&xiq));
6975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(xiq,xi,npoints*dim));
6985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureSetData(quadrature,dim,1,npoints,(const PetscReal*)xiq,NULL));
6995f80ce2aSJacob Faibussowitsch   CHKERRQ(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe));
7005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetQuadrature(fe,quadrature));
7015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetDimension(fe,&nbasis));
7025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetCellTabulation(fe, 1, &T));
70392e40656SDave May 
70492e40656SDave May   /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */
70592e40656SDave May   /* 0->cell, 1->edge, 2->vert */
7065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dmc,0,&ps,&pe));
70792e40656SDave May   nel = pe - ps;
70892e40656SDave May 
7095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(dm,npoints*nel,-1));
7105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
7115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
71292e40656SDave May 
7135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dmc,&coorlocal));
7145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dmc,&coordSection));
71592e40656SDave May 
71692e40656SDave May   pcnt = 0;
71792e40656SDave May   for (e=0; e<nel; e++) {
7185f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor));
71992e40656SDave May 
72092e40656SDave May     for (p=0; p<npoints; p++) {
72192e40656SDave May       for (d=0; d<dim; d++) {
72292e40656SDave May         swarm_coor[dim*pcnt+d] = 0.0;
72392e40656SDave May         for (k=0; k<nbasis; k++) {
724ef0bb6c7SMatthew G. Knepley           swarm_coor[dim*pcnt+d] += T->T[0][p*nbasis + k] * PetscRealPart(elcoor[dim*k+d]);
72592e40656SDave May         }
72692e40656SDave May       }
72792e40656SDave May       swarm_cellid[pcnt] = e;
72892e40656SDave May       pcnt++;
72992e40656SDave May     }
7305f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor));
73192e40656SDave May   }
7325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
7335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
73492e40656SDave May 
7355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&quadrature));
7365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
737431382f2SDave May   PetscFunctionReturn(0);
738431382f2SDave May }
739