xref: /petsc/src/dm/impls/swarm/swarmpic_da.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
10e2ec84fSDave May #include <petscdm.h>
20e2ec84fSDave May #include <petscdmda.h>
30e2ec84fSDave May #include <petscdmswarm.h>
41f52046fSDave May #include <petsc/private/dmswarmimpl.h>
5279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
60e2ec84fSDave May 
70e2ec84fSDave May PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim,PetscInt np[],PetscInt *_npoints,PetscReal **_xi)
80e2ec84fSDave May {
90e2ec84fSDave May   PetscReal      *xi;
10f9368e01SSatish Balay   PetscInt       d,npoints=0,cnt;
110e2ec84fSDave May   PetscReal      ds[] = {0.0,0.0,0.0};
120e2ec84fSDave May   PetscInt       ii,jj,kk;
130e2ec84fSDave May 
140e2ec84fSDave May   PetscFunctionBegin;
150e2ec84fSDave May   switch (dim) {
160e2ec84fSDave May     case 1:
170e2ec84fSDave May       npoints = np[0];
180e2ec84fSDave May       break;
190e2ec84fSDave May     case 2:
200e2ec84fSDave May       npoints = np[0]*np[1];
210e2ec84fSDave May       break;
220e2ec84fSDave May     case 3:
230e2ec84fSDave May       npoints = np[0]*np[1]*np[2];
240e2ec84fSDave May       break;
250e2ec84fSDave May   }
260e2ec84fSDave May   for (d=0; d<dim; d++) {
270e2ec84fSDave May     ds[d] = 2.0 / ((PetscReal)np[d]);
280e2ec84fSDave May   }
290e2ec84fSDave May 
309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim*npoints,&xi));
310e2ec84fSDave May   switch (dim) {
320e2ec84fSDave May     case 1:
330e2ec84fSDave May       cnt = 0;
340e2ec84fSDave May       for (ii=0; ii<np[0]; ii++) {
350e2ec84fSDave May         xi[dim*cnt+0] = -1.0 + 0.5*ds[d] + ii*ds[0];
360e2ec84fSDave May         cnt++;
370e2ec84fSDave May       }
380e2ec84fSDave May       break;
390e2ec84fSDave May 
400e2ec84fSDave May     case 2:
410e2ec84fSDave May       cnt = 0;
420e2ec84fSDave May       for (jj=0; jj<np[1]; jj++) {
430e2ec84fSDave May         for (ii=0; ii<np[0]; ii++) {
440e2ec84fSDave May           xi[dim*cnt+0] = -1.0 + 0.5*ds[0] + ii*ds[0];
450e2ec84fSDave May           xi[dim*cnt+1] = -1.0 + 0.5*ds[1] + jj*ds[1];
460e2ec84fSDave May           cnt++;
470e2ec84fSDave May         }
480e2ec84fSDave May       }
490e2ec84fSDave May       break;
500e2ec84fSDave May 
510e2ec84fSDave May     case 3:
520e2ec84fSDave May       cnt = 0;
530e2ec84fSDave May       for (kk=0; kk<np[2]; kk++) {
540e2ec84fSDave May         for (jj=0; jj<np[1]; jj++) {
550e2ec84fSDave May           for (ii=0; ii<np[0]; ii++) {
560e2ec84fSDave May             xi[dim*cnt+0] = -1.0 + 0.5*ds[0] + ii*ds[0];
570e2ec84fSDave May             xi[dim*cnt+1] = -1.0 + 0.5*ds[1] + jj*ds[1];
580e2ec84fSDave May             xi[dim*cnt+2] = -1.0 + 0.5*ds[2] + kk*ds[2];
590e2ec84fSDave May             cnt++;
600e2ec84fSDave May           }
610e2ec84fSDave May         }
620e2ec84fSDave May       }
630e2ec84fSDave May       break;
640e2ec84fSDave May   }
650e2ec84fSDave May   *_npoints = npoints;
660e2ec84fSDave May   *_xi = xi;
670e2ec84fSDave May   PetscFunctionReturn(0);
680e2ec84fSDave May }
690e2ec84fSDave May 
700e2ec84fSDave May PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim,PetscInt np_1d,PetscInt *_npoints,PetscReal **_xi)
710e2ec84fSDave May {
720e2ec84fSDave May   PetscQuadrature quadrature;
730e2ec84fSDave May   const PetscReal *quadrature_xi;
740e2ec84fSDave May   PetscReal       *xi;
750e2ec84fSDave May   PetscInt        d,q,npoints_q;
760e2ec84fSDave May 
770e2ec84fSDave May   PetscFunctionBegin;
789566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussTensorQuadrature(dim,1,np_1d,-1.0,1.0,&quadrature));
799566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadrature,NULL,NULL,&npoints_q,&quadrature_xi,NULL));
809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim*npoints_q,&xi));
810e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
820e2ec84fSDave May     for (d=0; d<dim; d++) {
830e2ec84fSDave May       xi[dim*q+d] = quadrature_xi[dim*q+d];
840e2ec84fSDave May     }
850e2ec84fSDave May   }
869566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quadrature));
870e2ec84fSDave May   *_npoints = npoints_q;
880e2ec84fSDave May   *_xi = xi;
890e2ec84fSDave May   PetscFunctionReturn(0);
900e2ec84fSDave May }
910e2ec84fSDave May 
92d3393ee1SDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm,DM dmc,PetscInt npoints,DMSwarmPICLayoutType layout)
930e2ec84fSDave May {
940e2ec84fSDave May   PetscInt          dim,npoints_q;
950e2ec84fSDave May   PetscInt          nel,npe,e,q,k,d;
960e2ec84fSDave May   const PetscInt    *element_list;
970e2ec84fSDave May   PetscReal         **basis;
980e2ec84fSDave May   PetscReal         *xi;
990e2ec84fSDave May   Vec               coor;
1000e2ec84fSDave May   const PetscScalar *_coor;
1010e2ec84fSDave May   PetscReal         *elcoor;
1020e2ec84fSDave May   PetscReal         *swarm_coor;
1030e2ec84fSDave May   PetscInt          *swarm_cellid;
1040e2ec84fSDave May   PetscInt          pcnt;
1050e2ec84fSDave May 
1060e2ec84fSDave May   PetscFunctionBegin;
1079566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
1080e2ec84fSDave May   switch (layout) {
1090e2ec84fSDave May     case DMSWARMPIC_LAYOUT_REGULAR:
1100e2ec84fSDave May     {
1110e2ec84fSDave May       PetscInt np_dir[3];
1120e2ec84fSDave May       np_dir[0] = np_dir[1] = np_dir[2] = npoints;
1139566063dSJacob Faibussowitsch       PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim,np_dir,&npoints_q,&xi));
1140e2ec84fSDave May     }
1150e2ec84fSDave May       break;
1160e2ec84fSDave May     case DMSWARMPIC_LAYOUT_GAUSS:
1179566063dSJacob Faibussowitsch       PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(dim,npoints,&npoints_q,&xi));
1180e2ec84fSDave May       break;
119d3393ee1SDave May 
120d3393ee1SDave May     case DMSWARMPIC_LAYOUT_SUBDIVISION:
121d3393ee1SDave May     {
122d3393ee1SDave May       PetscInt s,nsub;
123d3393ee1SDave May       PetscInt np_dir[3];
124d3393ee1SDave May       nsub = npoints;
125d3393ee1SDave May       np_dir[0] = 1;
126d3393ee1SDave May       for (s=0; s<nsub; s++) {
127d3393ee1SDave May         np_dir[0] *= 2;
128d3393ee1SDave May       }
129d3393ee1SDave May       np_dir[1] = np_dir[0];
130d3393ee1SDave May       np_dir[2] = np_dir[0];
1319566063dSJacob Faibussowitsch       PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim,np_dir,&npoints_q,&xi));
132d3393ee1SDave May     }
133d3393ee1SDave May       break;
1340e2ec84fSDave May     default:
1350e2ec84fSDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"A valid DMSwarmPIC layout must be provided");
1360e2ec84fSDave May   }
1370e2ec84fSDave May 
1389566063dSJacob Faibussowitsch   PetscCall(DMDAGetElements(dmc,&nel,&npe,&element_list));
1399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim*npe,&elcoor));
1409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints_q,&basis));
1410e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
1429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npe,&basis[q]));
1430e2ec84fSDave May 
1440e2ec84fSDave May     switch (dim) {
1450e2ec84fSDave May       case 1:
1460e2ec84fSDave May         basis[q][0] = 0.5*(1.0 - xi[dim*q+0]);
1470e2ec84fSDave May         basis[q][1] = 0.5*(1.0 + xi[dim*q+0]);
1480e2ec84fSDave May         break;
1490e2ec84fSDave May       case 2:
1500e2ec84fSDave May         basis[q][0] = 0.25*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1]);
1510e2ec84fSDave May         basis[q][1] = 0.25*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1]);
1520e2ec84fSDave May         basis[q][2] = 0.25*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1]);
1530e2ec84fSDave May         basis[q][3] = 0.25*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1]);
1540e2ec84fSDave May         break;
1550e2ec84fSDave May 
1560e2ec84fSDave May       case 3:
1570e2ec84fSDave May         basis[q][0] = 0.125*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 - xi[dim*q+2]);
1580e2ec84fSDave May         basis[q][1] = 0.125*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 - xi[dim*q+2]);
1590e2ec84fSDave May         basis[q][2] = 0.125*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 - xi[dim*q+2]);
1600e2ec84fSDave May         basis[q][3] = 0.125*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 - xi[dim*q+2]);
1610e2ec84fSDave May         basis[q][4] = 0.125*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 + xi[dim*q+2]);
1620e2ec84fSDave May         basis[q][5] = 0.125*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 + xi[dim*q+2]);
1630e2ec84fSDave May         basis[q][6] = 0.125*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 + xi[dim*q+2]);
1640e2ec84fSDave May         basis[q][7] = 0.125*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 + xi[dim*q+2]);
1650e2ec84fSDave May         break;
1660e2ec84fSDave May     }
1670e2ec84fSDave May   }
1680e2ec84fSDave May 
1699566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm,npoints_q*nel,-1));
1709566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
1719566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
1720e2ec84fSDave May 
1739566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dmc,&coor));
1749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coor,&_coor));
1750e2ec84fSDave May   pcnt = 0;
1760e2ec84fSDave May   for (e=0; e<nel; e++) {
1770e2ec84fSDave May     const PetscInt *element = &element_list[npe*e];
1780e2ec84fSDave May 
1790e2ec84fSDave May     for (k=0; k<npe; k++) {
1800e2ec84fSDave May       for (d=0; d<dim; d++) {
1810ca99263SDave May         elcoor[dim*k+d] = PetscRealPart(_coor[ dim*element[k] + d ]);
1820e2ec84fSDave May       }
1830e2ec84fSDave May     }
1840e2ec84fSDave May 
1850e2ec84fSDave May     for (q=0; q<npoints_q; q++) {
1860e2ec84fSDave May       for (d=0; d<dim; d++) {
1870e2ec84fSDave May         swarm_coor[dim*pcnt+d] = 0.0;
1880e2ec84fSDave May       }
1890e2ec84fSDave May       for (k=0; k<npe; k++) {
1900e2ec84fSDave May         for (d=0; d<dim; d++) {
1910e2ec84fSDave May           swarm_coor[dim*pcnt+d] += basis[q][k] * elcoor[dim*k+d];
1920e2ec84fSDave May         }
1930e2ec84fSDave May       }
1940e2ec84fSDave May       swarm_cellid[pcnt] = e;
1950e2ec84fSDave May       pcnt++;
1960e2ec84fSDave May     }
1970e2ec84fSDave May   }
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coor,&_coor));
1999566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid));
2009566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor));
2019566063dSJacob Faibussowitsch   PetscCall(DMDARestoreElements(dmc,&nel,&npe,&element_list));
2020e2ec84fSDave May 
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(xi));
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree(elcoor));
2050e2ec84fSDave May   for (q=0; q<npoints_q; q++) {
2069566063dSJacob Faibussowitsch     PetscCall(PetscFree(basis[q]));
2070e2ec84fSDave May   }
2089566063dSJacob Faibussowitsch   PetscCall(PetscFree(basis));
2090e2ec84fSDave May   PetscFunctionReturn(0);
2100e2ec84fSDave May }
2110e2ec84fSDave May 
2120e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param)
2130e2ec84fSDave May {
2140e2ec84fSDave May   DMDAElementType etype;
2150e2ec84fSDave May   PetscInt        dim;
2160e2ec84fSDave May 
2170e2ec84fSDave May   PetscFunctionBegin;
2189566063dSJacob Faibussowitsch   PetscCall(DMDAGetElementType(celldm,&etype));
2199566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(celldm,&dim));
2200e2ec84fSDave May   switch (etype) {
2210e2ec84fSDave May     case DMDA_ELEMENT_P1:
2220e2ec84fSDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DA support is not currently available for DMDA_ELEMENT_P1");
2230e2ec84fSDave May     case DMDA_ELEMENT_Q1:
224*08401ef6SPierre Jolivet       PetscCheck(dim != 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Support only available for dim = 2, 3");
2259566063dSJacob Faibussowitsch       PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA_Q1(dm,celldm,layout_param,layout));
2260e2ec84fSDave May       break;
2270e2ec84fSDave May   }
2280e2ec84fSDave May   PetscFunctionReturn(0);
2290e2ec84fSDave May }
2301f52046fSDave May 
2311f52046fSDave May PetscErrorCode DMSwarmProjectField_ApproxQ1_DA_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field)
2321f52046fSDave May {
2331f52046fSDave May   Vec               v_field_l,denom_l,coor_l,denom;
234a5f152d1SDave May   PetscScalar       *_field_l,*_denom_l;
2351f52046fSDave May   PetscInt          k,p,e,npoints,nel,npe;
2361f52046fSDave May   PetscInt          *mpfield_cell;
2371f52046fSDave May   PetscReal         *mpfield_coor;
2381f52046fSDave May   const PetscInt    *element_list;
2391f52046fSDave May   const PetscInt    *element;
240a5f152d1SDave May   PetscScalar       xi_p[2],Ni[4];
241a5f152d1SDave May   const PetscScalar *_coor;
2421f52046fSDave May 
2431f52046fSDave May   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(v_field));
2451f52046fSDave May 
2469566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&v_field_l));
2479566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm,&denom));
2489566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&denom_l));
2499566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(v_field_l));
2509566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(denom));
2519566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(denom_l));
2521f52046fSDave May 
2539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v_field_l,&_field_l));
2549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(denom_l,&_denom_l));
2551f52046fSDave May 
2569566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm,&coor_l));
2579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coor_l,&_coor));
2581f52046fSDave May 
2599566063dSJacob Faibussowitsch   PetscCall(DMDAGetElements(dm,&nel,&npe,&element_list));
2609566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(swarm,&npoints));
2619566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor));
2629566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell));
2631f52046fSDave May 
2641f52046fSDave May   for (p=0; p<npoints; p++) {
2651f52046fSDave May     PetscReal         *coor_p;
266a5f152d1SDave May     const PetscScalar *x0;
267a5f152d1SDave May     const PetscScalar *x2;
268a5f152d1SDave May     PetscScalar       dx[2];
2691f52046fSDave May 
2701f52046fSDave May     e = mpfield_cell[p];
2711f52046fSDave May     coor_p = &mpfield_coor[2*p];
2721f52046fSDave May     element = &element_list[npe*e];
2731f52046fSDave May 
2741f52046fSDave May     /* compute local coordinates: (xp-x0)/dx = (xip+1)/2 */
2751f52046fSDave May     x0 = &_coor[2*element[0]];
2761f52046fSDave May     x2 = &_coor[2*element[2]];
2771f52046fSDave May 
2781f52046fSDave May     dx[0] = x2[0] - x0[0];
2791f52046fSDave May     dx[1] = x2[1] - x0[1];
2801f52046fSDave May 
2811f52046fSDave May     xi_p[0] = 2.0 * (coor_p[0] - x0[0])/dx[0] - 1.0;
2821f52046fSDave May     xi_p[1] = 2.0 * (coor_p[1] - x0[1])/dx[1] - 1.0;
2831f52046fSDave May 
2841f52046fSDave May     /* evaluate basis functions */
2851f52046fSDave May     Ni[0] = 0.25*(1.0 - xi_p[0])*(1.0 - xi_p[1]);
2861f52046fSDave May     Ni[1] = 0.25*(1.0 + xi_p[0])*(1.0 - xi_p[1]);
2871f52046fSDave May     Ni[2] = 0.25*(1.0 + xi_p[0])*(1.0 + xi_p[1]);
2881f52046fSDave May     Ni[3] = 0.25*(1.0 - xi_p[0])*(1.0 + xi_p[1]);
2891f52046fSDave May 
2901f52046fSDave May     for (k=0; k<npe; k++) {
2911f52046fSDave May       _field_l[ element[k] ] += Ni[k] * swarm_field[p];
2921f52046fSDave May       _denom_l[ element[k] ] += Ni[k];
2931f52046fSDave May     }
2941f52046fSDave May   }
2951f52046fSDave May 
2969566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell));
2979566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor));
2989566063dSJacob Faibussowitsch   PetscCall(DMDARestoreElements(dm,&nel,&npe,&element_list));
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coor_l,&_coor));
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v_field_l,&_field_l));
3019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(denom_l,&_denom_l));
3021f52046fSDave May 
3039566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field));
3049566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field));
3059566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom));
3069566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom));
3071f52046fSDave May 
3089566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(v_field,v_field,denom));
3091f52046fSDave May 
3109566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&v_field_l));
3119566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&denom_l));
3129566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm,&denom));
3131f52046fSDave May   PetscFunctionReturn(0);
3141f52046fSDave May }
3151f52046fSDave May 
31677048351SPatrick Sanan PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[])
3171f52046fSDave May {
3181f52046fSDave May   PetscInt        f,dim;
3191f52046fSDave May   DMDAElementType etype;
3201f52046fSDave May 
3211f52046fSDave May   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(DMDAGetElementType(celldm,&etype));
323*08401ef6SPierre Jolivet   PetscCheck(etype != DMDA_ELEMENT_P1,PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"Only Q1 DMDA supported");
3241f52046fSDave May 
3259566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(swarm,&dim));
3261f52046fSDave May   switch (dim) {
3271f52046fSDave May     case 2:
3281f52046fSDave May       for (f=0; f<nfields; f++) {
3291f52046fSDave May         PetscReal *swarm_field;
3301f52046fSDave May 
3319566063dSJacob Faibussowitsch         PetscCall(DMSwarmDataFieldGetEntries(dfield[f],(void**)&swarm_field));
3329566063dSJacob Faibussowitsch         PetscCall(DMSwarmProjectField_ApproxQ1_DA_2D(swarm,swarm_field,celldm,vecs[f]));
3331f52046fSDave May       }
3341f52046fSDave May       break;
3351f52046fSDave May     case 3:
3361f52046fSDave May       SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D");
3371f52046fSDave May     default:
3381f52046fSDave May       break;
3391f52046fSDave May   }
3401f52046fSDave May   PetscFunctionReturn(0);
3411f52046fSDave May }
342