10e2ec84fSDave May #include <petscdm.h> 20e2ec84fSDave May #include <petscdmplex.h> 30e2ec84fSDave May #include <petscdmswarm.h> 4e3cd5995SDave May #include "data_bucket.h" 50e2ec84fSDave May 60e2ec84fSDave May PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[3],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np) 70e2ec84fSDave May { 80e2ec84fSDave May PetscReal v12[2],v23[2],v31[2]; 90e2ec84fSDave May PetscInt i; 100e2ec84fSDave May PetscErrorCode ierr; 110e2ec84fSDave May 120e2ec84fSDave May PetscFunctionBegin; 130e2ec84fSDave May if (depth == max) { 140e2ec84fSDave May PetscReal cx[2]; 150e2ec84fSDave May 160e2ec84fSDave May cx[0] = (v1[0] + v2[0] + v3[0])/3.0; 170e2ec84fSDave May cx[1] = (v1[1] + v2[1] + v3[1])/3.0; 180e2ec84fSDave May 190e2ec84fSDave May xi[2*(*np)+0] = cx[0]; 200e2ec84fSDave May xi[2*(*np)+1] = cx[1]; 210e2ec84fSDave May (*np)++; 220e2ec84fSDave May PetscFunctionReturn(0); 230e2ec84fSDave May } 240e2ec84fSDave May 250e2ec84fSDave May /* calculate midpoints of each side */ 260e2ec84fSDave May for (i = 0; i < 2; i++) { 270e2ec84fSDave May v12[i] = (v1[i]+v2[i])/2.0; 280e2ec84fSDave May v23[i] = (v2[i]+v3[i])/2.0; 290e2ec84fSDave May v31[i] = (v3[i]+v1[i])/2.0; 300e2ec84fSDave May } 310e2ec84fSDave May 320e2ec84fSDave May /* recursively subdivide new triangles */ 330e2ec84fSDave May ierr = subdivide_triangle(v1,v12,v31,depth+1,max,xi,np);CHKERRQ(ierr); 340e2ec84fSDave May ierr = subdivide_triangle(v2,v23,v12,depth+1,max,xi,np);CHKERRQ(ierr); 350e2ec84fSDave May ierr = subdivide_triangle(v3,v31,v23,depth+1,max,xi,np);CHKERRQ(ierr); 360e2ec84fSDave May ierr = subdivide_triangle(v12,v23,v31,depth+1,max,xi,np);CHKERRQ(ierr); 370e2ec84fSDave May PetscFunctionReturn(0); 380e2ec84fSDave May } 390e2ec84fSDave May 400e2ec84fSDave May 410e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub) 420e2ec84fSDave May { 430e2ec84fSDave May PetscErrorCode ierr; 440e2ec84fSDave May const PetscInt dim = 2; 450e2ec84fSDave May PetscInt q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth; 460e2ec84fSDave May PetscReal *xi; 470e2ec84fSDave May PetscReal **basis; 480e2ec84fSDave May Vec coorlocal; 490e2ec84fSDave May PetscSection coordSection; 500e2ec84fSDave May PetscScalar *elcoor = NULL; 510e2ec84fSDave May PetscReal *swarm_coor; 520e2ec84fSDave May PetscInt *swarm_cellid; 530e2ec84fSDave May PetscReal v1[2],v2[2],v3[2]; 540e2ec84fSDave May 550e2ec84fSDave May PetscFunctionBegin; 560e2ec84fSDave May 570e2ec84fSDave May npoints_q = 1; 580e2ec84fSDave May for (d=0; d<nsub; d++) { npoints_q *= 4; } 590e2ec84fSDave May ierr = PetscMalloc1(dim*npoints_q,&xi);CHKERRQ(ierr); 600e2ec84fSDave May 610e2ec84fSDave May v1[0] = 0.0; v1[1] = 0.0; 620e2ec84fSDave May v2[0] = 1.0; v2[1] = 0.0; 630e2ec84fSDave May v3[0] = 0.0; v3[1] = 1.0; 640e2ec84fSDave May depth = 0; 650e2ec84fSDave May pcnt = 0; 660e2ec84fSDave May ierr = subdivide_triangle(v1,v2,v3,depth,nsub,xi,&pcnt);CHKERRQ(ierr); 670e2ec84fSDave May 680e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 690e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 700e2ec84fSDave May for (q=0; q<npoints_q; q++) { 710e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 720e2ec84fSDave May 730e2ec84fSDave May basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1]; 740e2ec84fSDave May basis[q][1] = xi[dim*q+0]; 750e2ec84fSDave May basis[q][2] = xi[dim*q+1]; 760e2ec84fSDave May } 770e2ec84fSDave May 780e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 790e2ec84fSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 800e2ec84fSDave May nel = pe - ps; 810e2ec84fSDave May 820e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 830e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 840e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 850e2ec84fSDave May 860e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 870e2ec84fSDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 880e2ec84fSDave May 890e2ec84fSDave May pcnt = 0; 900e2ec84fSDave May for (e=0; e<nel; e++) { 910e2ec84fSDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 920e2ec84fSDave May 930e2ec84fSDave May for (q=0; q<npoints_q; q++) { 940e2ec84fSDave May for (d=0; d<dim; d++) { 950e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 960e2ec84fSDave May for (k=0; k<npe; k++) { 970ca99263SDave May swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]); 980e2ec84fSDave May } 990e2ec84fSDave May } 1000e2ec84fSDave May swarm_cellid[pcnt] = e; 1010e2ec84fSDave May pcnt++; 1020e2ec84fSDave May } 1030e2ec84fSDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 1040e2ec84fSDave May } 1050e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1060e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1070e2ec84fSDave May 1080e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 1090e2ec84fSDave May for (q=0; q<npoints_q; q++) { 1100e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 1110e2ec84fSDave May } 1120e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 1130e2ec84fSDave May 1140e2ec84fSDave May PetscFunctionReturn(0); 1150e2ec84fSDave May } 1160e2ec84fSDave May 11778185223SDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm,DM dmc,PetscInt nsub) 11878185223SDave May { 11978185223SDave May PetscErrorCode ierr; 12078185223SDave May PetscInt dim,nfaces,nbasis; 12178185223SDave May PetscInt q,npoints_q,e,nel,pcnt,ps,pe,d,k,r; 12278185223SDave May PetscReal *B; 12378185223SDave May Vec coorlocal; 12478185223SDave May PetscSection coordSection; 12578185223SDave May PetscScalar *elcoor = NULL; 12678185223SDave May PetscReal *swarm_coor; 12778185223SDave May PetscInt *swarm_cellid; 12878185223SDave May const PetscReal *xiq; 12978185223SDave May PetscQuadrature quadrature; 13078185223SDave May PetscFE fe,feRef; 13178185223SDave May PetscBool is_simplex; 13278185223SDave May 13378185223SDave May PetscFunctionBegin; 13478185223SDave May 13578185223SDave May ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr); 13678185223SDave May 13778185223SDave May is_simplex = PETSC_FALSE; 13878185223SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 13978185223SDave May ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr); 14078185223SDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 14178185223SDave May 14278185223SDave May ierr = PetscFECreateDefault(dmc, dim, 1, is_simplex, NULL, 0, &fe);CHKERRQ(ierr); 14378185223SDave May 14478185223SDave May for (r=0; r<nsub; r++) { 14578185223SDave May ierr = PetscFERefine(fe,&feRef);CHKERRQ(ierr); 14678185223SDave May ierr = PetscFEGetQuadrature(feRef,&quadrature);CHKERRQ(ierr); 14778185223SDave May ierr = PetscFESetQuadrature(fe,quadrature);CHKERRQ(ierr); 14878185223SDave May ierr = PetscFEDestroy(&feRef);CHKERRQ(ierr); 14978185223SDave May } 15078185223SDave May 15178185223SDave May ierr = PetscFEGetQuadrature(fe,&quadrature);CHKERRQ(ierr); 15278185223SDave May ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL);CHKERRQ(ierr); 15378185223SDave May ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr); 15478185223SDave May ierr = PetscFEGetDefaultTabulation(fe, &B, NULL, NULL);CHKERRQ(ierr); 15578185223SDave May 15678185223SDave May /* 0->cell, 1->edge, 2->vert */ 15778185223SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 15878185223SDave May nel = pe - ps; 15978185223SDave May 16078185223SDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 16178185223SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 16278185223SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 16378185223SDave May 16478185223SDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 16578185223SDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 16678185223SDave May 16778185223SDave May pcnt = 0; 16878185223SDave May for (e=0; e<nel; e++) { 16978185223SDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 17078185223SDave May 17178185223SDave May for (q=0; q<npoints_q; q++) { 17278185223SDave May for (d=0; d<dim; d++) { 17378185223SDave May swarm_coor[dim*pcnt+d] = 0.0; 17478185223SDave May for (k=0; k<nbasis; k++) { 17578185223SDave May swarm_coor[dim*pcnt+d] += B[q*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); 17678185223SDave May } 17778185223SDave May } 17878185223SDave May swarm_cellid[pcnt] = e; 17978185223SDave May pcnt++; 18078185223SDave May } 181*92e40656SDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 18278185223SDave May } 18378185223SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 18478185223SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 18578185223SDave May 18678185223SDave May ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 18778185223SDave May PetscFunctionReturn(0); 18878185223SDave May } 18978185223SDave May 1900e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints) 1910e2ec84fSDave May { 1920e2ec84fSDave May PetscErrorCode ierr; 1930e2ec84fSDave May const PetscInt dim = 2; 1940e2ec84fSDave May PetscInt ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k; 1950e2ec84fSDave May PetscReal *xi,ds,ds2; 1960e2ec84fSDave May PetscReal **basis; 1970e2ec84fSDave May Vec coorlocal; 1980e2ec84fSDave May PetscSection coordSection; 1990e2ec84fSDave May PetscScalar *elcoor = NULL; 2000e2ec84fSDave May PetscReal *swarm_coor; 2010e2ec84fSDave May PetscInt *swarm_cellid; 2020e2ec84fSDave May 2030e2ec84fSDave May PetscFunctionBegin; 2040e2ec84fSDave May ierr = PetscMalloc1(dim*npoints*npoints,&xi);CHKERRQ(ierr); 2050e2ec84fSDave May pcnt = 0; 2060e2ec84fSDave May ds = 1.0/((PetscReal)(npoints-1)); 2070e2ec84fSDave May ds2 = 1.0/((PetscReal)(npoints)); 2080e2ec84fSDave May for (jj = 0; jj<npoints; jj++) { 2090e2ec84fSDave May for (ii=0; ii<npoints-jj; ii++) { 2100e2ec84fSDave May xi[dim*pcnt+0] = ii * ds; 2110e2ec84fSDave May xi[dim*pcnt+1] = jj * ds; 2120e2ec84fSDave May 2130e2ec84fSDave May xi[dim*pcnt+0] *= (1.0 - 1.2*ds2); 2140e2ec84fSDave May xi[dim*pcnt+1] *= (1.0 - 1.2*ds2); 2150e2ec84fSDave May 2160e2ec84fSDave May xi[dim*pcnt+0] += 0.35*ds2; 2170e2ec84fSDave May xi[dim*pcnt+1] += 0.35*ds2; 2180e2ec84fSDave May 2190e2ec84fSDave May pcnt++; 2200e2ec84fSDave May } 2210e2ec84fSDave May } 2220e2ec84fSDave May npoints_q = pcnt; 2230e2ec84fSDave May 2240e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 2250e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 2260e2ec84fSDave May for (q=0; q<npoints_q; q++) { 2270e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 2280e2ec84fSDave May 2290e2ec84fSDave May basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1]; 2300e2ec84fSDave May basis[q][1] = xi[dim*q+0]; 2310e2ec84fSDave May basis[q][2] = xi[dim*q+1]; 2320e2ec84fSDave May } 2330e2ec84fSDave May 2340e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 2350e2ec84fSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 2360e2ec84fSDave May nel = pe - ps; 2370e2ec84fSDave May 2380e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 2390e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 2400e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 2410e2ec84fSDave May 2420e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 2430e2ec84fSDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 2440e2ec84fSDave May 2450e2ec84fSDave May pcnt = 0; 2460e2ec84fSDave May for (e=0; e<nel; e++) { 2470e2ec84fSDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 2480e2ec84fSDave May 2490e2ec84fSDave May for (q=0; q<npoints_q; q++) { 2500e2ec84fSDave May for (d=0; d<dim; d++) { 2510e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 2520e2ec84fSDave May for (k=0; k<npe; k++) { 2530ca99263SDave May swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]); 2540e2ec84fSDave May } 2550e2ec84fSDave May } 2560e2ec84fSDave May swarm_cellid[pcnt] = e; 2570e2ec84fSDave May pcnt++; 2580e2ec84fSDave May } 2590e2ec84fSDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 2600e2ec84fSDave May } 2610e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 2620e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 2630e2ec84fSDave May 2640e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 2650e2ec84fSDave May for (q=0; q<npoints_q; q++) { 2660e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 2670e2ec84fSDave May } 2680e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 2690e2ec84fSDave May 2700e2ec84fSDave May PetscFunctionReturn(0); 2710e2ec84fSDave May } 2720e2ec84fSDave May 2730e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param) 2740e2ec84fSDave May { 2750e2ec84fSDave May PetscErrorCode ierr; 2760e2ec84fSDave May PetscInt dim; 2770e2ec84fSDave May 2780e2ec84fSDave May PetscFunctionBegin; 2790e2ec84fSDave May ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr); 2800e2ec84fSDave May switch (layout) { 2810e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 28278185223SDave May if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No 3D support for REGULAR+PLEX"); 2830e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm,celldm,layout_param);CHKERRQ(ierr); 2840e2ec84fSDave May break; 2850e2ec84fSDave May case DMSWARMPIC_LAYOUT_GAUSS: 2860e2ec84fSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Gauss layout not supported for PLEX"); 2870e2ec84fSDave May break; 2880e2ec84fSDave May case DMSWARMPIC_LAYOUT_SUBDIVISION: 28978185223SDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm,celldm,layout_param);CHKERRQ(ierr); 2900e2ec84fSDave May break; 2910e2ec84fSDave May } 2920e2ec84fSDave May PetscFunctionReturn(0); 2930e2ec84fSDave May } 294e3cd5995SDave May 295e3cd5995SDave May /* 296e3cd5995SDave May typedef struct { 297e3cd5995SDave May PetscReal x,y; 298e3cd5995SDave May } Point2d; 299e3cd5995SDave May 300e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s) 301e3cd5995SDave May { 302e3cd5995SDave May PetscFunctionBegin; 303e3cd5995SDave May *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y); 304e3cd5995SDave May PetscFunctionReturn(0); 305e3cd5995SDave May } 306e3cd5995SDave May */ 307e3cd5995SDave May /* 308e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v) 309e3cd5995SDave May { 310e3cd5995SDave May PetscReal s1,s2,s3; 311e3cd5995SDave May PetscBool b1, b2, b3; 312e3cd5995SDave May 313e3cd5995SDave May PetscFunctionBegin; 314e3cd5995SDave May signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f; 315e3cd5995SDave May signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f; 316e3cd5995SDave May signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f; 317e3cd5995SDave May 318e3cd5995SDave May *v = ((b1 == b2) && (b2 == b3)); 319e3cd5995SDave May PetscFunctionReturn(0); 320e3cd5995SDave May } 321e3cd5995SDave May */ 322e3cd5995SDave May /* 323e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ) 324e3cd5995SDave May { 325e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 326e3cd5995SDave May PetscReal c,b[2],A[2][2],inv[2][2],detJ,od; 327e3cd5995SDave May 328e3cd5995SDave May PetscFunctionBegin; 329e3cd5995SDave May x1 = coords[2*0+0]; 330e3cd5995SDave May x2 = coords[2*1+0]; 331e3cd5995SDave May x3 = coords[2*2+0]; 332e3cd5995SDave May 333e3cd5995SDave May y1 = coords[2*0+1]; 334e3cd5995SDave May y2 = coords[2*1+1]; 335e3cd5995SDave May y3 = coords[2*2+1]; 336e3cd5995SDave May 337e3cd5995SDave May c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3; 338e3cd5995SDave May b[0] = xp[0] - c; 339e3cd5995SDave May c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3; 340e3cd5995SDave May b[1] = xp[1] - c; 341e3cd5995SDave May 342e3cd5995SDave May A[0][0] = -0.5*x1 + 0.5*x2; A[0][1] = -0.5*x1 + 0.5*x3; 343e3cd5995SDave May A[1][0] = -0.5*y1 + 0.5*y2; A[1][1] = -0.5*y1 + 0.5*y3; 344e3cd5995SDave May 345e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 346e3cd5995SDave May *dJ = PetscAbsReal(detJ); 347e3cd5995SDave May od = 1.0/detJ; 348e3cd5995SDave May 349e3cd5995SDave May inv[0][0] = A[1][1] * od; 350e3cd5995SDave May inv[0][1] = -A[0][1] * od; 351e3cd5995SDave May inv[1][0] = -A[1][0] * od; 352e3cd5995SDave May inv[1][1] = A[0][0] * od; 353e3cd5995SDave May 354e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 355e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 356e3cd5995SDave May PetscFunctionReturn(0); 357e3cd5995SDave May } 358e3cd5995SDave May */ 359e3cd5995SDave May 360a5f152d1SDave May static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscScalar coords[],PetscReal xip[],PetscReal *dJ) 361e3cd5995SDave May { 362e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 363e3cd5995SDave May PetscReal b[2],A[2][2],inv[2][2],detJ,od; 364e3cd5995SDave May 365e3cd5995SDave May PetscFunctionBegin; 366a5f152d1SDave May x1 = PetscRealPart(coords[2*0+0]); 367a5f152d1SDave May x2 = PetscRealPart(coords[2*1+0]); 368a5f152d1SDave May x3 = PetscRealPart(coords[2*2+0]); 369e3cd5995SDave May 370a5f152d1SDave May y1 = PetscRealPart(coords[2*0+1]); 371a5f152d1SDave May y2 = PetscRealPart(coords[2*1+1]); 372a5f152d1SDave May y3 = PetscRealPart(coords[2*2+1]); 373e3cd5995SDave May 374e3cd5995SDave May b[0] = xp[0] - x1; 375e3cd5995SDave May b[1] = xp[1] - y1; 376e3cd5995SDave May 377e3cd5995SDave May A[0][0] = x2-x1; A[0][1] = x3-x1; 378e3cd5995SDave May A[1][0] = y2-y1; A[1][1] = y3-y1; 379e3cd5995SDave May 380e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 381e3cd5995SDave May *dJ = PetscAbsReal(detJ); 382e3cd5995SDave May od = 1.0/detJ; 383e3cd5995SDave May 384e3cd5995SDave May inv[0][0] = A[1][1] * od; 385e3cd5995SDave May inv[0][1] = -A[0][1] * od; 386e3cd5995SDave May inv[1][0] = -A[1][0] * od; 387e3cd5995SDave May inv[1][1] = A[0][0] * od; 388e3cd5995SDave May 389e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 390e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 391e3cd5995SDave May PetscFunctionReturn(0); 392e3cd5995SDave May } 393e3cd5995SDave May 394e3cd5995SDave May PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field) 395e3cd5995SDave May { 396e3cd5995SDave May PetscErrorCode ierr; 397e3cd5995SDave May const PetscReal PLEX_C_EPS = 1.0e-8; 398e3cd5995SDave May Vec v_field_l,denom_l,coor_l,denom; 399e3cd5995SDave May PetscInt k,p,e,npoints; 400e3cd5995SDave May PetscInt *mpfield_cell; 401e3cd5995SDave May PetscReal *mpfield_coor; 402a5f152d1SDave May PetscReal xi_p[2]; 403a5f152d1SDave May PetscScalar Ni[3]; 404e3cd5995SDave May PetscSection coordSection; 405e3cd5995SDave May PetscScalar *elcoor = NULL; 406e3cd5995SDave May 407e3cd5995SDave May PetscFunctionBegin; 408e3cd5995SDave May ierr = VecZeroEntries(v_field);CHKERRQ(ierr); 409e3cd5995SDave May 410e3cd5995SDave May ierr = DMGetLocalVector(dm,&v_field_l);CHKERRQ(ierr); 411e3cd5995SDave May ierr = DMGetGlobalVector(dm,&denom);CHKERRQ(ierr); 412e3cd5995SDave May ierr = DMGetLocalVector(dm,&denom_l);CHKERRQ(ierr); 413e3cd5995SDave May ierr = VecZeroEntries(v_field_l);CHKERRQ(ierr); 414e3cd5995SDave May ierr = VecZeroEntries(denom);CHKERRQ(ierr); 415e3cd5995SDave May ierr = VecZeroEntries(denom_l);CHKERRQ(ierr); 416e3cd5995SDave May 417e3cd5995SDave May ierr = DMGetCoordinatesLocal(dm,&coor_l);CHKERRQ(ierr); 418e3cd5995SDave May ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 419e3cd5995SDave May 420e3cd5995SDave May ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr); 421e3cd5995SDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 422e3cd5995SDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 423e3cd5995SDave May 424e3cd5995SDave May for (p=0; p<npoints; p++) { 425a5f152d1SDave May PetscReal *coor_p,dJ; 426a5f152d1SDave May PetscScalar elfield[3]; 427e3cd5995SDave May PetscBool point_located; 428e3cd5995SDave May 429e3cd5995SDave May e = mpfield_cell[p]; 430e3cd5995SDave May coor_p = &mpfield_coor[2*p]; 431e3cd5995SDave May 432e3cd5995SDave May ierr = DMPlexVecGetClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr); 433e3cd5995SDave May 434e3cd5995SDave May /* 435e3cd5995SDave May while (!point_located && (failed_counter < 25)) { 436e3cd5995SDave May ierr = PointInTriangle(point, coords[0], coords[1], coords[2], &point_located);CHKERRQ(ierr); 437e3cd5995SDave May point.x = coor_p[0]; 438e3cd5995SDave May point.y = coor_p[1]; 439e3cd5995SDave May point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 440e3cd5995SDave May point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 441e3cd5995SDave May failed_counter++; 442e3cd5995SDave May } 443e3cd5995SDave May 444e3cd5995SDave May if (!point_located){ 445e3cd5995SDave 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); 446e3cd5995SDave May } 447e3cd5995SDave May 448e3cd5995SDave 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); 449e3cd5995SDave May else { 450e3cd5995SDave May ierr = _ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr); 451e3cd5995SDave May xi_p[0] = 0.5*(xi_p[0] + 1.0); 452e3cd5995SDave May xi_p[1] = 0.5*(xi_p[1] + 1.0); 453e3cd5995SDave May 454e3cd5995SDave 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]); 455e3cd5995SDave May 456e3cd5995SDave May } 457e3cd5995SDave May */ 458e3cd5995SDave May 459e3cd5995SDave May ierr = ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr); 460e3cd5995SDave May /* 461e3cd5995SDave 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]); 462e3cd5995SDave May */ 463e3cd5995SDave May /* 464e3cd5995SDave May point_located = PETSC_TRUE; 465e3cd5995SDave May if (xi_p[0] < 0.0) { 466e3cd5995SDave May if (xi_p[0] > -PLEX_C_EPS) { 467e3cd5995SDave May xi_p[0] = 0.0; 468e3cd5995SDave May } else { 469e3cd5995SDave May point_located = PETSC_FALSE; 470e3cd5995SDave May } 471e3cd5995SDave May } 472e3cd5995SDave May if (xi_p[1] < 0.0) { 473e3cd5995SDave May if (xi_p[1] > -PLEX_C_EPS) { 474e3cd5995SDave May xi_p[1] = 0.0; 475e3cd5995SDave May } else { 476e3cd5995SDave May point_located = PETSC_FALSE; 477e3cd5995SDave May } 478e3cd5995SDave May } 479e3cd5995SDave May if (xi_p[1] > (1.0-xi_p[0])) { 480e3cd5995SDave May if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { 481e3cd5995SDave May xi_p[1] = 1.0 - xi_p[0]; 482e3cd5995SDave May } else { 483e3cd5995SDave May point_located = PETSC_FALSE; 484e3cd5995SDave May } 485e3cd5995SDave May } 486e3cd5995SDave May if (!point_located){ 487e3cd5995SDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); 488e3cd5995SDave 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]); 489e3cd5995SDave May } 490e3cd5995SDave 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); 491e3cd5995SDave May */ 492e3cd5995SDave May 493e3cd5995SDave May Ni[0] = 1.0 - xi_p[0] - xi_p[1]; 494e3cd5995SDave May Ni[1] = xi_p[0]; 495e3cd5995SDave May Ni[2] = xi_p[1]; 496e3cd5995SDave May 497e3cd5995SDave May point_located = PETSC_TRUE; 498e3cd5995SDave May for (k=0; k<3; k++) { 499a5f152d1SDave May if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; 500a5f152d1SDave May if (PetscRealPart(Ni[k]) > (1.0+PLEX_C_EPS)) point_located = PETSC_FALSE; 501e3cd5995SDave May } 502e3cd5995SDave May if (!point_located){ 503d81390bbSDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",(double)xi_p[0],(double)xi_p[1]); 504d81390bbSDave 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])); 505e3cd5995SDave May } 506d81390bbSDave 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); 507e3cd5995SDave May 508e3cd5995SDave May for (k=0; k<3; k++) { 509e3cd5995SDave May Ni[k] = Ni[k] * dJ; 510e3cd5995SDave May elfield[k] = Ni[k] * swarm_field[p]; 511e3cd5995SDave May } 512e3cd5995SDave May 513e3cd5995SDave May ierr = DMPlexVecRestoreClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr); 514e3cd5995SDave May 515e3cd5995SDave May ierr = DMPlexVecSetClosure(dm, NULL,v_field_l, e, elfield, ADD_VALUES);CHKERRQ(ierr); 516e3cd5995SDave May ierr = DMPlexVecSetClosure(dm, NULL,denom_l, e, Ni, ADD_VALUES);CHKERRQ(ierr); 517e3cd5995SDave May } 518e3cd5995SDave May 519e3cd5995SDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 520e3cd5995SDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 521e3cd5995SDave May 522e3cd5995SDave May ierr = DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 523e3cd5995SDave May ierr = DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 524e3cd5995SDave May ierr = DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 525e3cd5995SDave May ierr = DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 526e3cd5995SDave May 527e3cd5995SDave May ierr = VecPointwiseDivide(v_field,v_field,denom);CHKERRQ(ierr); 528e3cd5995SDave May 529e3cd5995SDave May ierr = DMRestoreLocalVector(dm,&v_field_l);CHKERRQ(ierr); 530e3cd5995SDave May ierr = DMRestoreLocalVector(dm,&denom_l);CHKERRQ(ierr); 531e3cd5995SDave May ierr = DMRestoreGlobalVector(dm,&denom);CHKERRQ(ierr); 532e3cd5995SDave May 533e3cd5995SDave May PetscFunctionReturn(0); 534e3cd5995SDave May } 535e3cd5995SDave May 536e3cd5995SDave May PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DataField dfield[],Vec vecs[]) 537e3cd5995SDave May { 538e3cd5995SDave May PetscErrorCode ierr; 539e3cd5995SDave May PetscInt f,dim; 540e3cd5995SDave May 541e3cd5995SDave May PetscFunctionBegin; 542e3cd5995SDave May ierr = DMGetDimension(swarm,&dim);CHKERRQ(ierr); 543e3cd5995SDave May switch (dim) { 544e3cd5995SDave May case 2: 545e3cd5995SDave May for (f=0; f<nfields; f++) { 546e3cd5995SDave May PetscReal *swarm_field; 547e3cd5995SDave May 548e3cd5995SDave May ierr = DataFieldGetEntries(dfield[f],(void**)&swarm_field);CHKERRQ(ierr); 549e3cd5995SDave May ierr = DMSwarmProjectField_ApproxP1_PLEX_2D(swarm,swarm_field,celldm,vecs[f]);CHKERRQ(ierr); 550e3cd5995SDave May } 551e3cd5995SDave May break; 552e3cd5995SDave May case 3: 553e3cd5995SDave May SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D"); 554e3cd5995SDave May break; 555e3cd5995SDave May default: 556e3cd5995SDave May break; 557e3cd5995SDave May } 558e3cd5995SDave May 559e3cd5995SDave May PetscFunctionReturn(0); 560e3cd5995SDave May } 561431382f2SDave May 562431382f2SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm,DM dmc,PetscInt npoints,PetscReal xi[]) 563431382f2SDave May { 564431382f2SDave May PetscBool is_simplex,is_tensorcell; 565431382f2SDave May PetscErrorCode ierr; 566*92e40656SDave May PetscInt dim,nfaces,ps,pe,p,d,nbasis,pcnt,e,k,nel; 567*92e40656SDave May PetscFE fe; 568*92e40656SDave May PetscQuadrature quadrature; 569*92e40656SDave May PetscReal *B; 570*92e40656SDave May Vec coorlocal; 571*92e40656SDave May PetscSection coordSection; 572*92e40656SDave May PetscScalar *elcoor = NULL; 573*92e40656SDave May PetscReal *swarm_coor; 574*92e40656SDave May PetscInt *swarm_cellid; 575431382f2SDave May 576*92e40656SDave May PetscFunctionBegin; 577431382f2SDave May ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr); 578431382f2SDave May 579431382f2SDave May is_simplex = PETSC_FALSE; 580431382f2SDave May is_tensorcell = PETSC_FALSE; 581431382f2SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 582431382f2SDave May ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr); 583431382f2SDave May 584431382f2SDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 585431382f2SDave May 586431382f2SDave May switch (dim) { 587431382f2SDave May case 2: 588431382f2SDave May if (nfaces == 4) { is_tensorcell = PETSC_TRUE; } 589431382f2SDave May break; 590431382f2SDave May case 3: 591431382f2SDave May if (nfaces == 6) { is_tensorcell = PETSC_TRUE; } 592431382f2SDave May break; 593431382f2SDave May default: 594431382f2SDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for 2D, 3D"); 595431382f2SDave May break; 596431382f2SDave May } 597431382f2SDave May 598*92e40656SDave May /* check points provided fail inside the reference cell */ 599431382f2SDave May if (is_simplex) { 600*92e40656SDave May for (p=0; p<npoints; p++) { 601*92e40656SDave May PetscReal sum; 602*92e40656SDave May for (d=0; d<dim; d++) { 603*92e40656SDave May if (xi[dim*p+d] < -1.0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain"); 604*92e40656SDave May } 605*92e40656SDave May sum = 0.0; 606*92e40656SDave May for (d=0; d<dim; d++) { 607*92e40656SDave May sum += xi[dim*p+d]; 608*92e40656SDave May } 609*92e40656SDave May if (sum > 0.0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain"); 610*92e40656SDave May } 611431382f2SDave May } else if (is_tensorcell) { 612*92e40656SDave May for (p=0; p<npoints; p++) { 613*92e40656SDave May for (d=0; d<dim; d++) { 614*92e40656SDave 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"); 615*92e40656SDave May } 616*92e40656SDave May } 617431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for d-simplex and d-tensorcell"); 618431382f2SDave May 619*92e40656SDave May ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)dm),&quadrature);CHKERRQ(ierr); 620*92e40656SDave May ierr = PetscQuadratureSetData(quadrature,dim,1,npoints,(const PetscReal*)xi,NULL);CHKERRQ(ierr); 621*92e40656SDave May ierr = PetscFECreateDefault(dmc, dim, 1, is_simplex, NULL, 0, &fe);CHKERRQ(ierr); 622*92e40656SDave May ierr = PetscFESetQuadrature(fe,quadrature);CHKERRQ(ierr); 623*92e40656SDave May ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr); 624*92e40656SDave May ierr = PetscFEGetDefaultTabulation(fe, &B, NULL, NULL);CHKERRQ(ierr); 625*92e40656SDave May 626*92e40656SDave May /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ 627*92e40656SDave May /* 0->cell, 1->edge, 2->vert */ 628*92e40656SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 629*92e40656SDave May nel = pe - ps; 630*92e40656SDave May 631*92e40656SDave May ierr = DMSwarmSetLocalSizes(dm,npoints*nel,-1);CHKERRQ(ierr); 632*92e40656SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 633*92e40656SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 634*92e40656SDave May 635*92e40656SDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 636*92e40656SDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 637*92e40656SDave May 638*92e40656SDave May pcnt = 0; 639*92e40656SDave May for (e=0; e<nel; e++) { 640*92e40656SDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 641*92e40656SDave May 642*92e40656SDave May for (p=0; p<npoints; p++) { 643*92e40656SDave May for (d=0; d<dim; d++) { 644*92e40656SDave May swarm_coor[dim*pcnt+d] = 0.0; 645*92e40656SDave May for (k=0; k<nbasis; k++) { 646*92e40656SDave May swarm_coor[dim*pcnt+d] += B[p*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); 647*92e40656SDave May } 648*92e40656SDave May } 649*92e40656SDave May swarm_cellid[pcnt] = e; 650*92e40656SDave May pcnt++; 651*92e40656SDave May } 652*92e40656SDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 653*92e40656SDave May } 654*92e40656SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 655*92e40656SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 656*92e40656SDave May 657*92e40656SDave May ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr); 658*92e40656SDave May ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 659*92e40656SDave May 660431382f2SDave May PetscFunctionReturn(0); 661431382f2SDave May } 662