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 7b12d2206SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM,DM,PetscInt,PetscReal*xi); 8b12d2206SDave May 92cf64d79SDave May static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem) 102cf64d79SDave May { 112cf64d79SDave May const PetscInt Nc = 1; 122cf64d79SDave May PetscQuadrature q, fq; 132cf64d79SDave May DM K; 142cf64d79SDave May PetscSpace P; 152cf64d79SDave May PetscDualSpace Q; 162cf64d79SDave May PetscInt order, quadPointsPerEdge; 172cf64d79SDave May PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 182cf64d79SDave May PetscErrorCode ierr; 192cf64d79SDave May 202cf64d79SDave May PetscFunctionBegin; 212cf64d79SDave May /* Create space */ 222cf64d79SDave May ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P);CHKERRQ(ierr); 232cf64d79SDave May /*ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);*/ 242cf64d79SDave May ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr); 252cf64d79SDave May /*ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);*/ 262cf64d79SDave May ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 27d39dd5f5SToby Isaac ierr = PetscSpaceSetDegree(P,1,PETSC_DETERMINE);CHKERRQ(ierr); 282cf64d79SDave May ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr); 29157782e2SToby Isaac ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 302cf64d79SDave May ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 31cdcc6362SToby Isaac ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr); 322cf64d79SDave May ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr); 332cf64d79SDave May /* Create dual space */ 342cf64d79SDave May ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q);CHKERRQ(ierr); 352cf64d79SDave May ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 362cf64d79SDave May /*ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);*/ 372cf64d79SDave May ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 382cf64d79SDave May ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 392cf64d79SDave May ierr = DMDestroy(&K);CHKERRQ(ierr); 402cf64d79SDave May ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr); 412cf64d79SDave May ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr); 422cf64d79SDave May ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr); 432cf64d79SDave May /*ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);*/ 442cf64d79SDave May ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 452cf64d79SDave May ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 462cf64d79SDave May /* Create element */ 472cf64d79SDave May ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem);CHKERRQ(ierr); 482cf64d79SDave May /*ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);*/ 492cf64d79SDave May /*ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);*/ 502cf64d79SDave May ierr = PetscFESetType(*fem,PETSCFEBASIC);CHKERRQ(ierr); 512cf64d79SDave May ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr); 522cf64d79SDave May ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr); 532cf64d79SDave May ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr); 542cf64d79SDave May ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 552cf64d79SDave May ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 562cf64d79SDave May ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 572cf64d79SDave May /* Create quadrature (with specified order if given) */ 582cf64d79SDave May qorder = qorder >= 0 ? qorder : order; 592cf64d79SDave May quadPointsPerEdge = PetscMax(qorder + 1,1); 602cf64d79SDave May if (isSimplex) { 612cf64d79SDave May ierr = PetscDTGaussJacobiQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 622cf64d79SDave May ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 632cf64d79SDave May } 642cf64d79SDave May else { 652cf64d79SDave May ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 662cf64d79SDave May ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 672cf64d79SDave May } 682cf64d79SDave May ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr); 692cf64d79SDave May ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr); 702cf64d79SDave May ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 712cf64d79SDave May ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 722cf64d79SDave May PetscFunctionReturn(0); 732cf64d79SDave May } 742cf64d79SDave May 750e2ec84fSDave May PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[3],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np) 760e2ec84fSDave May { 770e2ec84fSDave May PetscReal v12[2],v23[2],v31[2]; 780e2ec84fSDave May PetscInt i; 790e2ec84fSDave May PetscErrorCode ierr; 800e2ec84fSDave May 810e2ec84fSDave May PetscFunctionBegin; 820e2ec84fSDave May if (depth == max) { 830e2ec84fSDave May PetscReal cx[2]; 840e2ec84fSDave May 850e2ec84fSDave May cx[0] = (v1[0] + v2[0] + v3[0])/3.0; 860e2ec84fSDave May cx[1] = (v1[1] + v2[1] + v3[1])/3.0; 870e2ec84fSDave May 880e2ec84fSDave May xi[2*(*np)+0] = cx[0]; 890e2ec84fSDave May xi[2*(*np)+1] = cx[1]; 900e2ec84fSDave May (*np)++; 910e2ec84fSDave May PetscFunctionReturn(0); 920e2ec84fSDave May } 930e2ec84fSDave May 940e2ec84fSDave May /* calculate midpoints of each side */ 950e2ec84fSDave May for (i = 0; i < 2; i++) { 960e2ec84fSDave May v12[i] = (v1[i]+v2[i])/2.0; 970e2ec84fSDave May v23[i] = (v2[i]+v3[i])/2.0; 980e2ec84fSDave May v31[i] = (v3[i]+v1[i])/2.0; 990e2ec84fSDave May } 1000e2ec84fSDave May 1010e2ec84fSDave May /* recursively subdivide new triangles */ 1020e2ec84fSDave May ierr = subdivide_triangle(v1,v12,v31,depth+1,max,xi,np);CHKERRQ(ierr); 1030e2ec84fSDave May ierr = subdivide_triangle(v2,v23,v12,depth+1,max,xi,np);CHKERRQ(ierr); 1040e2ec84fSDave May ierr = subdivide_triangle(v3,v31,v23,depth+1,max,xi,np);CHKERRQ(ierr); 1050e2ec84fSDave May ierr = subdivide_triangle(v12,v23,v31,depth+1,max,xi,np);CHKERRQ(ierr); 1060e2ec84fSDave May PetscFunctionReturn(0); 1070e2ec84fSDave May } 1080e2ec84fSDave May 1090e2ec84fSDave May 1100e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub) 1110e2ec84fSDave May { 1120e2ec84fSDave May PetscErrorCode ierr; 1130e2ec84fSDave May const PetscInt dim = 2; 1140e2ec84fSDave May PetscInt q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth; 1150e2ec84fSDave May PetscReal *xi; 1160e2ec84fSDave May PetscReal **basis; 1170e2ec84fSDave May Vec coorlocal; 1180e2ec84fSDave May PetscSection coordSection; 1190e2ec84fSDave May PetscScalar *elcoor = NULL; 1200e2ec84fSDave May PetscReal *swarm_coor; 1210e2ec84fSDave May PetscInt *swarm_cellid; 1220e2ec84fSDave May PetscReal v1[2],v2[2],v3[2]; 1230e2ec84fSDave May 1240e2ec84fSDave May PetscFunctionBegin; 1250e2ec84fSDave May 1260e2ec84fSDave May npoints_q = 1; 1270e2ec84fSDave May for (d=0; d<nsub; d++) { npoints_q *= 4; } 1280e2ec84fSDave May ierr = PetscMalloc1(dim*npoints_q,&xi);CHKERRQ(ierr); 1290e2ec84fSDave May 1300e2ec84fSDave May v1[0] = 0.0; v1[1] = 0.0; 1310e2ec84fSDave May v2[0] = 1.0; v2[1] = 0.0; 1320e2ec84fSDave May v3[0] = 0.0; v3[1] = 1.0; 1330e2ec84fSDave May depth = 0; 1340e2ec84fSDave May pcnt = 0; 1350e2ec84fSDave May ierr = subdivide_triangle(v1,v2,v3,depth,nsub,xi,&pcnt);CHKERRQ(ierr); 1360e2ec84fSDave May 1370e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 1380e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 1390e2ec84fSDave May for (q=0; q<npoints_q; q++) { 1400e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 1410e2ec84fSDave May 1420e2ec84fSDave May basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1]; 1430e2ec84fSDave May basis[q][1] = xi[dim*q+0]; 1440e2ec84fSDave May basis[q][2] = xi[dim*q+1]; 1450e2ec84fSDave May } 1460e2ec84fSDave May 1470e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 1480e2ec84fSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 1490e2ec84fSDave May nel = pe - ps; 1500e2ec84fSDave May 1510e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 1520e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1530e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1540e2ec84fSDave May 1550e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 1560e2ec84fSDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 1570e2ec84fSDave May 1580e2ec84fSDave May pcnt = 0; 1590e2ec84fSDave May for (e=0; e<nel; e++) { 1600e2ec84fSDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 1610e2ec84fSDave May 1620e2ec84fSDave May for (q=0; q<npoints_q; q++) { 1630e2ec84fSDave May for (d=0; d<dim; d++) { 1640e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 1650e2ec84fSDave May for (k=0; k<npe; k++) { 1660ca99263SDave May swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]); 1670e2ec84fSDave May } 1680e2ec84fSDave May } 1690e2ec84fSDave May swarm_cellid[pcnt] = e; 1700e2ec84fSDave May pcnt++; 1710e2ec84fSDave May } 1720e2ec84fSDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 1730e2ec84fSDave May } 1740e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1750e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1760e2ec84fSDave May 1770e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 1780e2ec84fSDave May for (q=0; q<npoints_q; q++) { 1790e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 1800e2ec84fSDave May } 1810e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 1820e2ec84fSDave May 1830e2ec84fSDave May PetscFunctionReturn(0); 1840e2ec84fSDave May } 1850e2ec84fSDave May 18678185223SDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm,DM dmc,PetscInt nsub) 18778185223SDave May { 18878185223SDave May PetscErrorCode ierr; 18978185223SDave May PetscInt dim,nfaces,nbasis; 19078185223SDave May PetscInt q,npoints_q,e,nel,pcnt,ps,pe,d,k,r; 191*ef0bb6c7SMatthew G. Knepley PetscTabulation T; 19278185223SDave May Vec coorlocal; 19378185223SDave May PetscSection coordSection; 19478185223SDave May PetscScalar *elcoor = NULL; 19578185223SDave May PetscReal *swarm_coor; 19678185223SDave May PetscInt *swarm_cellid; 19778185223SDave May const PetscReal *xiq; 19878185223SDave May PetscQuadrature quadrature; 19978185223SDave May PetscFE fe,feRef; 20078185223SDave May PetscBool is_simplex; 20178185223SDave May 20278185223SDave May PetscFunctionBegin; 20378185223SDave May 20478185223SDave May ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr); 20578185223SDave May 20678185223SDave May is_simplex = PETSC_FALSE; 20778185223SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 20878185223SDave May ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr); 20978185223SDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 21078185223SDave May 2112cf64d79SDave May ierr = private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe);CHKERRQ(ierr); 21278185223SDave May 21378185223SDave May for (r=0; r<nsub; r++) { 21478185223SDave May ierr = PetscFERefine(fe,&feRef);CHKERRQ(ierr); 2155dc5c000SMatthew G. Knepley ierr = PetscFECopyQuadrature(feRef,fe);CHKERRQ(ierr); 21678185223SDave May ierr = PetscFEDestroy(&feRef);CHKERRQ(ierr); 21778185223SDave May } 21878185223SDave May 21978185223SDave May ierr = PetscFEGetQuadrature(fe,&quadrature);CHKERRQ(ierr); 22078185223SDave May ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL);CHKERRQ(ierr); 22178185223SDave May ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr); 222*ef0bb6c7SMatthew G. Knepley ierr = PetscFEGetCellTabulation(fe, &T);CHKERRQ(ierr); 22378185223SDave May 22478185223SDave May /* 0->cell, 1->edge, 2->vert */ 22578185223SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 22678185223SDave May nel = pe - ps; 22778185223SDave May 22878185223SDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 22978185223SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 23078185223SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 23178185223SDave May 23278185223SDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 23378185223SDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 23478185223SDave May 23578185223SDave May pcnt = 0; 23678185223SDave May for (e=0; e<nel; e++) { 23778185223SDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 23878185223SDave May 23978185223SDave May for (q=0; q<npoints_q; q++) { 24078185223SDave May for (d=0; d<dim; d++) { 24178185223SDave May swarm_coor[dim*pcnt+d] = 0.0; 24278185223SDave May for (k=0; k<nbasis; k++) { 243*ef0bb6c7SMatthew G. Knepley swarm_coor[dim*pcnt+d] += T->T[0][q*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); 24478185223SDave May } 24578185223SDave May } 24678185223SDave May swarm_cellid[pcnt] = e; 24778185223SDave May pcnt++; 24878185223SDave May } 24992e40656SDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 25078185223SDave May } 25178185223SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 25278185223SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 25378185223SDave May 25478185223SDave May ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 25578185223SDave May PetscFunctionReturn(0); 25678185223SDave May } 25778185223SDave May 2580e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints) 2590e2ec84fSDave May { 2600e2ec84fSDave May PetscErrorCode ierr; 261bc53056eSDave May PetscInt dim; 262bc53056eSDave May PetscInt ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,nfaces; 2630e2ec84fSDave May PetscReal *xi,ds,ds2; 2640e2ec84fSDave May PetscReal **basis; 2650e2ec84fSDave May Vec coorlocal; 2660e2ec84fSDave May PetscSection coordSection; 2670e2ec84fSDave May PetscScalar *elcoor = NULL; 2680e2ec84fSDave May PetscReal *swarm_coor; 2690e2ec84fSDave May PetscInt *swarm_cellid; 270bc53056eSDave May PetscBool is_simplex; 2710e2ec84fSDave May 2720e2ec84fSDave May PetscFunctionBegin; 273bc53056eSDave May ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr); 274bc53056eSDave May if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only 2D is supported"); 275bc53056eSDave May is_simplex = PETSC_FALSE; 276bc53056eSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 277bc53056eSDave May ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr); 278bc53056eSDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 279bc53056eSDave May if (!is_simplex) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only the simplex is supported"); 280bc53056eSDave May 2810e2ec84fSDave May ierr = PetscMalloc1(dim*npoints*npoints,&xi);CHKERRQ(ierr); 2820e2ec84fSDave May pcnt = 0; 2830e2ec84fSDave May ds = 1.0/((PetscReal)(npoints-1)); 2840e2ec84fSDave May ds2 = 1.0/((PetscReal)(npoints)); 2850e2ec84fSDave May for (jj = 0; jj<npoints; jj++) { 2860e2ec84fSDave May for (ii=0; ii<npoints-jj; ii++) { 2870e2ec84fSDave May xi[dim*pcnt+0] = ii * ds; 2880e2ec84fSDave May xi[dim*pcnt+1] = jj * ds; 2890e2ec84fSDave May 2900e2ec84fSDave May xi[dim*pcnt+0] *= (1.0 - 1.2*ds2); 2910e2ec84fSDave May xi[dim*pcnt+1] *= (1.0 - 1.2*ds2); 2920e2ec84fSDave May 2930e2ec84fSDave May xi[dim*pcnt+0] += 0.35*ds2; 2940e2ec84fSDave May xi[dim*pcnt+1] += 0.35*ds2; 2950e2ec84fSDave May 2960e2ec84fSDave May pcnt++; 2970e2ec84fSDave May } 2980e2ec84fSDave May } 2990e2ec84fSDave May npoints_q = pcnt; 3000e2ec84fSDave May 3010e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 3020e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 3030e2ec84fSDave May for (q=0; q<npoints_q; q++) { 3040e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 3050e2ec84fSDave May 3060e2ec84fSDave May basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1]; 3070e2ec84fSDave May basis[q][1] = xi[dim*q+0]; 3080e2ec84fSDave May basis[q][2] = xi[dim*q+1]; 3090e2ec84fSDave May } 3100e2ec84fSDave May 3110e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 3120e2ec84fSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 3130e2ec84fSDave May nel = pe - ps; 3140e2ec84fSDave May 3150e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 3160e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3170e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3180e2ec84fSDave May 3190e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 3200e2ec84fSDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 3210e2ec84fSDave May 3220e2ec84fSDave May pcnt = 0; 3230e2ec84fSDave May for (e=0; e<nel; e++) { 3240e2ec84fSDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 3250e2ec84fSDave May 3260e2ec84fSDave May for (q=0; q<npoints_q; q++) { 3270e2ec84fSDave May for (d=0; d<dim; d++) { 3280e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 3290e2ec84fSDave May for (k=0; k<npe; k++) { 3300ca99263SDave May swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]); 3310e2ec84fSDave May } 3320e2ec84fSDave May } 3330e2ec84fSDave May swarm_cellid[pcnt] = e; 3340e2ec84fSDave May pcnt++; 3350e2ec84fSDave May } 3360e2ec84fSDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 3370e2ec84fSDave May } 3380e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3390e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3400e2ec84fSDave May 3410e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 3420e2ec84fSDave May for (q=0; q<npoints_q; q++) { 3430e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 3440e2ec84fSDave May } 3450e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 3460e2ec84fSDave May 3470e2ec84fSDave May PetscFunctionReturn(0); 3480e2ec84fSDave May } 3490e2ec84fSDave May 3500e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param) 3510e2ec84fSDave May { 3520e2ec84fSDave May PetscErrorCode ierr; 3530e2ec84fSDave May PetscInt dim; 3540e2ec84fSDave May 3550e2ec84fSDave May PetscFunctionBegin; 3560e2ec84fSDave May ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr); 3570e2ec84fSDave May switch (layout) { 3580e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 35978185223SDave May if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No 3D support for REGULAR+PLEX"); 3600e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm,celldm,layout_param);CHKERRQ(ierr); 3610e2ec84fSDave May break; 3620e2ec84fSDave May case DMSWARMPIC_LAYOUT_GAUSS: 363b12d2206SDave May { 364b12d2206SDave May PetscInt npoints,npoints1,ps,pe,nfaces; 365b12d2206SDave May const PetscReal *xi; 366b12d2206SDave May PetscBool is_simplex; 367b12d2206SDave May PetscQuadrature quadrature; 368b12d2206SDave May 369b12d2206SDave May is_simplex = PETSC_FALSE; 370b12d2206SDave May ierr = DMPlexGetHeightStratum(celldm,0,&ps,&pe);CHKERRQ(ierr); 371b12d2206SDave May ierr = DMPlexGetConeSize(celldm, ps, &nfaces);CHKERRQ(ierr); 372b12d2206SDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 373b12d2206SDave May 374b12d2206SDave May npoints1 = layout_param; 375b12d2206SDave May if (is_simplex) { 376b12d2206SDave May ierr = PetscDTGaussJacobiQuadrature(dim,1,npoints1,-1.0,1.0,&quadrature);CHKERRQ(ierr); 377b12d2206SDave May } else { 378b12d2206SDave May ierr = PetscDTGaussTensorQuadrature(dim,1,npoints1,-1.0,1.0,&quadrature);CHKERRQ(ierr); 379b12d2206SDave May } 380b12d2206SDave May ierr = PetscQuadratureGetData(quadrature,NULL,NULL,&npoints,&xi,NULL);CHKERRQ(ierr); 381b12d2206SDave May ierr = private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm,celldm,npoints,(PetscReal*)xi);CHKERRQ(ierr); 382b12d2206SDave May ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr); 383b12d2206SDave May } 3840e2ec84fSDave May break; 3850e2ec84fSDave May case DMSWARMPIC_LAYOUT_SUBDIVISION: 38678185223SDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm,celldm,layout_param);CHKERRQ(ierr); 3870e2ec84fSDave May break; 3880e2ec84fSDave May } 3890e2ec84fSDave May PetscFunctionReturn(0); 3900e2ec84fSDave May } 391e3cd5995SDave May 392e3cd5995SDave May /* 393e3cd5995SDave May typedef struct { 394e3cd5995SDave May PetscReal x,y; 395e3cd5995SDave May } Point2d; 396e3cd5995SDave May 397e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s) 398e3cd5995SDave May { 399e3cd5995SDave May PetscFunctionBegin; 400e3cd5995SDave May *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y); 401e3cd5995SDave May PetscFunctionReturn(0); 402e3cd5995SDave May } 403e3cd5995SDave May */ 404e3cd5995SDave May /* 405e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v) 406e3cd5995SDave May { 407e3cd5995SDave May PetscReal s1,s2,s3; 408e3cd5995SDave May PetscBool b1, b2, b3; 409e3cd5995SDave May 410e3cd5995SDave May PetscFunctionBegin; 411e3cd5995SDave May signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f; 412e3cd5995SDave May signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f; 413e3cd5995SDave May signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f; 414e3cd5995SDave May 415e3cd5995SDave May *v = ((b1 == b2) && (b2 == b3)); 416e3cd5995SDave May PetscFunctionReturn(0); 417e3cd5995SDave May } 418e3cd5995SDave May */ 419e3cd5995SDave May /* 420e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ) 421e3cd5995SDave May { 422e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 423e3cd5995SDave May PetscReal c,b[2],A[2][2],inv[2][2],detJ,od; 424e3cd5995SDave May 425e3cd5995SDave May PetscFunctionBegin; 426e3cd5995SDave May x1 = coords[2*0+0]; 427e3cd5995SDave May x2 = coords[2*1+0]; 428e3cd5995SDave May x3 = coords[2*2+0]; 429e3cd5995SDave May 430e3cd5995SDave May y1 = coords[2*0+1]; 431e3cd5995SDave May y2 = coords[2*1+1]; 432e3cd5995SDave May y3 = coords[2*2+1]; 433e3cd5995SDave May 434e3cd5995SDave May c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3; 435e3cd5995SDave May b[0] = xp[0] - c; 436e3cd5995SDave May c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3; 437e3cd5995SDave May b[1] = xp[1] - c; 438e3cd5995SDave May 439e3cd5995SDave May A[0][0] = -0.5*x1 + 0.5*x2; A[0][1] = -0.5*x1 + 0.5*x3; 440e3cd5995SDave May A[1][0] = -0.5*y1 + 0.5*y2; A[1][1] = -0.5*y1 + 0.5*y3; 441e3cd5995SDave May 442e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 443e3cd5995SDave May *dJ = PetscAbsReal(detJ); 444e3cd5995SDave May od = 1.0/detJ; 445e3cd5995SDave May 446e3cd5995SDave May inv[0][0] = A[1][1] * od; 447e3cd5995SDave May inv[0][1] = -A[0][1] * od; 448e3cd5995SDave May inv[1][0] = -A[1][0] * od; 449e3cd5995SDave May inv[1][1] = A[0][0] * od; 450e3cd5995SDave May 451e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 452e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 453e3cd5995SDave May PetscFunctionReturn(0); 454e3cd5995SDave May } 455e3cd5995SDave May */ 456e3cd5995SDave May 457a5f152d1SDave May static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscScalar coords[],PetscReal xip[],PetscReal *dJ) 458e3cd5995SDave May { 459e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 460e3cd5995SDave May PetscReal b[2],A[2][2],inv[2][2],detJ,od; 461e3cd5995SDave May 462e3cd5995SDave May PetscFunctionBegin; 463a5f152d1SDave May x1 = PetscRealPart(coords[2*0+0]); 464a5f152d1SDave May x2 = PetscRealPart(coords[2*1+0]); 465a5f152d1SDave May x3 = PetscRealPart(coords[2*2+0]); 466e3cd5995SDave May 467a5f152d1SDave May y1 = PetscRealPart(coords[2*0+1]); 468a5f152d1SDave May y2 = PetscRealPart(coords[2*1+1]); 469a5f152d1SDave May y3 = PetscRealPart(coords[2*2+1]); 470e3cd5995SDave May 471e3cd5995SDave May b[0] = xp[0] - x1; 472e3cd5995SDave May b[1] = xp[1] - y1; 473e3cd5995SDave May 474e3cd5995SDave May A[0][0] = x2-x1; A[0][1] = x3-x1; 475e3cd5995SDave May A[1][0] = y2-y1; A[1][1] = y3-y1; 476e3cd5995SDave May 477e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 478e3cd5995SDave May *dJ = PetscAbsReal(detJ); 479e3cd5995SDave May od = 1.0/detJ; 480e3cd5995SDave May 481e3cd5995SDave May inv[0][0] = A[1][1] * od; 482e3cd5995SDave May inv[0][1] = -A[0][1] * od; 483e3cd5995SDave May inv[1][0] = -A[1][0] * od; 484e3cd5995SDave May inv[1][1] = A[0][0] * od; 485e3cd5995SDave May 486e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 487e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 488e3cd5995SDave May PetscFunctionReturn(0); 489e3cd5995SDave May } 490e3cd5995SDave May 491e3cd5995SDave May PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field) 492e3cd5995SDave May { 493e3cd5995SDave May PetscErrorCode ierr; 494e3cd5995SDave May const PetscReal PLEX_C_EPS = 1.0e-8; 495e3cd5995SDave May Vec v_field_l,denom_l,coor_l,denom; 496e3cd5995SDave May PetscInt k,p,e,npoints; 497e3cd5995SDave May PetscInt *mpfield_cell; 498e3cd5995SDave May PetscReal *mpfield_coor; 499a5f152d1SDave May PetscReal xi_p[2]; 500a5f152d1SDave May PetscScalar Ni[3]; 501e3cd5995SDave May PetscSection coordSection; 502e3cd5995SDave May PetscScalar *elcoor = NULL; 503e3cd5995SDave May 504e3cd5995SDave May PetscFunctionBegin; 505e3cd5995SDave May ierr = VecZeroEntries(v_field);CHKERRQ(ierr); 506e3cd5995SDave May 507e3cd5995SDave May ierr = DMGetLocalVector(dm,&v_field_l);CHKERRQ(ierr); 508e3cd5995SDave May ierr = DMGetGlobalVector(dm,&denom);CHKERRQ(ierr); 509e3cd5995SDave May ierr = DMGetLocalVector(dm,&denom_l);CHKERRQ(ierr); 510e3cd5995SDave May ierr = VecZeroEntries(v_field_l);CHKERRQ(ierr); 511e3cd5995SDave May ierr = VecZeroEntries(denom);CHKERRQ(ierr); 512e3cd5995SDave May ierr = VecZeroEntries(denom_l);CHKERRQ(ierr); 513e3cd5995SDave May 514e3cd5995SDave May ierr = DMGetCoordinatesLocal(dm,&coor_l);CHKERRQ(ierr); 515e3cd5995SDave May ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 516e3cd5995SDave May 517e3cd5995SDave May ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr); 518e3cd5995SDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 519e3cd5995SDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 520e3cd5995SDave May 521e3cd5995SDave May for (p=0; p<npoints; p++) { 522a5f152d1SDave May PetscReal *coor_p,dJ; 523a5f152d1SDave May PetscScalar elfield[3]; 524e3cd5995SDave May PetscBool point_located; 525e3cd5995SDave May 526e3cd5995SDave May e = mpfield_cell[p]; 527e3cd5995SDave May coor_p = &mpfield_coor[2*p]; 528e3cd5995SDave May 529e3cd5995SDave May ierr = DMPlexVecGetClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr); 530e3cd5995SDave May 531e3cd5995SDave May /* 532e3cd5995SDave May while (!point_located && (failed_counter < 25)) { 533e3cd5995SDave May ierr = PointInTriangle(point, coords[0], coords[1], coords[2], &point_located);CHKERRQ(ierr); 534e3cd5995SDave May point.x = coor_p[0]; 535e3cd5995SDave May point.y = coor_p[1]; 536e3cd5995SDave May point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 537e3cd5995SDave May point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 538e3cd5995SDave May failed_counter++; 539e3cd5995SDave May } 540e3cd5995SDave May 541e3cd5995SDave May if (!point_located){ 542e3cd5995SDave 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); 543e3cd5995SDave May } 544e3cd5995SDave May 545e3cd5995SDave May if (!point_located) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)\n",point.x,point.y,e); 546e3cd5995SDave May else { 547e3cd5995SDave May ierr = _ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr); 548e3cd5995SDave May xi_p[0] = 0.5*(xi_p[0] + 1.0); 549e3cd5995SDave May xi_p[1] = 0.5*(xi_p[1] + 1.0); 550e3cd5995SDave May 551e3cd5995SDave 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]); 552e3cd5995SDave May 553e3cd5995SDave May } 554e3cd5995SDave May */ 555e3cd5995SDave May 556e3cd5995SDave May ierr = ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr); 557e3cd5995SDave May /* 558e3cd5995SDave 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]); 559e3cd5995SDave May */ 560e3cd5995SDave May /* 561e3cd5995SDave May point_located = PETSC_TRUE; 562e3cd5995SDave May if (xi_p[0] < 0.0) { 563e3cd5995SDave May if (xi_p[0] > -PLEX_C_EPS) { 564e3cd5995SDave May xi_p[0] = 0.0; 565e3cd5995SDave May } else { 566e3cd5995SDave May point_located = PETSC_FALSE; 567e3cd5995SDave May } 568e3cd5995SDave May } 569e3cd5995SDave May if (xi_p[1] < 0.0) { 570e3cd5995SDave May if (xi_p[1] > -PLEX_C_EPS) { 571e3cd5995SDave May xi_p[1] = 0.0; 572e3cd5995SDave May } else { 573e3cd5995SDave May point_located = PETSC_FALSE; 574e3cd5995SDave May } 575e3cd5995SDave May } 576e3cd5995SDave May if (xi_p[1] > (1.0-xi_p[0])) { 577e3cd5995SDave May if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { 578e3cd5995SDave May xi_p[1] = 1.0 - xi_p[0]; 579e3cd5995SDave May } else { 580e3cd5995SDave May point_located = PETSC_FALSE; 581e3cd5995SDave May } 582e3cd5995SDave May } 583e3cd5995SDave May if (!point_located){ 584e3cd5995SDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); 585e3cd5995SDave 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]); 586e3cd5995SDave May } 587e3cd5995SDave May if (!point_located) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)\n",coor_p[0],coor_p[1],e); 588e3cd5995SDave May */ 589e3cd5995SDave May 590e3cd5995SDave May Ni[0] = 1.0 - xi_p[0] - xi_p[1]; 591e3cd5995SDave May Ni[1] = xi_p[0]; 592e3cd5995SDave May Ni[2] = xi_p[1]; 593e3cd5995SDave May 594e3cd5995SDave May point_located = PETSC_TRUE; 595e3cd5995SDave May for (k=0; k<3; k++) { 596a5f152d1SDave May if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; 597a5f152d1SDave May if (PetscRealPart(Ni[k]) > (1.0+PLEX_C_EPS)) point_located = PETSC_FALSE; 598e3cd5995SDave May } 599e3cd5995SDave May if (!point_located){ 600d81390bbSDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",(double)xi_p[0],(double)xi_p[1]); 601d81390bbSDave 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",(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])); 602e3cd5995SDave May } 603d81390bbSDave May if (!point_located) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)\n",(double)coor_p[0],(double)coor_p[1],e); 604e3cd5995SDave May 605e3cd5995SDave May for (k=0; k<3; k++) { 606e3cd5995SDave May Ni[k] = Ni[k] * dJ; 607e3cd5995SDave May elfield[k] = Ni[k] * swarm_field[p]; 608e3cd5995SDave May } 609e3cd5995SDave May 610e3cd5995SDave May ierr = DMPlexVecRestoreClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr); 611e3cd5995SDave May 612e3cd5995SDave May ierr = DMPlexVecSetClosure(dm, NULL,v_field_l, e, elfield, ADD_VALUES);CHKERRQ(ierr); 613e3cd5995SDave May ierr = DMPlexVecSetClosure(dm, NULL,denom_l, e, Ni, ADD_VALUES);CHKERRQ(ierr); 614e3cd5995SDave May } 615e3cd5995SDave May 616e3cd5995SDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 617e3cd5995SDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 618e3cd5995SDave May 619e3cd5995SDave May ierr = DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 620e3cd5995SDave May ierr = DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 621e3cd5995SDave May ierr = DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 622e3cd5995SDave May ierr = DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 623e3cd5995SDave May 624e3cd5995SDave May ierr = VecPointwiseDivide(v_field,v_field,denom);CHKERRQ(ierr); 625e3cd5995SDave May 626e3cd5995SDave May ierr = DMRestoreLocalVector(dm,&v_field_l);CHKERRQ(ierr); 627e3cd5995SDave May ierr = DMRestoreLocalVector(dm,&denom_l);CHKERRQ(ierr); 628e3cd5995SDave May ierr = DMRestoreGlobalVector(dm,&denom);CHKERRQ(ierr); 629e3cd5995SDave May 630e3cd5995SDave May PetscFunctionReturn(0); 631e3cd5995SDave May } 632e3cd5995SDave May 63377048351SPatrick Sanan PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[]) 634e3cd5995SDave May { 635e3cd5995SDave May PetscErrorCode ierr; 636e3cd5995SDave May PetscInt f,dim; 637e3cd5995SDave May 638e3cd5995SDave May PetscFunctionBegin; 639e3cd5995SDave May ierr = DMGetDimension(swarm,&dim);CHKERRQ(ierr); 640e3cd5995SDave May switch (dim) { 641e3cd5995SDave May case 2: 642e3cd5995SDave May for (f=0; f<nfields; f++) { 643e3cd5995SDave May PetscReal *swarm_field; 644e3cd5995SDave May 64577048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(dfield[f],(void**)&swarm_field);CHKERRQ(ierr); 646e3cd5995SDave May ierr = DMSwarmProjectField_ApproxP1_PLEX_2D(swarm,swarm_field,celldm,vecs[f]);CHKERRQ(ierr); 647e3cd5995SDave May } 648e3cd5995SDave May break; 649e3cd5995SDave May case 3: 650e3cd5995SDave May SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D"); 651e3cd5995SDave May break; 652e3cd5995SDave May default: 653e3cd5995SDave May break; 654e3cd5995SDave May } 655e3cd5995SDave May 656e3cd5995SDave May PetscFunctionReturn(0); 657e3cd5995SDave May } 658431382f2SDave May 659431382f2SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm,DM dmc,PetscInt npoints,PetscReal xi[]) 660431382f2SDave May { 661431382f2SDave May PetscBool is_simplex,is_tensorcell; 662431382f2SDave May PetscErrorCode ierr; 66392e40656SDave May PetscInt dim,nfaces,ps,pe,p,d,nbasis,pcnt,e,k,nel; 66492e40656SDave May PetscFE fe; 66592e40656SDave May PetscQuadrature quadrature; 666*ef0bb6c7SMatthew G. Knepley PetscTabulation T; 667*ef0bb6c7SMatthew G. Knepley PetscReal *xiq; 66892e40656SDave May Vec coorlocal; 66992e40656SDave May PetscSection coordSection; 67092e40656SDave May PetscScalar *elcoor = NULL; 67192e40656SDave May PetscReal *swarm_coor; 67292e40656SDave May PetscInt *swarm_cellid; 673431382f2SDave May 67492e40656SDave May PetscFunctionBegin; 675431382f2SDave May ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr); 676431382f2SDave May 677431382f2SDave May is_simplex = PETSC_FALSE; 678431382f2SDave May is_tensorcell = PETSC_FALSE; 679431382f2SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 680431382f2SDave May ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr); 681431382f2SDave May 682431382f2SDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 683431382f2SDave May 684431382f2SDave May switch (dim) { 685431382f2SDave May case 2: 686431382f2SDave May if (nfaces == 4) { is_tensorcell = PETSC_TRUE; } 687431382f2SDave May break; 688431382f2SDave May case 3: 689431382f2SDave May if (nfaces == 6) { is_tensorcell = PETSC_TRUE; } 690431382f2SDave May break; 691431382f2SDave May default: 692431382f2SDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for 2D, 3D"); 693431382f2SDave May break; 694431382f2SDave May } 695431382f2SDave May 69692e40656SDave May /* check points provided fail inside the reference cell */ 697431382f2SDave May if (is_simplex) { 69892e40656SDave May for (p=0; p<npoints; p++) { 69992e40656SDave May PetscReal sum; 70092e40656SDave May for (d=0; d<dim; d++) { 70192e40656SDave May if (xi[dim*p+d] < -1.0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain"); 70292e40656SDave May } 70392e40656SDave May sum = 0.0; 70492e40656SDave May for (d=0; d<dim; d++) { 70592e40656SDave May sum += xi[dim*p+d]; 70692e40656SDave May } 70792e40656SDave May if (sum > 0.0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain"); 70892e40656SDave May } 709431382f2SDave May } else if (is_tensorcell) { 71092e40656SDave May for (p=0; p<npoints; p++) { 71192e40656SDave May for (d=0; d<dim; d++) { 71292e40656SDave May if (PetscAbsReal(xi[dim*p+d]) > 1.0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the tensor domain [-1,1]^d"); 71392e40656SDave May } 71492e40656SDave May } 715431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for d-simplex and d-tensorcell"); 716431382f2SDave May 71792e40656SDave May ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)dm),&quadrature);CHKERRQ(ierr); 718b12d2206SDave May ierr = PetscMalloc1(npoints*dim,&xiq);CHKERRQ(ierr); 719580bdb30SBarry Smith ierr = PetscArraycpy(xiq,xi,npoints*dim);CHKERRQ(ierr); 720b12d2206SDave May ierr = PetscQuadratureSetData(quadrature,dim,1,npoints,(const PetscReal*)xiq,NULL);CHKERRQ(ierr); 7212cf64d79SDave May ierr = private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe);CHKERRQ(ierr); 72292e40656SDave May ierr = PetscFESetQuadrature(fe,quadrature);CHKERRQ(ierr); 72392e40656SDave May ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr); 724*ef0bb6c7SMatthew G. Knepley ierr = PetscFEGetCellTabulation(fe, &T);CHKERRQ(ierr); 72592e40656SDave May 72692e40656SDave May /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ 72792e40656SDave May /* 0->cell, 1->edge, 2->vert */ 72892e40656SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 72992e40656SDave May nel = pe - ps; 73092e40656SDave May 73192e40656SDave May ierr = DMSwarmSetLocalSizes(dm,npoints*nel,-1);CHKERRQ(ierr); 73292e40656SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 73392e40656SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 73492e40656SDave May 73592e40656SDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 73692e40656SDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 73792e40656SDave May 73892e40656SDave May pcnt = 0; 73992e40656SDave May for (e=0; e<nel; e++) { 74092e40656SDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 74192e40656SDave May 74292e40656SDave May for (p=0; p<npoints; p++) { 74392e40656SDave May for (d=0; d<dim; d++) { 74492e40656SDave May swarm_coor[dim*pcnt+d] = 0.0; 74592e40656SDave May for (k=0; k<nbasis; k++) { 746*ef0bb6c7SMatthew G. Knepley swarm_coor[dim*pcnt+d] += T->T[0][p*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); 74792e40656SDave May } 74892e40656SDave May } 74992e40656SDave May swarm_cellid[pcnt] = e; 75092e40656SDave May pcnt++; 75192e40656SDave May } 75292e40656SDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 75392e40656SDave May } 75492e40656SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 75592e40656SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 75692e40656SDave May 75792e40656SDave May ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr); 75892e40656SDave May ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 75992e40656SDave May 760431382f2SDave May PetscFunctionReturn(0); 761431382f2SDave May } 762