10e2ec84fSDave May #include <petscdm.h> 20e2ec84fSDave May #include <petscdmplex.h> 30e2ec84fSDave May #include <petscdmswarm.h> 4e3cd5995SDave May #include "data_bucket.h" 50e2ec84fSDave May 6*2cf64d79SDave May static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem) 7*2cf64d79SDave May { 8*2cf64d79SDave May const PetscInt Nc = 1; 9*2cf64d79SDave May PetscQuadrature q, fq; 10*2cf64d79SDave May DM K; 11*2cf64d79SDave May PetscSpace P; 12*2cf64d79SDave May PetscDualSpace Q; 13*2cf64d79SDave May PetscInt order, quadPointsPerEdge; 14*2cf64d79SDave May PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 15*2cf64d79SDave May PetscErrorCode ierr; 16*2cf64d79SDave May 17*2cf64d79SDave May PetscFunctionBegin; 18*2cf64d79SDave May /* Create space */ 19*2cf64d79SDave May ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P);CHKERRQ(ierr); 20*2cf64d79SDave May /*ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);*/ 21*2cf64d79SDave May ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr); 22*2cf64d79SDave May /*ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);*/ 23*2cf64d79SDave May ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 24*2cf64d79SDave May ierr = PetscSpaceSetOrder(P,1);CHKERRQ(ierr); 25*2cf64d79SDave May ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr); 26*2cf64d79SDave May ierr = PetscSpacePolynomialSetNumVariables(P, dim);CHKERRQ(ierr); 27*2cf64d79SDave May ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 28*2cf64d79SDave May ierr = PetscSpaceGetOrder(P, &order);CHKERRQ(ierr); 29*2cf64d79SDave May ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr); 30*2cf64d79SDave May /* Create dual space */ 31*2cf64d79SDave May ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q);CHKERRQ(ierr); 32*2cf64d79SDave May ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 33*2cf64d79SDave May /*ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);*/ 34*2cf64d79SDave May ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 35*2cf64d79SDave May ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 36*2cf64d79SDave May ierr = DMDestroy(&K);CHKERRQ(ierr); 37*2cf64d79SDave May ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr); 38*2cf64d79SDave May ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr); 39*2cf64d79SDave May ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr); 40*2cf64d79SDave May /*ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);*/ 41*2cf64d79SDave May ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 42*2cf64d79SDave May ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 43*2cf64d79SDave May /* Create element */ 44*2cf64d79SDave May ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem);CHKERRQ(ierr); 45*2cf64d79SDave May /*ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);*/ 46*2cf64d79SDave May /*ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);*/ 47*2cf64d79SDave May ierr = PetscFESetType(*fem,PETSCFEBASIC);CHKERRQ(ierr); 48*2cf64d79SDave May ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr); 49*2cf64d79SDave May ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr); 50*2cf64d79SDave May ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr); 51*2cf64d79SDave May ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 52*2cf64d79SDave May ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 53*2cf64d79SDave May ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 54*2cf64d79SDave May /* Create quadrature (with specified order if given) */ 55*2cf64d79SDave May qorder = qorder >= 0 ? qorder : order; 56*2cf64d79SDave May quadPointsPerEdge = PetscMax(qorder + 1,1); 57*2cf64d79SDave May if (isSimplex) { 58*2cf64d79SDave May ierr = PetscDTGaussJacobiQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 59*2cf64d79SDave May ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 60*2cf64d79SDave May } 61*2cf64d79SDave May else { 62*2cf64d79SDave May ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 63*2cf64d79SDave May ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 64*2cf64d79SDave May } 65*2cf64d79SDave May ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr); 66*2cf64d79SDave May ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr); 67*2cf64d79SDave May ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 68*2cf64d79SDave May ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 69*2cf64d79SDave May PetscFunctionReturn(0); 70*2cf64d79SDave May } 71*2cf64d79SDave May 720e2ec84fSDave May PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[3],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np) 730e2ec84fSDave May { 740e2ec84fSDave May PetscReal v12[2],v23[2],v31[2]; 750e2ec84fSDave May PetscInt i; 760e2ec84fSDave May PetscErrorCode ierr; 770e2ec84fSDave May 780e2ec84fSDave May PetscFunctionBegin; 790e2ec84fSDave May if (depth == max) { 800e2ec84fSDave May PetscReal cx[2]; 810e2ec84fSDave May 820e2ec84fSDave May cx[0] = (v1[0] + v2[0] + v3[0])/3.0; 830e2ec84fSDave May cx[1] = (v1[1] + v2[1] + v3[1])/3.0; 840e2ec84fSDave May 850e2ec84fSDave May xi[2*(*np)+0] = cx[0]; 860e2ec84fSDave May xi[2*(*np)+1] = cx[1]; 870e2ec84fSDave May (*np)++; 880e2ec84fSDave May PetscFunctionReturn(0); 890e2ec84fSDave May } 900e2ec84fSDave May 910e2ec84fSDave May /* calculate midpoints of each side */ 920e2ec84fSDave May for (i = 0; i < 2; i++) { 930e2ec84fSDave May v12[i] = (v1[i]+v2[i])/2.0; 940e2ec84fSDave May v23[i] = (v2[i]+v3[i])/2.0; 950e2ec84fSDave May v31[i] = (v3[i]+v1[i])/2.0; 960e2ec84fSDave May } 970e2ec84fSDave May 980e2ec84fSDave May /* recursively subdivide new triangles */ 990e2ec84fSDave May ierr = subdivide_triangle(v1,v12,v31,depth+1,max,xi,np);CHKERRQ(ierr); 1000e2ec84fSDave May ierr = subdivide_triangle(v2,v23,v12,depth+1,max,xi,np);CHKERRQ(ierr); 1010e2ec84fSDave May ierr = subdivide_triangle(v3,v31,v23,depth+1,max,xi,np);CHKERRQ(ierr); 1020e2ec84fSDave May ierr = subdivide_triangle(v12,v23,v31,depth+1,max,xi,np);CHKERRQ(ierr); 1030e2ec84fSDave May PetscFunctionReturn(0); 1040e2ec84fSDave May } 1050e2ec84fSDave May 1060e2ec84fSDave May 1070e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub) 1080e2ec84fSDave May { 1090e2ec84fSDave May PetscErrorCode ierr; 1100e2ec84fSDave May const PetscInt dim = 2; 1110e2ec84fSDave May PetscInt q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth; 1120e2ec84fSDave May PetscReal *xi; 1130e2ec84fSDave May PetscReal **basis; 1140e2ec84fSDave May Vec coorlocal; 1150e2ec84fSDave May PetscSection coordSection; 1160e2ec84fSDave May PetscScalar *elcoor = NULL; 1170e2ec84fSDave May PetscReal *swarm_coor; 1180e2ec84fSDave May PetscInt *swarm_cellid; 1190e2ec84fSDave May PetscReal v1[2],v2[2],v3[2]; 1200e2ec84fSDave May 1210e2ec84fSDave May PetscFunctionBegin; 1220e2ec84fSDave May 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 1800e2ec84fSDave May PetscFunctionReturn(0); 1810e2ec84fSDave May } 1820e2ec84fSDave May 18378185223SDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm,DM dmc,PetscInt nsub) 18478185223SDave May { 18578185223SDave May PetscErrorCode ierr; 18678185223SDave May PetscInt dim,nfaces,nbasis; 18778185223SDave May PetscInt q,npoints_q,e,nel,pcnt,ps,pe,d,k,r; 18878185223SDave May PetscReal *B; 18978185223SDave May Vec coorlocal; 19078185223SDave May PetscSection coordSection; 19178185223SDave May PetscScalar *elcoor = NULL; 19278185223SDave May PetscReal *swarm_coor; 19378185223SDave May PetscInt *swarm_cellid; 19478185223SDave May const PetscReal *xiq; 19578185223SDave May PetscQuadrature quadrature; 19678185223SDave May PetscFE fe,feRef; 19778185223SDave May PetscBool is_simplex; 19878185223SDave May 19978185223SDave May PetscFunctionBegin; 20078185223SDave May 20178185223SDave May ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr); 20278185223SDave May 20378185223SDave May is_simplex = PETSC_FALSE; 20478185223SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 20578185223SDave May ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr); 20678185223SDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 20778185223SDave May 208*2cf64d79SDave May ierr = private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe);CHKERRQ(ierr); 20978185223SDave May 21078185223SDave May for (r=0; r<nsub; r++) { 21178185223SDave May ierr = PetscFERefine(fe,&feRef);CHKERRQ(ierr); 21278185223SDave May ierr = PetscFEGetQuadrature(feRef,&quadrature);CHKERRQ(ierr); 21378185223SDave May ierr = PetscFESetQuadrature(fe,quadrature);CHKERRQ(ierr); 21478185223SDave May ierr = PetscFEDestroy(&feRef);CHKERRQ(ierr); 21578185223SDave May } 21678185223SDave May 21778185223SDave May ierr = PetscFEGetQuadrature(fe,&quadrature);CHKERRQ(ierr); 21878185223SDave May ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL);CHKERRQ(ierr); 21978185223SDave May ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr); 22078185223SDave May ierr = PetscFEGetDefaultTabulation(fe, &B, NULL, NULL);CHKERRQ(ierr); 22178185223SDave May 22278185223SDave May /* 0->cell, 1->edge, 2->vert */ 22378185223SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 22478185223SDave May nel = pe - ps; 22578185223SDave May 22678185223SDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 22778185223SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 22878185223SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 22978185223SDave May 23078185223SDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 23178185223SDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 23278185223SDave May 23378185223SDave May pcnt = 0; 23478185223SDave May for (e=0; e<nel; e++) { 23578185223SDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 23678185223SDave May 23778185223SDave May for (q=0; q<npoints_q; q++) { 23878185223SDave May for (d=0; d<dim; d++) { 23978185223SDave May swarm_coor[dim*pcnt+d] = 0.0; 24078185223SDave May for (k=0; k<nbasis; k++) { 24178185223SDave May swarm_coor[dim*pcnt+d] += B[q*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); 24278185223SDave May } 24378185223SDave May } 24478185223SDave May swarm_cellid[pcnt] = e; 24578185223SDave May pcnt++; 24678185223SDave May } 24792e40656SDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 24878185223SDave May } 24978185223SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 25078185223SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 25178185223SDave May 25278185223SDave May ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 25378185223SDave May PetscFunctionReturn(0); 25478185223SDave May } 25578185223SDave May 2560e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints) 2570e2ec84fSDave May { 2580e2ec84fSDave May PetscErrorCode ierr; 2590e2ec84fSDave May const PetscInt dim = 2; 2600e2ec84fSDave May PetscInt ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k; 2610e2ec84fSDave May PetscReal *xi,ds,ds2; 2620e2ec84fSDave May PetscReal **basis; 2630e2ec84fSDave May Vec coorlocal; 2640e2ec84fSDave May PetscSection coordSection; 2650e2ec84fSDave May PetscScalar *elcoor = NULL; 2660e2ec84fSDave May PetscReal *swarm_coor; 2670e2ec84fSDave May PetscInt *swarm_cellid; 2680e2ec84fSDave May 2690e2ec84fSDave May PetscFunctionBegin; 2700e2ec84fSDave May ierr = PetscMalloc1(dim*npoints*npoints,&xi);CHKERRQ(ierr); 2710e2ec84fSDave May pcnt = 0; 2720e2ec84fSDave May ds = 1.0/((PetscReal)(npoints-1)); 2730e2ec84fSDave May ds2 = 1.0/((PetscReal)(npoints)); 2740e2ec84fSDave May for (jj = 0; jj<npoints; jj++) { 2750e2ec84fSDave May for (ii=0; ii<npoints-jj; ii++) { 2760e2ec84fSDave May xi[dim*pcnt+0] = ii * ds; 2770e2ec84fSDave May xi[dim*pcnt+1] = jj * ds; 2780e2ec84fSDave May 2790e2ec84fSDave May xi[dim*pcnt+0] *= (1.0 - 1.2*ds2); 2800e2ec84fSDave May xi[dim*pcnt+1] *= (1.0 - 1.2*ds2); 2810e2ec84fSDave May 2820e2ec84fSDave May xi[dim*pcnt+0] += 0.35*ds2; 2830e2ec84fSDave May xi[dim*pcnt+1] += 0.35*ds2; 2840e2ec84fSDave May 2850e2ec84fSDave May pcnt++; 2860e2ec84fSDave May } 2870e2ec84fSDave May } 2880e2ec84fSDave May npoints_q = pcnt; 2890e2ec84fSDave May 2900e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 2910e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 2920e2ec84fSDave May for (q=0; q<npoints_q; q++) { 2930e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 2940e2ec84fSDave May 2950e2ec84fSDave May basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1]; 2960e2ec84fSDave May basis[q][1] = xi[dim*q+0]; 2970e2ec84fSDave May basis[q][2] = xi[dim*q+1]; 2980e2ec84fSDave May } 2990e2ec84fSDave May 3000e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 3010e2ec84fSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 3020e2ec84fSDave May nel = pe - ps; 3030e2ec84fSDave May 3040e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 3050e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3060e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3070e2ec84fSDave May 3080e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 3090e2ec84fSDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 3100e2ec84fSDave May 3110e2ec84fSDave May pcnt = 0; 3120e2ec84fSDave May for (e=0; e<nel; e++) { 3130e2ec84fSDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 3140e2ec84fSDave May 3150e2ec84fSDave May for (q=0; q<npoints_q; q++) { 3160e2ec84fSDave May for (d=0; d<dim; d++) { 3170e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 3180e2ec84fSDave May for (k=0; k<npe; k++) { 3190ca99263SDave May swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]); 3200e2ec84fSDave May } 3210e2ec84fSDave May } 3220e2ec84fSDave May swarm_cellid[pcnt] = e; 3230e2ec84fSDave May pcnt++; 3240e2ec84fSDave May } 3250e2ec84fSDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 3260e2ec84fSDave May } 3270e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3280e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3290e2ec84fSDave May 3300e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 3310e2ec84fSDave May for (q=0; q<npoints_q; q++) { 3320e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 3330e2ec84fSDave May } 3340e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 3350e2ec84fSDave May 3360e2ec84fSDave May PetscFunctionReturn(0); 3370e2ec84fSDave May } 3380e2ec84fSDave May 3390e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param) 3400e2ec84fSDave May { 3410e2ec84fSDave May PetscErrorCode ierr; 3420e2ec84fSDave May PetscInt dim; 3430e2ec84fSDave May 3440e2ec84fSDave May PetscFunctionBegin; 3450e2ec84fSDave May ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr); 3460e2ec84fSDave May switch (layout) { 3470e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 34878185223SDave May if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No 3D support for REGULAR+PLEX"); 3490e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm,celldm,layout_param);CHKERRQ(ierr); 3500e2ec84fSDave May break; 3510e2ec84fSDave May case DMSWARMPIC_LAYOUT_GAUSS: 3520e2ec84fSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Gauss layout not supported for PLEX"); 3530e2ec84fSDave May break; 3540e2ec84fSDave May case DMSWARMPIC_LAYOUT_SUBDIVISION: 35578185223SDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm,celldm,layout_param);CHKERRQ(ierr); 3560e2ec84fSDave May break; 3570e2ec84fSDave May } 3580e2ec84fSDave May PetscFunctionReturn(0); 3590e2ec84fSDave May } 360e3cd5995SDave May 361e3cd5995SDave May /* 362e3cd5995SDave May typedef struct { 363e3cd5995SDave May PetscReal x,y; 364e3cd5995SDave May } Point2d; 365e3cd5995SDave May 366e3cd5995SDave May static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s) 367e3cd5995SDave May { 368e3cd5995SDave May PetscFunctionBegin; 369e3cd5995SDave May *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y); 370e3cd5995SDave May PetscFunctionReturn(0); 371e3cd5995SDave May } 372e3cd5995SDave May */ 373e3cd5995SDave May /* 374e3cd5995SDave May static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v) 375e3cd5995SDave May { 376e3cd5995SDave May PetscReal s1,s2,s3; 377e3cd5995SDave May PetscBool b1, b2, b3; 378e3cd5995SDave May 379e3cd5995SDave May PetscFunctionBegin; 380e3cd5995SDave May signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f; 381e3cd5995SDave May signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f; 382e3cd5995SDave May signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f; 383e3cd5995SDave May 384e3cd5995SDave May *v = ((b1 == b2) && (b2 == b3)); 385e3cd5995SDave May PetscFunctionReturn(0); 386e3cd5995SDave May } 387e3cd5995SDave May */ 388e3cd5995SDave May /* 389e3cd5995SDave May static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ) 390e3cd5995SDave May { 391e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 392e3cd5995SDave May PetscReal c,b[2],A[2][2],inv[2][2],detJ,od; 393e3cd5995SDave May 394e3cd5995SDave May PetscFunctionBegin; 395e3cd5995SDave May x1 = coords[2*0+0]; 396e3cd5995SDave May x2 = coords[2*1+0]; 397e3cd5995SDave May x3 = coords[2*2+0]; 398e3cd5995SDave May 399e3cd5995SDave May y1 = coords[2*0+1]; 400e3cd5995SDave May y2 = coords[2*1+1]; 401e3cd5995SDave May y3 = coords[2*2+1]; 402e3cd5995SDave May 403e3cd5995SDave May c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3; 404e3cd5995SDave May b[0] = xp[0] - c; 405e3cd5995SDave May c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3; 406e3cd5995SDave May b[1] = xp[1] - c; 407e3cd5995SDave May 408e3cd5995SDave May A[0][0] = -0.5*x1 + 0.5*x2; A[0][1] = -0.5*x1 + 0.5*x3; 409e3cd5995SDave May A[1][0] = -0.5*y1 + 0.5*y2; A[1][1] = -0.5*y1 + 0.5*y3; 410e3cd5995SDave May 411e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 412e3cd5995SDave May *dJ = PetscAbsReal(detJ); 413e3cd5995SDave May od = 1.0/detJ; 414e3cd5995SDave May 415e3cd5995SDave May inv[0][0] = A[1][1] * od; 416e3cd5995SDave May inv[0][1] = -A[0][1] * od; 417e3cd5995SDave May inv[1][0] = -A[1][0] * od; 418e3cd5995SDave May inv[1][1] = A[0][0] * od; 419e3cd5995SDave May 420e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 421e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 422e3cd5995SDave May PetscFunctionReturn(0); 423e3cd5995SDave May } 424e3cd5995SDave May */ 425e3cd5995SDave May 426a5f152d1SDave May static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscScalar coords[],PetscReal xip[],PetscReal *dJ) 427e3cd5995SDave May { 428e3cd5995SDave May PetscReal x1,y1,x2,y2,x3,y3; 429e3cd5995SDave May PetscReal b[2],A[2][2],inv[2][2],detJ,od; 430e3cd5995SDave May 431e3cd5995SDave May PetscFunctionBegin; 432a5f152d1SDave May x1 = PetscRealPart(coords[2*0+0]); 433a5f152d1SDave May x2 = PetscRealPart(coords[2*1+0]); 434a5f152d1SDave May x3 = PetscRealPart(coords[2*2+0]); 435e3cd5995SDave May 436a5f152d1SDave May y1 = PetscRealPart(coords[2*0+1]); 437a5f152d1SDave May y2 = PetscRealPart(coords[2*1+1]); 438a5f152d1SDave May y3 = PetscRealPart(coords[2*2+1]); 439e3cd5995SDave May 440e3cd5995SDave May b[0] = xp[0] - x1; 441e3cd5995SDave May b[1] = xp[1] - y1; 442e3cd5995SDave May 443e3cd5995SDave May A[0][0] = x2-x1; A[0][1] = x3-x1; 444e3cd5995SDave May A[1][0] = y2-y1; A[1][1] = y3-y1; 445e3cd5995SDave May 446e3cd5995SDave May detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0]; 447e3cd5995SDave May *dJ = PetscAbsReal(detJ); 448e3cd5995SDave May od = 1.0/detJ; 449e3cd5995SDave May 450e3cd5995SDave May inv[0][0] = A[1][1] * od; 451e3cd5995SDave May inv[0][1] = -A[0][1] * od; 452e3cd5995SDave May inv[1][0] = -A[1][0] * od; 453e3cd5995SDave May inv[1][1] = A[0][0] * od; 454e3cd5995SDave May 455e3cd5995SDave May xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1]; 456e3cd5995SDave May xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1]; 457e3cd5995SDave May PetscFunctionReturn(0); 458e3cd5995SDave May } 459e3cd5995SDave May 460e3cd5995SDave May PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field) 461e3cd5995SDave May { 462e3cd5995SDave May PetscErrorCode ierr; 463e3cd5995SDave May const PetscReal PLEX_C_EPS = 1.0e-8; 464e3cd5995SDave May Vec v_field_l,denom_l,coor_l,denom; 465e3cd5995SDave May PetscInt k,p,e,npoints; 466e3cd5995SDave May PetscInt *mpfield_cell; 467e3cd5995SDave May PetscReal *mpfield_coor; 468a5f152d1SDave May PetscReal xi_p[2]; 469a5f152d1SDave May PetscScalar Ni[3]; 470e3cd5995SDave May PetscSection coordSection; 471e3cd5995SDave May PetscScalar *elcoor = NULL; 472e3cd5995SDave May 473e3cd5995SDave May PetscFunctionBegin; 474e3cd5995SDave May ierr = VecZeroEntries(v_field);CHKERRQ(ierr); 475e3cd5995SDave May 476e3cd5995SDave May ierr = DMGetLocalVector(dm,&v_field_l);CHKERRQ(ierr); 477e3cd5995SDave May ierr = DMGetGlobalVector(dm,&denom);CHKERRQ(ierr); 478e3cd5995SDave May ierr = DMGetLocalVector(dm,&denom_l);CHKERRQ(ierr); 479e3cd5995SDave May ierr = VecZeroEntries(v_field_l);CHKERRQ(ierr); 480e3cd5995SDave May ierr = VecZeroEntries(denom);CHKERRQ(ierr); 481e3cd5995SDave May ierr = VecZeroEntries(denom_l);CHKERRQ(ierr); 482e3cd5995SDave May 483e3cd5995SDave May ierr = DMGetCoordinatesLocal(dm,&coor_l);CHKERRQ(ierr); 484e3cd5995SDave May ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 485e3cd5995SDave May 486e3cd5995SDave May ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr); 487e3cd5995SDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 488e3cd5995SDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 489e3cd5995SDave May 490e3cd5995SDave May for (p=0; p<npoints; p++) { 491a5f152d1SDave May PetscReal *coor_p,dJ; 492a5f152d1SDave May PetscScalar elfield[3]; 493e3cd5995SDave May PetscBool point_located; 494e3cd5995SDave May 495e3cd5995SDave May e = mpfield_cell[p]; 496e3cd5995SDave May coor_p = &mpfield_coor[2*p]; 497e3cd5995SDave May 498e3cd5995SDave May ierr = DMPlexVecGetClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr); 499e3cd5995SDave May 500e3cd5995SDave May /* 501e3cd5995SDave May while (!point_located && (failed_counter < 25)) { 502e3cd5995SDave May ierr = PointInTriangle(point, coords[0], coords[1], coords[2], &point_located);CHKERRQ(ierr); 503e3cd5995SDave May point.x = coor_p[0]; 504e3cd5995SDave May point.y = coor_p[1]; 505e3cd5995SDave May point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 506e3cd5995SDave May point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0); 507e3cd5995SDave May failed_counter++; 508e3cd5995SDave May } 509e3cd5995SDave May 510e3cd5995SDave May if (!point_located){ 511e3cd5995SDave 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); 512e3cd5995SDave May } 513e3cd5995SDave May 514e3cd5995SDave 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); 515e3cd5995SDave May else { 516e3cd5995SDave May ierr = _ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr); 517e3cd5995SDave May xi_p[0] = 0.5*(xi_p[0] + 1.0); 518e3cd5995SDave May xi_p[1] = 0.5*(xi_p[1] + 1.0); 519e3cd5995SDave May 520e3cd5995SDave 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]); 521e3cd5995SDave May 522e3cd5995SDave May } 523e3cd5995SDave May */ 524e3cd5995SDave May 525e3cd5995SDave May ierr = ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr); 526e3cd5995SDave May /* 527e3cd5995SDave 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]); 528e3cd5995SDave May */ 529e3cd5995SDave May /* 530e3cd5995SDave May point_located = PETSC_TRUE; 531e3cd5995SDave May if (xi_p[0] < 0.0) { 532e3cd5995SDave May if (xi_p[0] > -PLEX_C_EPS) { 533e3cd5995SDave May xi_p[0] = 0.0; 534e3cd5995SDave May } else { 535e3cd5995SDave May point_located = PETSC_FALSE; 536e3cd5995SDave May } 537e3cd5995SDave May } 538e3cd5995SDave May if (xi_p[1] < 0.0) { 539e3cd5995SDave May if (xi_p[1] > -PLEX_C_EPS) { 540e3cd5995SDave May xi_p[1] = 0.0; 541e3cd5995SDave May } else { 542e3cd5995SDave May point_located = PETSC_FALSE; 543e3cd5995SDave May } 544e3cd5995SDave May } 545e3cd5995SDave May if (xi_p[1] > (1.0-xi_p[0])) { 546e3cd5995SDave May if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { 547e3cd5995SDave May xi_p[1] = 1.0 - xi_p[0]; 548e3cd5995SDave May } else { 549e3cd5995SDave May point_located = PETSC_FALSE; 550e3cd5995SDave May } 551e3cd5995SDave May } 552e3cd5995SDave May if (!point_located){ 553e3cd5995SDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); 554e3cd5995SDave 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]); 555e3cd5995SDave May } 556e3cd5995SDave 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); 557e3cd5995SDave May */ 558e3cd5995SDave May 559e3cd5995SDave May Ni[0] = 1.0 - xi_p[0] - xi_p[1]; 560e3cd5995SDave May Ni[1] = xi_p[0]; 561e3cd5995SDave May Ni[2] = xi_p[1]; 562e3cd5995SDave May 563e3cd5995SDave May point_located = PETSC_TRUE; 564e3cd5995SDave May for (k=0; k<3; k++) { 565a5f152d1SDave May if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; 566a5f152d1SDave May if (PetscRealPart(Ni[k]) > (1.0+PLEX_C_EPS)) point_located = PETSC_FALSE; 567e3cd5995SDave May } 568e3cd5995SDave May if (!point_located){ 569d81390bbSDave May PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",(double)xi_p[0],(double)xi_p[1]); 570d81390bbSDave 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])); 571e3cd5995SDave May } 572d81390bbSDave 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); 573e3cd5995SDave May 574e3cd5995SDave May for (k=0; k<3; k++) { 575e3cd5995SDave May Ni[k] = Ni[k] * dJ; 576e3cd5995SDave May elfield[k] = Ni[k] * swarm_field[p]; 577e3cd5995SDave May } 578e3cd5995SDave May 579e3cd5995SDave May ierr = DMPlexVecRestoreClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr); 580e3cd5995SDave May 581e3cd5995SDave May ierr = DMPlexVecSetClosure(dm, NULL,v_field_l, e, elfield, ADD_VALUES);CHKERRQ(ierr); 582e3cd5995SDave May ierr = DMPlexVecSetClosure(dm, NULL,denom_l, e, Ni, ADD_VALUES);CHKERRQ(ierr); 583e3cd5995SDave May } 584e3cd5995SDave May 585e3cd5995SDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 586e3cd5995SDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 587e3cd5995SDave May 588e3cd5995SDave May ierr = DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 589e3cd5995SDave May ierr = DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 590e3cd5995SDave May ierr = DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 591e3cd5995SDave May ierr = DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 592e3cd5995SDave May 593e3cd5995SDave May ierr = VecPointwiseDivide(v_field,v_field,denom);CHKERRQ(ierr); 594e3cd5995SDave May 595e3cd5995SDave May ierr = DMRestoreLocalVector(dm,&v_field_l);CHKERRQ(ierr); 596e3cd5995SDave May ierr = DMRestoreLocalVector(dm,&denom_l);CHKERRQ(ierr); 597e3cd5995SDave May ierr = DMRestoreGlobalVector(dm,&denom);CHKERRQ(ierr); 598e3cd5995SDave May 599e3cd5995SDave May PetscFunctionReturn(0); 600e3cd5995SDave May } 601e3cd5995SDave May 602e3cd5995SDave May PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DataField dfield[],Vec vecs[]) 603e3cd5995SDave May { 604e3cd5995SDave May PetscErrorCode ierr; 605e3cd5995SDave May PetscInt f,dim; 606e3cd5995SDave May 607e3cd5995SDave May PetscFunctionBegin; 608e3cd5995SDave May ierr = DMGetDimension(swarm,&dim);CHKERRQ(ierr); 609e3cd5995SDave May switch (dim) { 610e3cd5995SDave May case 2: 611e3cd5995SDave May for (f=0; f<nfields; f++) { 612e3cd5995SDave May PetscReal *swarm_field; 613e3cd5995SDave May 614e3cd5995SDave May ierr = DataFieldGetEntries(dfield[f],(void**)&swarm_field);CHKERRQ(ierr); 615e3cd5995SDave May ierr = DMSwarmProjectField_ApproxP1_PLEX_2D(swarm,swarm_field,celldm,vecs[f]);CHKERRQ(ierr); 616e3cd5995SDave May } 617e3cd5995SDave May break; 618e3cd5995SDave May case 3: 619e3cd5995SDave May SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D"); 620e3cd5995SDave May break; 621e3cd5995SDave May default: 622e3cd5995SDave May break; 623e3cd5995SDave May } 624e3cd5995SDave May 625e3cd5995SDave May PetscFunctionReturn(0); 626e3cd5995SDave May } 627431382f2SDave May 628431382f2SDave May PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm,DM dmc,PetscInt npoints,PetscReal xi[]) 629431382f2SDave May { 630431382f2SDave May PetscBool is_simplex,is_tensorcell; 631431382f2SDave May PetscErrorCode ierr; 63292e40656SDave May PetscInt dim,nfaces,ps,pe,p,d,nbasis,pcnt,e,k,nel; 63392e40656SDave May PetscFE fe; 63492e40656SDave May PetscQuadrature quadrature; 63592e40656SDave May PetscReal *B; 63692e40656SDave May Vec coorlocal; 63792e40656SDave May PetscSection coordSection; 63892e40656SDave May PetscScalar *elcoor = NULL; 63992e40656SDave May PetscReal *swarm_coor; 64092e40656SDave May PetscInt *swarm_cellid; 641431382f2SDave May 64292e40656SDave May PetscFunctionBegin; 643431382f2SDave May ierr = DMGetDimension(dmc,&dim);CHKERRQ(ierr); 644431382f2SDave May 645431382f2SDave May is_simplex = PETSC_FALSE; 646431382f2SDave May is_tensorcell = PETSC_FALSE; 647431382f2SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 648431382f2SDave May ierr = DMPlexGetConeSize(dmc, ps, &nfaces);CHKERRQ(ierr); 649431382f2SDave May 650431382f2SDave May if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } 651431382f2SDave May 652431382f2SDave May switch (dim) { 653431382f2SDave May case 2: 654431382f2SDave May if (nfaces == 4) { is_tensorcell = PETSC_TRUE; } 655431382f2SDave May break; 656431382f2SDave May case 3: 657431382f2SDave May if (nfaces == 6) { is_tensorcell = PETSC_TRUE; } 658431382f2SDave May break; 659431382f2SDave May default: 660431382f2SDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for 2D, 3D"); 661431382f2SDave May break; 662431382f2SDave May } 663431382f2SDave May 66492e40656SDave May /* check points provided fail inside the reference cell */ 665431382f2SDave May if (is_simplex) { 66692e40656SDave May for (p=0; p<npoints; p++) { 66792e40656SDave May PetscReal sum; 66892e40656SDave May for (d=0; d<dim; d++) { 66992e40656SDave May if (xi[dim*p+d] < -1.0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain"); 67092e40656SDave May } 67192e40656SDave May sum = 0.0; 67292e40656SDave May for (d=0; d<dim; d++) { 67392e40656SDave May sum += xi[dim*p+d]; 67492e40656SDave May } 67592e40656SDave May if (sum > 0.0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain"); 67692e40656SDave May } 677431382f2SDave May } else if (is_tensorcell) { 67892e40656SDave May for (p=0; p<npoints; p++) { 67992e40656SDave May for (d=0; d<dim; d++) { 68092e40656SDave 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"); 68192e40656SDave May } 68292e40656SDave May } 683431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for d-simplex and d-tensorcell"); 684431382f2SDave May 68592e40656SDave May ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)dm),&quadrature);CHKERRQ(ierr); 68692e40656SDave May ierr = PetscQuadratureSetData(quadrature,dim,1,npoints,(const PetscReal*)xi,NULL);CHKERRQ(ierr); 687*2cf64d79SDave May ierr = private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe);CHKERRQ(ierr); 68892e40656SDave May ierr = PetscFESetQuadrature(fe,quadrature);CHKERRQ(ierr); 68992e40656SDave May ierr = PetscFEGetDimension(fe,&nbasis);CHKERRQ(ierr); 69092e40656SDave May ierr = PetscFEGetDefaultTabulation(fe, &B, NULL, NULL);CHKERRQ(ierr); 69192e40656SDave May 69292e40656SDave May /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ 69392e40656SDave May /* 0->cell, 1->edge, 2->vert */ 69492e40656SDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 69592e40656SDave May nel = pe - ps; 69692e40656SDave May 69792e40656SDave May ierr = DMSwarmSetLocalSizes(dm,npoints*nel,-1);CHKERRQ(ierr); 69892e40656SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 69992e40656SDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 70092e40656SDave May 70192e40656SDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 70292e40656SDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 70392e40656SDave May 70492e40656SDave May pcnt = 0; 70592e40656SDave May for (e=0; e<nel; e++) { 70692e40656SDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 70792e40656SDave May 70892e40656SDave May for (p=0; p<npoints; p++) { 70992e40656SDave May for (d=0; d<dim; d++) { 71092e40656SDave May swarm_coor[dim*pcnt+d] = 0.0; 71192e40656SDave May for (k=0; k<nbasis; k++) { 71292e40656SDave May swarm_coor[dim*pcnt+d] += B[p*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); 71392e40656SDave May } 71492e40656SDave May } 71592e40656SDave May swarm_cellid[pcnt] = e; 71692e40656SDave May pcnt++; 71792e40656SDave May } 71892e40656SDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor);CHKERRQ(ierr); 71992e40656SDave May } 72092e40656SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 72192e40656SDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 72292e40656SDave May 72392e40656SDave May ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr); 72492e40656SDave May ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 72592e40656SDave May 726431382f2SDave May PetscFunctionReturn(0); 727431382f2SDave May } 728