xref: /petsc/src/dm/impls/swarm/swarmpic_plex.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode  ierr;
182cf64d79SDave May 
192cf64d79SDave May   PetscFunctionBegin;
202cf64d79SDave May   /* Create space */
212cf64d79SDave May   ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P);CHKERRQ(ierr);
222cf64d79SDave May   /*ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);*/
232cf64d79SDave May   ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr);
242cf64d79SDave May   /*ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);*/
252cf64d79SDave May   ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
26d39dd5f5SToby Isaac   ierr = PetscSpaceSetDegree(P,1,PETSC_DETERMINE);CHKERRQ(ierr);
272cf64d79SDave May   ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr);
28157782e2SToby Isaac   ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr);
292cf64d79SDave May   ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
30cdcc6362SToby Isaac   ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr);
312cf64d79SDave May   ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr);
322cf64d79SDave May   /* Create dual space */
332cf64d79SDave May   ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q);CHKERRQ(ierr);
342cf64d79SDave May   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
352cf64d79SDave May   /*ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);*/
362cf64d79SDave May   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr);
372cf64d79SDave May   ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr);
382cf64d79SDave May   ierr = DMDestroy(&K);CHKERRQ(ierr);
392cf64d79SDave May   ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr);
402cf64d79SDave May   ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr);
412cf64d79SDave May   ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr);
422cf64d79SDave May   /*ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);*/
432cf64d79SDave May   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
442cf64d79SDave May   ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
452cf64d79SDave May   /* Create element */
462cf64d79SDave May   ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem);CHKERRQ(ierr);
472cf64d79SDave May   /*ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);*/
482cf64d79SDave May   /*ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);*/
492cf64d79SDave May   ierr = PetscFESetType(*fem,PETSCFEBASIC);CHKERRQ(ierr);
502cf64d79SDave May   ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr);
512cf64d79SDave May   ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr);
522cf64d79SDave May   ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr);
532cf64d79SDave May   ierr = PetscFESetUp(*fem);CHKERRQ(ierr);
542cf64d79SDave May   ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
552cf64d79SDave May   ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
562cf64d79SDave May   /* Create quadrature (with specified order if given) */
572cf64d79SDave May   qorder = qorder >= 0 ? qorder : order;
582cf64d79SDave May   quadPointsPerEdge = PetscMax(qorder + 1,1);
592cf64d79SDave May   if (isSimplex) {
60e6a796c3SToby Isaac     ierr = PetscDTStroudConicalQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
61e6a796c3SToby Isaac     ierr = PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
622cf64d79SDave May   }
632cf64d79SDave May   else {
642cf64d79SDave May     ierr = PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
652cf64d79SDave May     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
662cf64d79SDave May   }
672cf64d79SDave May   ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr);
682cf64d79SDave May   ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr);
692cf64d79SDave May   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
702cf64d79SDave May   ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr);
712cf64d79SDave May   PetscFunctionReturn(0);
722cf64d79SDave May }
732cf64d79SDave May 
746335e310SSatish Balay PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[2],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np)
750e2ec84fSDave May {
760e2ec84fSDave May   PetscReal      v12[2],v23[2],v31[2];
770e2ec84fSDave May   PetscInt       i;
780e2ec84fSDave May   PetscErrorCode ierr;
790e2ec84fSDave May 
800e2ec84fSDave May   PetscFunctionBegin;
810e2ec84fSDave May   if (depth == max) {
820e2ec84fSDave May     PetscReal cx[2];
830e2ec84fSDave May 
840e2ec84fSDave May     cx[0] = (v1[0] + v2[0] + v3[0])/3.0;
850e2ec84fSDave May     cx[1] = (v1[1] + v2[1] + v3[1])/3.0;
860e2ec84fSDave May 
870e2ec84fSDave May     xi[2*(*np)+0] = cx[0];
880e2ec84fSDave May     xi[2*(*np)+1] = cx[1];
890e2ec84fSDave May     (*np)++;
900e2ec84fSDave May     PetscFunctionReturn(0);
910e2ec84fSDave May   }
920e2ec84fSDave May 
930e2ec84fSDave May   /* calculate midpoints of each side */
940e2ec84fSDave May   for (i = 0; i < 2; i++) {
950e2ec84fSDave May     v12[i] = (v1[i]+v2[i])/2.0;
960e2ec84fSDave May     v23[i] = (v2[i]+v3[i])/2.0;
970e2ec84fSDave May     v31[i] = (v3[i]+v1[i])/2.0;
980e2ec84fSDave May   }
990e2ec84fSDave May 
1000e2ec84fSDave May   /* recursively subdivide new triangles */
1010e2ec84fSDave May   ierr = subdivide_triangle(v1,v12,v31,depth+1,max,xi,np);CHKERRQ(ierr);
1020e2ec84fSDave May   ierr = subdivide_triangle(v2,v23,v12,depth+1,max,xi,np);CHKERRQ(ierr);
1030e2ec84fSDave May   ierr = subdivide_triangle(v3,v31,v23,depth+1,max,xi,np);CHKERRQ(ierr);
1040e2ec84fSDave May   ierr = subdivide_triangle(v12,v23,v31,depth+1,max,xi,np);CHKERRQ(ierr);
1050e2ec84fSDave May   PetscFunctionReturn(0);
1060e2ec84fSDave May }
1070e2ec84fSDave May 
1080e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub)
1090e2ec84fSDave May {
1100e2ec84fSDave May   PetscErrorCode ierr;
1110e2ec84fSDave May   const PetscInt dim = 2;
1120e2ec84fSDave May   PetscInt       q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth;
1130e2ec84fSDave May   PetscReal      *xi;
1140e2ec84fSDave May   PetscReal      **basis;
1150e2ec84fSDave May   Vec            coorlocal;
1160e2ec84fSDave May   PetscSection   coordSection;
1170e2ec84fSDave May   PetscScalar    *elcoor = NULL;
1180e2ec84fSDave May   PetscReal      *swarm_coor;
1190e2ec84fSDave May   PetscInt       *swarm_cellid;
1200e2ec84fSDave May   PetscReal      v1[2],v2[2],v3[2];
1210e2ec84fSDave May 
1220e2ec84fSDave May   PetscFunctionBegin;
1230e2ec84fSDave May   npoints_q = 1;
1240e2ec84fSDave May   for (d=0; d<nsub; d++) { npoints_q *= 4; }
1250e2ec84fSDave May   ierr = PetscMalloc1(dim*npoints_q,&xi);CHKERRQ(ierr);
1260e2ec84fSDave May 
1270e2ec84fSDave May   v1[0] = 0.0;  v1[1] = 0.0;
1280e2ec84fSDave May   v2[0] = 1.0;  v2[1] = 0.0;
1290e2ec84fSDave May   v3[0] = 0.0;  v3[1] = 1.0;
1300e2ec84fSDave May   depth = 0;
1310e2ec84fSDave May   pcnt = 0;
1320e2ec84fSDave May   ierr = subdivide_triangle(v1,v2,v3,depth,nsub,xi,&pcnt);CHKERRQ(ierr);
1330e2ec84fSDave May 
1340e2ec84fSDave May   npe = 3; /* nodes per element (triangle) */
1350e2ec84fSDave May   ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr);
1360e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
1370e2ec84fSDave May     ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr);
1380e2ec84fSDave May 
1390e2ec84fSDave May     basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1];
1400e2ec84fSDave May     basis[q][1] = xi[dim*q+0];
1410e2ec84fSDave May     basis[q][2] = xi[dim*q+1];
1420e2ec84fSDave May   }
1430e2ec84fSDave May 
1440e2ec84fSDave May   /* 0->cell, 1->edge, 2->vert */
1450e2ec84fSDave May   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
1460e2ec84fSDave May   nel = pe - ps;
1470e2ec84fSDave May 
1480e2ec84fSDave May   ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr);
1490e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
1500e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
1510e2ec84fSDave May 
1520e2ec84fSDave May   ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr);
1530e2ec84fSDave May   ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr);
1540e2ec84fSDave May 
1550e2ec84fSDave May   pcnt = 0;
1560e2ec84fSDave May   for (e=0; e<nel; e++) {
1570e2ec84fSDave May     ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
1580e2ec84fSDave May 
1590e2ec84fSDave May     for (q=0; q<npoints_q; q++) {
1600e2ec84fSDave May       for (d=0; d<dim; d++) {
1610e2ec84fSDave May         swarm_coor[dim*pcnt+d] = 0.0;
1620e2ec84fSDave May         for (k=0; k<npe; k++) {
1630ca99263SDave May           swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]);
1640e2ec84fSDave May         }
1650e2ec84fSDave May       }
1660e2ec84fSDave May       swarm_cellid[pcnt] = e;
1670e2ec84fSDave May       pcnt++;
1680e2ec84fSDave May     }
1690e2ec84fSDave May     ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
1700e2ec84fSDave May   }
1710e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
1720e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
1730e2ec84fSDave May 
1740e2ec84fSDave May   ierr = PetscFree(xi);CHKERRQ(ierr);
1750e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
1760e2ec84fSDave May     ierr = PetscFree(basis[q]);CHKERRQ(ierr);
1770e2ec84fSDave May   }
1780e2ec84fSDave May   ierr = PetscFree(basis);CHKERRQ(ierr);
1790e2ec84fSDave May   PetscFunctionReturn(0);
1800e2ec84fSDave May }
1810e2ec84fSDave May 
18278185223SDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm,DM dmc,PetscInt nsub)
18378185223SDave May {
18478185223SDave May   PetscErrorCode  ierr;
18578185223SDave May   PetscInt        dim,nfaces,nbasis;
18678185223SDave May   PetscInt        q,npoints_q,e,nel,pcnt,ps,pe,d,k,r;
187ef0bb6c7SMatthew G. Knepley   PetscTabulation T;
18878185223SDave May   Vec             coorlocal;
18978185223SDave May   PetscSection    coordSection;
19078185223SDave May   PetscScalar     *elcoor = NULL;
19178185223SDave May   PetscReal       *swarm_coor;
19278185223SDave May   PetscInt        *swarm_cellid;
19378185223SDave May   const PetscReal *xiq;
19478185223SDave May   PetscQuadrature quadrature;
19578185223SDave May   PetscFE         fe,feRef;
19678185223SDave May   PetscBool       is_simplex;
19778185223SDave May 
19878185223SDave May   PetscFunctionBegin;
19978185223SDave May   ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr);
20078185223SDave May   is_simplex = PETSC_FALSE;
20178185223SDave May   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
20278185223SDave May   ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr);
20378185223SDave May   if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
20478185223SDave May 
2052cf64d79SDave May   ierr = private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe);CHKERRQ(ierr);
20678185223SDave May 
20778185223SDave May   for (r=0; r<nsub; r++) {
20878185223SDave May     ierr = PetscFERefine(fe,&feRef);CHKERRQ(ierr);
2095dc5c000SMatthew G. Knepley     ierr = PetscFECopyQuadrature(feRef,fe);CHKERRQ(ierr);
21078185223SDave May     ierr = PetscFEDestroy(&feRef);CHKERRQ(ierr);
21178185223SDave May   }
21278185223SDave May 
21378185223SDave May   ierr = PetscFEGetQuadrature(fe,&quadrature);CHKERRQ(ierr);
21478185223SDave May   ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL);CHKERRQ(ierr);
21578185223SDave May   ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr);
216f9244615SMatthew G. Knepley   ierr = PetscFEGetCellTabulation(fe, 1, &T);CHKERRQ(ierr);
21778185223SDave May 
21878185223SDave May   /* 0->cell, 1->edge, 2->vert */
21978185223SDave May   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
22078185223SDave May   nel = pe - ps;
22178185223SDave May 
22278185223SDave May   ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr);
22378185223SDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
22478185223SDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
22578185223SDave May 
22678185223SDave May   ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr);
22778185223SDave May   ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr);
22878185223SDave May 
22978185223SDave May   pcnt = 0;
23078185223SDave May   for (e=0; e<nel; e++) {
23178185223SDave May     ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr);
23278185223SDave May 
23378185223SDave May     for (q=0; q<npoints_q; q++) {
23478185223SDave May       for (d=0; d<dim; d++) {
23578185223SDave May         swarm_coor[dim*pcnt+d] = 0.0;
23678185223SDave May         for (k=0; k<nbasis; k++) {
237ef0bb6c7SMatthew G. Knepley           swarm_coor[dim*pcnt+d] += T->T[0][q*nbasis + k] * PetscRealPart(elcoor[dim*k+d]);
23878185223SDave May         }
23978185223SDave May       }
24078185223SDave May       swarm_cellid[pcnt] = e;
24178185223SDave May       pcnt++;
24278185223SDave May     }
24392e40656SDave May     ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr);
24478185223SDave May   }
24578185223SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
24678185223SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
24778185223SDave May 
24878185223SDave May   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
24978185223SDave May   PetscFunctionReturn(0);
25078185223SDave May }
25178185223SDave May 
2520e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints)
2530e2ec84fSDave May {
2540e2ec84fSDave May   PetscErrorCode ierr;
255bc53056eSDave May   PetscInt       dim;
256bc53056eSDave May   PetscInt       ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,nfaces;
2570e2ec84fSDave May   PetscReal      *xi,ds,ds2;
2580e2ec84fSDave May   PetscReal      **basis;
2590e2ec84fSDave May   Vec            coorlocal;
2600e2ec84fSDave May   PetscSection   coordSection;
2610e2ec84fSDave May   PetscScalar    *elcoor = NULL;
2620e2ec84fSDave May   PetscReal      *swarm_coor;
2630e2ec84fSDave May   PetscInt       *swarm_cellid;
264bc53056eSDave May   PetscBool      is_simplex;
2650e2ec84fSDave May 
2660e2ec84fSDave May   PetscFunctionBegin;
267bc53056eSDave May   ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr);
268*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dim != 2,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only 2D is supported");
269bc53056eSDave May   is_simplex = PETSC_FALSE;
270bc53056eSDave May   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
271bc53056eSDave May   ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr);
272bc53056eSDave May   if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
273*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!is_simplex,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only the simplex is supported");
274bc53056eSDave May 
2750e2ec84fSDave May   ierr = PetscMalloc1(dim*npoints*npoints,&xi);CHKERRQ(ierr);
2760e2ec84fSDave May   pcnt = 0;
2770e2ec84fSDave May   ds = 1.0/((PetscReal)(npoints-1));
2780e2ec84fSDave May   ds2 = 1.0/((PetscReal)(npoints));
2790e2ec84fSDave May   for (jj = 0; jj<npoints; jj++) {
2800e2ec84fSDave May     for (ii=0; ii<npoints-jj; ii++) {
2810e2ec84fSDave May       xi[dim*pcnt+0] = ii * ds;
2820e2ec84fSDave May       xi[dim*pcnt+1] = jj * ds;
2830e2ec84fSDave May 
2840e2ec84fSDave May       xi[dim*pcnt+0] *= (1.0 - 1.2*ds2);
2850e2ec84fSDave May       xi[dim*pcnt+1] *= (1.0 - 1.2*ds2);
2860e2ec84fSDave May 
2870e2ec84fSDave May       xi[dim*pcnt+0] += 0.35*ds2;
2880e2ec84fSDave May       xi[dim*pcnt+1] += 0.35*ds2;
2890e2ec84fSDave May       pcnt++;
2900e2ec84fSDave May     }
2910e2ec84fSDave May   }
2920e2ec84fSDave May   npoints_q = pcnt;
2930e2ec84fSDave May 
2940e2ec84fSDave May   npe = 3; /* nodes per element (triangle) */
2950e2ec84fSDave May   ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr);
2960e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
2970e2ec84fSDave May     ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr);
2980e2ec84fSDave May 
2990e2ec84fSDave May     basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1];
3000e2ec84fSDave May     basis[q][1] = xi[dim*q+0];
3010e2ec84fSDave May     basis[q][2] = xi[dim*q+1];
3020e2ec84fSDave May   }
3030e2ec84fSDave May 
3040e2ec84fSDave May   /* 0->cell, 1->edge, 2->vert */
3050e2ec84fSDave May   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
3060e2ec84fSDave May   nel = pe - ps;
3070e2ec84fSDave May 
3080e2ec84fSDave May   ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr);
3090e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
3100e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
3110e2ec84fSDave May 
3120e2ec84fSDave May   ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr);
3130e2ec84fSDave May   ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr);
3140e2ec84fSDave May 
3150e2ec84fSDave May   pcnt = 0;
3160e2ec84fSDave May   for (e=0; e<nel; e++) {
3170e2ec84fSDave May     ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
3180e2ec84fSDave May 
3190e2ec84fSDave May     for (q=0; q<npoints_q; q++) {
3200e2ec84fSDave May       for (d=0; d<dim; d++) {
3210e2ec84fSDave May         swarm_coor[dim*pcnt+d] = 0.0;
3220e2ec84fSDave May         for (k=0; k<npe; k++) {
3230ca99263SDave May           swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]);
3240e2ec84fSDave May         }
3250e2ec84fSDave May       }
3260e2ec84fSDave May       swarm_cellid[pcnt] = e;
3270e2ec84fSDave May       pcnt++;
3280e2ec84fSDave May     }
3290e2ec84fSDave May     ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
3300e2ec84fSDave May   }
3310e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
3320e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
3330e2ec84fSDave May 
3340e2ec84fSDave May   ierr = PetscFree(xi);CHKERRQ(ierr);
3350e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
3360e2ec84fSDave May     ierr = PetscFree(basis[q]);CHKERRQ(ierr);
3370e2ec84fSDave May   }
3380e2ec84fSDave May   ierr = PetscFree(basis);CHKERRQ(ierr);
3390e2ec84fSDave May   PetscFunctionReturn(0);
3400e2ec84fSDave May }
3410e2ec84fSDave May 
3420e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param)
3430e2ec84fSDave May {
3440e2ec84fSDave May   PetscErrorCode ierr;
3450e2ec84fSDave May   PetscInt       dim;
3460e2ec84fSDave May 
3470e2ec84fSDave May   PetscFunctionBegin;
3480e2ec84fSDave May   ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr);
3490e2ec84fSDave May   switch (layout) {
3500e2ec84fSDave May     case DMSWARMPIC_LAYOUT_REGULAR:
351*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(dim == 3,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No 3D support for REGULAR+PLEX");
3520e2ec84fSDave May       ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm,celldm,layout_param);CHKERRQ(ierr);
3530e2ec84fSDave May       break;
3540e2ec84fSDave May     case DMSWARMPIC_LAYOUT_GAUSS:
355b12d2206SDave May     {
356b12d2206SDave May       PetscInt npoints,npoints1,ps,pe,nfaces;
357b12d2206SDave May       const PetscReal *xi;
358b12d2206SDave May       PetscBool is_simplex;
359b12d2206SDave May       PetscQuadrature quadrature;
360b12d2206SDave May 
361b12d2206SDave May       is_simplex = PETSC_FALSE;
362b12d2206SDave May       ierr = DMPlexGetHeightStratum(celldm,0,&ps,&pe);CHKERRQ(ierr);
363b12d2206SDave May       ierr = DMPlexGetConeSize(celldm, ps, &nfaces);CHKERRQ(ierr);
364b12d2206SDave May       if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
365b12d2206SDave May 
366b12d2206SDave May       npoints1 = layout_param;
367b12d2206SDave May       if (is_simplex) {
368e6a796c3SToby Isaac         ierr = PetscDTStroudConicalQuadrature(dim,1,npoints1,-1.0,1.0,&quadrature);CHKERRQ(ierr);
369b12d2206SDave May       } else {
370b12d2206SDave May         ierr = PetscDTGaussTensorQuadrature(dim,1,npoints1,-1.0,1.0,&quadrature);CHKERRQ(ierr);
371b12d2206SDave May       }
372b12d2206SDave May       ierr = PetscQuadratureGetData(quadrature,NULL,NULL,&npoints,&xi,NULL);CHKERRQ(ierr);
373b12d2206SDave May       ierr = private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm,celldm,npoints,(PetscReal*)xi);CHKERRQ(ierr);
374b12d2206SDave May       ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr);
375b12d2206SDave May     }
3760e2ec84fSDave May       break;
3770e2ec84fSDave May     case DMSWARMPIC_LAYOUT_SUBDIVISION:
37878185223SDave May       ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm,celldm,layout_param);CHKERRQ(ierr);
3790e2ec84fSDave May       break;
3800e2ec84fSDave May   }
3810e2ec84fSDave May   PetscFunctionReturn(0);
3820e2ec84fSDave May }
383e3cd5995SDave May 
384e3cd5995SDave May /*
385e3cd5995SDave May typedef struct {
386e3cd5995SDave May   PetscReal x,y;
387e3cd5995SDave May } Point2d;
388e3cd5995SDave May 
389e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s)
390e3cd5995SDave May {
391e3cd5995SDave May   PetscFunctionBegin;
392e3cd5995SDave May   *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y);
393e3cd5995SDave May   PetscFunctionReturn(0);
394e3cd5995SDave May }
395e3cd5995SDave May */
396e3cd5995SDave May /*
397e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v)
398e3cd5995SDave May {
399e3cd5995SDave May   PetscReal s1,s2,s3;
400e3cd5995SDave May   PetscBool b1, b2, b3;
401e3cd5995SDave May 
402e3cd5995SDave May   PetscFunctionBegin;
403e3cd5995SDave May   signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f;
404e3cd5995SDave May   signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f;
405e3cd5995SDave May   signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f;
406e3cd5995SDave May 
407e3cd5995SDave May   *v = ((b1 == b2) && (b2 == b3));
408e3cd5995SDave May   PetscFunctionReturn(0);
409e3cd5995SDave May }
410e3cd5995SDave May */
411e3cd5995SDave May /*
412e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ)
413e3cd5995SDave May {
414e3cd5995SDave May   PetscReal x1,y1,x2,y2,x3,y3;
415e3cd5995SDave May   PetscReal c,b[2],A[2][2],inv[2][2],detJ,od;
416e3cd5995SDave May 
417e3cd5995SDave May   PetscFunctionBegin;
418e3cd5995SDave May   x1 = coords[2*0+0];
419e3cd5995SDave May   x2 = coords[2*1+0];
420e3cd5995SDave May   x3 = coords[2*2+0];
421e3cd5995SDave May 
422e3cd5995SDave May   y1 = coords[2*0+1];
423e3cd5995SDave May   y2 = coords[2*1+1];
424e3cd5995SDave May   y3 = coords[2*2+1];
425e3cd5995SDave May 
426e3cd5995SDave May   c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3;
427e3cd5995SDave May   b[0] = xp[0] - c;
428e3cd5995SDave May   c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3;
429e3cd5995SDave May   b[1] = xp[1] - c;
430e3cd5995SDave May 
431e3cd5995SDave May   A[0][0] = -0.5*x1 + 0.5*x2;   A[0][1] = -0.5*x1 + 0.5*x3;
432e3cd5995SDave May   A[1][0] = -0.5*y1 + 0.5*y2;   A[1][1] = -0.5*y1 + 0.5*y3;
433e3cd5995SDave May 
434e3cd5995SDave May   detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0];
435e3cd5995SDave May   *dJ = PetscAbsReal(detJ);
436e3cd5995SDave May   od = 1.0/detJ;
437e3cd5995SDave May 
438e3cd5995SDave May   inv[0][0] =  A[1][1] * od;
439e3cd5995SDave May   inv[0][1] = -A[0][1] * od;
440e3cd5995SDave May   inv[1][0] = -A[1][0] * od;
441e3cd5995SDave May   inv[1][1] =  A[0][0] * od;
442e3cd5995SDave May 
443e3cd5995SDave May   xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1];
444e3cd5995SDave May   xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1];
445e3cd5995SDave May   PetscFunctionReturn(0);
446e3cd5995SDave May }
447e3cd5995SDave May */
448e3cd5995SDave May 
449a5f152d1SDave May static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscScalar coords[],PetscReal xip[],PetscReal *dJ)
450e3cd5995SDave May {
451e3cd5995SDave May   PetscReal x1,y1,x2,y2,x3,y3;
452e3cd5995SDave May   PetscReal b[2],A[2][2],inv[2][2],detJ,od;
453e3cd5995SDave May 
454e3cd5995SDave May   PetscFunctionBegin;
455a5f152d1SDave May   x1 = PetscRealPart(coords[2*0+0]);
456a5f152d1SDave May   x2 = PetscRealPart(coords[2*1+0]);
457a5f152d1SDave May   x3 = PetscRealPart(coords[2*2+0]);
458e3cd5995SDave May 
459a5f152d1SDave May   y1 = PetscRealPart(coords[2*0+1]);
460a5f152d1SDave May   y2 = PetscRealPart(coords[2*1+1]);
461a5f152d1SDave May   y3 = PetscRealPart(coords[2*2+1]);
462e3cd5995SDave May 
463e3cd5995SDave May   b[0] = xp[0] - x1;
464e3cd5995SDave May   b[1] = xp[1] - y1;
465e3cd5995SDave May 
466e3cd5995SDave May   A[0][0] = x2-x1;   A[0][1] = x3-x1;
467e3cd5995SDave May   A[1][0] = y2-y1;   A[1][1] = y3-y1;
468e3cd5995SDave May 
469e3cd5995SDave May   detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0];
470e3cd5995SDave May   *dJ = PetscAbsReal(detJ);
471e3cd5995SDave May   od = 1.0/detJ;
472e3cd5995SDave May 
473e3cd5995SDave May   inv[0][0] =  A[1][1] * od;
474e3cd5995SDave May   inv[0][1] = -A[0][1] * od;
475e3cd5995SDave May   inv[1][0] = -A[1][0] * od;
476e3cd5995SDave May   inv[1][1] =  A[0][0] * od;
477e3cd5995SDave May 
478e3cd5995SDave May   xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1];
479e3cd5995SDave May   xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1];
480e3cd5995SDave May   PetscFunctionReturn(0);
481e3cd5995SDave May }
482e3cd5995SDave May 
483e3cd5995SDave May PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field)
484e3cd5995SDave May {
485e3cd5995SDave May   PetscErrorCode  ierr;
486e3cd5995SDave May   const PetscReal PLEX_C_EPS = 1.0e-8;
487e3cd5995SDave May   Vec             v_field_l,denom_l,coor_l,denom;
488e3cd5995SDave May   PetscInt        k,p,e,npoints;
489e3cd5995SDave May   PetscInt        *mpfield_cell;
490e3cd5995SDave May   PetscReal       *mpfield_coor;
491a5f152d1SDave May   PetscReal       xi_p[2];
492a5f152d1SDave May   PetscScalar     Ni[3];
493e3cd5995SDave May   PetscSection    coordSection;
494e3cd5995SDave May   PetscScalar     *elcoor = NULL;
495e3cd5995SDave May 
496e3cd5995SDave May   PetscFunctionBegin;
497e3cd5995SDave May   ierr = VecZeroEntries(v_field);CHKERRQ(ierr);
498e3cd5995SDave May 
499e3cd5995SDave May   ierr = DMGetLocalVector(dm,&v_field_l);CHKERRQ(ierr);
500e3cd5995SDave May   ierr = DMGetGlobalVector(dm,&denom);CHKERRQ(ierr);
501e3cd5995SDave May   ierr = DMGetLocalVector(dm,&denom_l);CHKERRQ(ierr);
502e3cd5995SDave May   ierr = VecZeroEntries(v_field_l);CHKERRQ(ierr);
503e3cd5995SDave May   ierr = VecZeroEntries(denom);CHKERRQ(ierr);
504e3cd5995SDave May   ierr = VecZeroEntries(denom_l);CHKERRQ(ierr);
505e3cd5995SDave May 
506e3cd5995SDave May   ierr = DMGetCoordinatesLocal(dm,&coor_l);CHKERRQ(ierr);
507e3cd5995SDave May   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
508e3cd5995SDave May 
509e3cd5995SDave May   ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr);
510e3cd5995SDave May   ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr);
511e3cd5995SDave May   ierr = DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr);
512e3cd5995SDave May 
513e3cd5995SDave May   for (p=0; p<npoints; p++) {
514a5f152d1SDave May     PetscReal   *coor_p,dJ;
515a5f152d1SDave May     PetscScalar elfield[3];
516e3cd5995SDave May     PetscBool   point_located;
517e3cd5995SDave May 
518e3cd5995SDave May     e       = mpfield_cell[p];
519e3cd5995SDave May     coor_p  = &mpfield_coor[2*p];
520e3cd5995SDave May 
521e3cd5995SDave May     ierr = DMPlexVecGetClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr);
522e3cd5995SDave May 
523e3cd5995SDave May /*
524e3cd5995SDave May     while (!point_located && (failed_counter < 25)) {
525e3cd5995SDave May       ierr = PointInTriangle(point, coords[0], coords[1], coords[2], &point_located);CHKERRQ(ierr);
526e3cd5995SDave May       point.x = coor_p[0];
527e3cd5995SDave May       point.y = coor_p[1];
528e3cd5995SDave May       point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0);
529e3cd5995SDave May       point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0);
530e3cd5995SDave May       failed_counter++;
531e3cd5995SDave May     }
532e3cd5995SDave May 
533e3cd5995SDave May     if (!point_located) {
534e3cd5995SDave 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);
535e3cd5995SDave May     }
536e3cd5995SDave May 
537*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!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);
538e3cd5995SDave May     else {
539e3cd5995SDave May       ierr = _ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr);
540e3cd5995SDave May       xi_p[0] = 0.5*(xi_p[0] + 1.0);
541e3cd5995SDave May       xi_p[1] = 0.5*(xi_p[1] + 1.0);
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 */
547e3cd5995SDave May 
548e3cd5995SDave May     ierr = ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr);
549e3cd5995SDave May     /*
550e3cd5995SDave 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]);
551e3cd5995SDave May     */
552e3cd5995SDave May     /*
553e3cd5995SDave May      point_located = PETSC_TRUE;
554e3cd5995SDave May     if (xi_p[0] < 0.0) {
555e3cd5995SDave May       if (xi_p[0] > -PLEX_C_EPS) {
556e3cd5995SDave May         xi_p[0] = 0.0;
557e3cd5995SDave May       } else {
558e3cd5995SDave May         point_located = PETSC_FALSE;
559e3cd5995SDave May       }
560e3cd5995SDave May     }
561e3cd5995SDave May     if (xi_p[1] < 0.0) {
562e3cd5995SDave May       if (xi_p[1] > -PLEX_C_EPS) {
563e3cd5995SDave May         xi_p[1] = 0.0;
564e3cd5995SDave May       } else {
565e3cd5995SDave May         point_located = PETSC_FALSE;
566e3cd5995SDave May       }
567e3cd5995SDave May     }
568e3cd5995SDave May     if (xi_p[1] > (1.0-xi_p[0])) {
569e3cd5995SDave May       if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) {
570e3cd5995SDave May         xi_p[1] = 1.0 - xi_p[0];
571e3cd5995SDave May       } else {
572e3cd5995SDave May         point_located = PETSC_FALSE;
573e3cd5995SDave May       }
574e3cd5995SDave May     }
575e3cd5995SDave May     if (!point_located) {
576e3cd5995SDave May       PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]);
577e3cd5995SDave 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]);
578e3cd5995SDave May     }
579*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!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);
580e3cd5995SDave May     */
581e3cd5995SDave May 
582e3cd5995SDave May     Ni[0] = 1.0 - xi_p[0] - xi_p[1];
583e3cd5995SDave May     Ni[1] = xi_p[0];
584e3cd5995SDave May     Ni[2] = xi_p[1];
585e3cd5995SDave May 
586e3cd5995SDave May     point_located = PETSC_TRUE;
587e3cd5995SDave May     for (k=0; k<3; k++) {
588a5f152d1SDave May       if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE;
589a5f152d1SDave May       if (PetscRealPart(Ni[k]) > (1.0+PLEX_C_EPS)) point_located = PETSC_FALSE;
590e3cd5995SDave May     }
591e3cd5995SDave May     if (!point_located) {
592b5675b0fSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",(double)xi_p[0],(double)xi_p[1]);CHKERRQ(ierr);
593b5675b0fSBarry Smith       ierr = 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]));CHKERRQ(ierr);
594e3cd5995SDave May     }
595*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!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);
596e3cd5995SDave May 
597e3cd5995SDave May     for (k=0; k<3; k++) {
598e3cd5995SDave May       Ni[k] = Ni[k] * dJ;
599e3cd5995SDave May       elfield[k] = Ni[k] * swarm_field[p];
600e3cd5995SDave May     }
601e3cd5995SDave May     ierr = DMPlexVecRestoreClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr);
602e3cd5995SDave May 
603e3cd5995SDave May     ierr = DMPlexVecSetClosure(dm, NULL,v_field_l, e, elfield, ADD_VALUES);CHKERRQ(ierr);
604e3cd5995SDave May     ierr = DMPlexVecSetClosure(dm, NULL,denom_l, e, Ni, ADD_VALUES);CHKERRQ(ierr);
605e3cd5995SDave May   }
606e3cd5995SDave May 
607e3cd5995SDave May   ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr);
608e3cd5995SDave May   ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr);
609e3cd5995SDave May 
610e3cd5995SDave May   ierr = DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr);
611e3cd5995SDave May   ierr = DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr);
612e3cd5995SDave May   ierr = DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr);
613e3cd5995SDave May   ierr = DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr);
614e3cd5995SDave May 
615e3cd5995SDave May   ierr = VecPointwiseDivide(v_field,v_field,denom);CHKERRQ(ierr);
616e3cd5995SDave May 
617e3cd5995SDave May   ierr = DMRestoreLocalVector(dm,&v_field_l);CHKERRQ(ierr);
618e3cd5995SDave May   ierr = DMRestoreLocalVector(dm,&denom_l);CHKERRQ(ierr);
619e3cd5995SDave May   ierr = DMRestoreGlobalVector(dm,&denom);CHKERRQ(ierr);
620e3cd5995SDave May   PetscFunctionReturn(0);
621e3cd5995SDave May }
622e3cd5995SDave May 
62377048351SPatrick Sanan PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[])
624e3cd5995SDave May {
625e3cd5995SDave May   PetscErrorCode ierr;
626e3cd5995SDave May   PetscInt       f,dim;
627e3cd5995SDave May 
628e3cd5995SDave May   PetscFunctionBegin;
629e3cd5995SDave May   ierr = DMGetDimension(swarm,&dim);CHKERRQ(ierr);
630e3cd5995SDave May   switch (dim) {
631e3cd5995SDave May     case 2:
632e3cd5995SDave May       for (f=0; f<nfields; f++) {
633e3cd5995SDave May         PetscReal *swarm_field;
634e3cd5995SDave May 
63577048351SPatrick Sanan         ierr = DMSwarmDataFieldGetEntries(dfield[f],(void**)&swarm_field);CHKERRQ(ierr);
636e3cd5995SDave May         ierr = DMSwarmProjectField_ApproxP1_PLEX_2D(swarm,swarm_field,celldm,vecs[f]);CHKERRQ(ierr);
637e3cd5995SDave May       }
638e3cd5995SDave May       break;
639e3cd5995SDave May     case 3:
640e3cd5995SDave May       SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D");
641e3cd5995SDave May     default:
642e3cd5995SDave May       break;
643e3cd5995SDave May   }
644e3cd5995SDave May   PetscFunctionReturn(0);
645e3cd5995SDave May }
646431382f2SDave May 
647431382f2SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm,DM dmc,PetscInt npoints,PetscReal xi[])
648431382f2SDave May {
649431382f2SDave May   PetscBool       is_simplex,is_tensorcell;
650431382f2SDave May   PetscErrorCode  ierr;
65192e40656SDave May   PetscInt        dim,nfaces,ps,pe,p,d,nbasis,pcnt,e,k,nel;
65292e40656SDave May   PetscFE         fe;
65392e40656SDave May   PetscQuadrature quadrature;
654ef0bb6c7SMatthew G. Knepley   PetscTabulation T;
655ef0bb6c7SMatthew G. Knepley   PetscReal       *xiq;
65692e40656SDave May   Vec             coorlocal;
65792e40656SDave May   PetscSection    coordSection;
65892e40656SDave May   PetscScalar     *elcoor = NULL;
65992e40656SDave May   PetscReal       *swarm_coor;
66092e40656SDave May   PetscInt        *swarm_cellid;
661431382f2SDave May 
66292e40656SDave May   PetscFunctionBegin;
663431382f2SDave May   ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr);
664431382f2SDave May 
665431382f2SDave May   is_simplex = PETSC_FALSE;
666431382f2SDave May   is_tensorcell = PETSC_FALSE;
667431382f2SDave May   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
668431382f2SDave May   ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr);
669431382f2SDave May 
670431382f2SDave May   if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; }
671431382f2SDave May 
672431382f2SDave May   switch (dim) {
673431382f2SDave May     case 2:
674431382f2SDave May       if (nfaces == 4) { is_tensorcell = PETSC_TRUE; }
675431382f2SDave May       break;
676431382f2SDave May     case 3:
677431382f2SDave May       if (nfaces == 6) { is_tensorcell = PETSC_TRUE; }
678431382f2SDave May       break;
679431382f2SDave May     default:
680431382f2SDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for 2D, 3D");
681431382f2SDave May   }
682431382f2SDave May 
68392e40656SDave May   /* check points provided fail inside the reference cell */
684431382f2SDave May   if (is_simplex) {
68592e40656SDave May     for (p=0; p<npoints; p++) {
68692e40656SDave May       PetscReal sum;
68792e40656SDave May       for (d=0; d<dim; d++) {
688*2c71b3e2SJacob Faibussowitsch         PetscCheckFalse(xi[dim*p+d] < -1.0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain");
68992e40656SDave May       }
69092e40656SDave May       sum = 0.0;
69192e40656SDave May       for (d=0; d<dim; d++) {
69292e40656SDave May         sum += xi[dim*p+d];
69392e40656SDave May       }
694*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(sum > 0.0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain");
69592e40656SDave May     }
696431382f2SDave May   } else if (is_tensorcell) {
69792e40656SDave May     for (p=0; p<npoints; p++) {
69892e40656SDave May       for (d=0; d<dim; d++) {
699*2c71b3e2SJacob 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");
70092e40656SDave May       }
70192e40656SDave May     }
702431382f2SDave May   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for d-simplex and d-tensorcell");
703431382f2SDave May 
70492e40656SDave May   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)dm),&quadrature);CHKERRQ(ierr);
705b12d2206SDave May   ierr = PetscMalloc1(npoints*dim,&xiq);CHKERRQ(ierr);
706580bdb30SBarry Smith   ierr = PetscArraycpy(xiq,xi,npoints*dim);CHKERRQ(ierr);
707b12d2206SDave May   ierr = PetscQuadratureSetData(quadrature,dim,1,npoints,(const PetscReal*)xiq,NULL);CHKERRQ(ierr);
7082cf64d79SDave May   ierr = private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe);CHKERRQ(ierr);
70992e40656SDave May   ierr = PetscFESetQuadrature(fe,quadrature);CHKERRQ(ierr);
71092e40656SDave May   ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr);
711f9244615SMatthew G. Knepley   ierr = PetscFEGetCellTabulation(fe, 1, &T);CHKERRQ(ierr);
71292e40656SDave May 
71392e40656SDave May   /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */
71492e40656SDave May   /* 0->cell, 1->edge, 2->vert */
71592e40656SDave May   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
71692e40656SDave May   nel = pe - ps;
71792e40656SDave May 
71892e40656SDave May   ierr = DMSwarmSetLocalSizes(dm,npoints*nel,-1);CHKERRQ(ierr);
71992e40656SDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
72092e40656SDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
72192e40656SDave May 
72292e40656SDave May   ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr);
72392e40656SDave May   ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr);
72492e40656SDave May 
72592e40656SDave May   pcnt = 0;
72692e40656SDave May   for (e=0; e<nel; e++) {
72792e40656SDave May     ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr);
72892e40656SDave May 
72992e40656SDave May     for (p=0; p<npoints; p++) {
73092e40656SDave May       for (d=0; d<dim; d++) {
73192e40656SDave May         swarm_coor[dim*pcnt+d] = 0.0;
73292e40656SDave May         for (k=0; k<nbasis; k++) {
733ef0bb6c7SMatthew G. Knepley           swarm_coor[dim*pcnt+d] += T->T[0][p*nbasis + k] * PetscRealPart(elcoor[dim*k+d]);
73492e40656SDave May         }
73592e40656SDave May       }
73692e40656SDave May       swarm_cellid[pcnt] = e;
73792e40656SDave May       pcnt++;
73892e40656SDave May     }
73992e40656SDave May     ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr);
74092e40656SDave May   }
74192e40656SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
74292e40656SDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
74392e40656SDave May 
74492e40656SDave May   ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr);
74592e40656SDave May   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
746431382f2SDave May   PetscFunctionReturn(0);
747431382f2SDave May }
748