xref: /petsc/src/dm/impls/swarm/swarmpic_da.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
7*d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim, PetscInt np[], PetscInt *_npoints, PetscReal **_xi)
8*d71ae5a4SJacob Faibussowitsch {
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) {
16*d71ae5a4SJacob Faibussowitsch   case 1:
17*d71ae5a4SJacob Faibussowitsch     npoints = np[0];
18*d71ae5a4SJacob Faibussowitsch     break;
19*d71ae5a4SJacob Faibussowitsch   case 2:
20*d71ae5a4SJacob Faibussowitsch     npoints = np[0] * np[1];
21*d71ae5a4SJacob Faibussowitsch     break;
22*d71ae5a4SJacob Faibussowitsch   case 3:
23*d71ae5a4SJacob Faibussowitsch     npoints = np[0] * np[1] * np[2];
24*d71ae5a4SJacob Faibussowitsch     break;
250e2ec84fSDave May   }
26ad540459SPierre Jolivet   for (d = 0; d < dim; d++) ds[d] = 2.0 / ((PetscReal)np[d]);
270e2ec84fSDave May 
289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * npoints, &xi));
290e2ec84fSDave May   switch (dim) {
300e2ec84fSDave May   case 1:
310e2ec84fSDave May     cnt = 0;
320e2ec84fSDave May     for (ii = 0; ii < np[0]; ii++) {
330e2ec84fSDave May       xi[dim * cnt + 0] = -1.0 + 0.5 * ds[d] + ii * ds[0];
340e2ec84fSDave May       cnt++;
350e2ec84fSDave May     }
360e2ec84fSDave May     break;
370e2ec84fSDave May 
380e2ec84fSDave May   case 2:
390e2ec84fSDave May     cnt = 0;
400e2ec84fSDave May     for (jj = 0; jj < np[1]; jj++) {
410e2ec84fSDave May       for (ii = 0; ii < np[0]; ii++) {
420e2ec84fSDave May         xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0];
430e2ec84fSDave May         xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1];
440e2ec84fSDave May         cnt++;
450e2ec84fSDave May       }
460e2ec84fSDave May     }
470e2ec84fSDave May     break;
480e2ec84fSDave May 
490e2ec84fSDave May   case 3:
500e2ec84fSDave May     cnt = 0;
510e2ec84fSDave May     for (kk = 0; kk < np[2]; kk++) {
520e2ec84fSDave May       for (jj = 0; jj < np[1]; jj++) {
530e2ec84fSDave May         for (ii = 0; ii < np[0]; ii++) {
540e2ec84fSDave May           xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0];
550e2ec84fSDave May           xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1];
560e2ec84fSDave May           xi[dim * cnt + 2] = -1.0 + 0.5 * ds[2] + kk * ds[2];
570e2ec84fSDave May           cnt++;
580e2ec84fSDave May         }
590e2ec84fSDave May       }
600e2ec84fSDave May     }
610e2ec84fSDave May     break;
620e2ec84fSDave May   }
630e2ec84fSDave May   *_npoints = npoints;
640e2ec84fSDave May   *_xi      = xi;
650e2ec84fSDave May   PetscFunctionReturn(0);
660e2ec84fSDave May }
670e2ec84fSDave May 
68*d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim, PetscInt np_1d, PetscInt *_npoints, PetscReal **_xi)
69*d71ae5a4SJacob Faibussowitsch {
700e2ec84fSDave May   PetscQuadrature  quadrature;
710e2ec84fSDave May   const PetscReal *quadrature_xi;
720e2ec84fSDave May   PetscReal       *xi;
730e2ec84fSDave May   PetscInt         d, q, npoints_q;
740e2ec84fSDave May 
750e2ec84fSDave May   PetscFunctionBegin;
769566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussTensorQuadrature(dim, 1, np_1d, -1.0, 1.0, &quadrature));
779566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &quadrature_xi, NULL));
789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * npoints_q, &xi));
790e2ec84fSDave May   for (q = 0; q < npoints_q; q++) {
80ad540459SPierre Jolivet     for (d = 0; d < dim; d++) xi[dim * q + d] = quadrature_xi[dim * q + d];
810e2ec84fSDave May   }
829566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quadrature));
830e2ec84fSDave May   *_npoints = npoints_q;
840e2ec84fSDave May   *_xi      = xi;
850e2ec84fSDave May   PetscFunctionReturn(0);
860e2ec84fSDave May }
870e2ec84fSDave May 
88*d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm, DM dmc, PetscInt npoints, DMSwarmPICLayoutType layout)
89*d71ae5a4SJacob Faibussowitsch {
900e2ec84fSDave May   PetscInt           dim, npoints_q;
910e2ec84fSDave May   PetscInt           nel, npe, e, q, k, d;
920e2ec84fSDave May   const PetscInt    *element_list;
930e2ec84fSDave May   PetscReal        **basis;
940e2ec84fSDave May   PetscReal         *xi;
950e2ec84fSDave May   Vec                coor;
960e2ec84fSDave May   const PetscScalar *_coor;
970e2ec84fSDave May   PetscReal         *elcoor;
980e2ec84fSDave May   PetscReal         *swarm_coor;
990e2ec84fSDave May   PetscInt          *swarm_cellid;
1000e2ec84fSDave May   PetscInt           pcnt;
1010e2ec84fSDave May 
1020e2ec84fSDave May   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1040e2ec84fSDave May   switch (layout) {
1059371c9d4SSatish Balay   case DMSWARMPIC_LAYOUT_REGULAR: {
1060e2ec84fSDave May     PetscInt np_dir[3];
1070e2ec84fSDave May     np_dir[0] = np_dir[1] = np_dir[2] = npoints;
1089566063dSJacob Faibussowitsch     PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi));
1099371c9d4SSatish Balay   } break;
110*d71ae5a4SJacob Faibussowitsch   case DMSWARMPIC_LAYOUT_GAUSS:
111*d71ae5a4SJacob Faibussowitsch     PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(dim, npoints, &npoints_q, &xi));
112*d71ae5a4SJacob Faibussowitsch     break;
113d3393ee1SDave May 
1149371c9d4SSatish Balay   case DMSWARMPIC_LAYOUT_SUBDIVISION: {
115d3393ee1SDave May     PetscInt s, nsub;
116d3393ee1SDave May     PetscInt np_dir[3];
117d3393ee1SDave May     nsub      = npoints;
118d3393ee1SDave May     np_dir[0] = 1;
119ad540459SPierre Jolivet     for (s = 0; s < nsub; s++) np_dir[0] *= 2;
120d3393ee1SDave May     np_dir[1] = np_dir[0];
121d3393ee1SDave May     np_dir[2] = np_dir[0];
1229566063dSJacob Faibussowitsch     PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi));
1239371c9d4SSatish Balay   } break;
124*d71ae5a4SJacob Faibussowitsch   default:
125*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "A valid DMSwarmPIC layout must be provided");
1260e2ec84fSDave May   }
1270e2ec84fSDave May 
1289566063dSJacob Faibussowitsch   PetscCall(DMDAGetElements(dmc, &nel, &npe, &element_list));
1299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * npe, &elcoor));
1309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints_q, &basis));
1310e2ec84fSDave May   for (q = 0; q < npoints_q; q++) {
1329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(npe, &basis[q]));
1330e2ec84fSDave May 
1340e2ec84fSDave May     switch (dim) {
1350e2ec84fSDave May     case 1:
1360e2ec84fSDave May       basis[q][0] = 0.5 * (1.0 - xi[dim * q + 0]);
1370e2ec84fSDave May       basis[q][1] = 0.5 * (1.0 + xi[dim * q + 0]);
1380e2ec84fSDave May       break;
1390e2ec84fSDave May     case 2:
1400e2ec84fSDave May       basis[q][0] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]);
1410e2ec84fSDave May       basis[q][1] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]);
1420e2ec84fSDave May       basis[q][2] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]);
1430e2ec84fSDave May       basis[q][3] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]);
1440e2ec84fSDave May       break;
1450e2ec84fSDave May 
1460e2ec84fSDave May     case 3:
1470e2ec84fSDave May       basis[q][0] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1480e2ec84fSDave May       basis[q][1] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1490e2ec84fSDave May       basis[q][2] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1500e2ec84fSDave May       basis[q][3] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1510e2ec84fSDave May       basis[q][4] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1520e2ec84fSDave May       basis[q][5] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1530e2ec84fSDave May       basis[q][6] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1540e2ec84fSDave May       basis[q][7] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1550e2ec84fSDave May       break;
1560e2ec84fSDave May     }
1570e2ec84fSDave May   }
1580e2ec84fSDave May 
1599566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
1609566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
1619566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
1620e2ec84fSDave May 
1639566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dmc, &coor));
1649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coor, &_coor));
1650e2ec84fSDave May   pcnt = 0;
1660e2ec84fSDave May   for (e = 0; e < nel; e++) {
1670e2ec84fSDave May     const PetscInt *element = &element_list[npe * e];
1680e2ec84fSDave May 
1690e2ec84fSDave May     for (k = 0; k < npe; k++) {
170ad540459SPierre Jolivet       for (d = 0; d < dim; d++) elcoor[dim * k + d] = PetscRealPart(_coor[dim * element[k] + d]);
1710e2ec84fSDave May     }
1720e2ec84fSDave May 
1730e2ec84fSDave May     for (q = 0; q < npoints_q; q++) {
174ad540459SPierre Jolivet       for (d = 0; d < dim; d++) swarm_coor[dim * pcnt + d] = 0.0;
1750e2ec84fSDave May       for (k = 0; k < npe; k++) {
176ad540459SPierre Jolivet         for (d = 0; d < dim; d++) swarm_coor[dim * pcnt + d] += basis[q][k] * elcoor[dim * k + d];
1770e2ec84fSDave May       }
1780e2ec84fSDave May       swarm_cellid[pcnt] = e;
1790e2ec84fSDave May       pcnt++;
1800e2ec84fSDave May     }
1810e2ec84fSDave May   }
1829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coor, &_coor));
1839566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
1849566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor));
1859566063dSJacob Faibussowitsch   PetscCall(DMDARestoreElements(dmc, &nel, &npe, &element_list));
1860e2ec84fSDave May 
1879566063dSJacob Faibussowitsch   PetscCall(PetscFree(xi));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFree(elcoor));
18948a46eb9SPierre Jolivet   for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q]));
1909566063dSJacob Faibussowitsch   PetscCall(PetscFree(basis));
1910e2ec84fSDave May   PetscFunctionReturn(0);
1920e2ec84fSDave May }
1930e2ec84fSDave May 
194*d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param)
195*d71ae5a4SJacob Faibussowitsch {
1960e2ec84fSDave May   DMDAElementType etype;
1970e2ec84fSDave May   PetscInt        dim;
1980e2ec84fSDave May 
1990e2ec84fSDave May   PetscFunctionBegin;
2009566063dSJacob Faibussowitsch   PetscCall(DMDAGetElementType(celldm, &etype));
2019566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(celldm, &dim));
2020e2ec84fSDave May   switch (etype) {
203*d71ae5a4SJacob Faibussowitsch   case DMDA_ELEMENT_P1:
204*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DA support is not currently available for DMDA_ELEMENT_P1");
2050e2ec84fSDave May   case DMDA_ELEMENT_Q1:
20608401ef6SPierre Jolivet     PetscCheck(dim != 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support only available for dim = 2, 3");
2079566063dSJacob Faibussowitsch     PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA_Q1(dm, celldm, layout_param, layout));
2080e2ec84fSDave May     break;
2090e2ec84fSDave May   }
2100e2ec84fSDave May   PetscFunctionReturn(0);
2110e2ec84fSDave May }
2121f52046fSDave May 
213*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmProjectField_ApproxQ1_DA_2D(DM swarm, PetscReal *swarm_field, DM dm, Vec v_field)
214*d71ae5a4SJacob Faibussowitsch {
2151f52046fSDave May   Vec                v_field_l, denom_l, coor_l, denom;
216a5f152d1SDave May   PetscScalar       *_field_l, *_denom_l;
2171f52046fSDave May   PetscInt           k, p, e, npoints, nel, npe;
2181f52046fSDave May   PetscInt          *mpfield_cell;
2191f52046fSDave May   PetscReal         *mpfield_coor;
2201f52046fSDave May   const PetscInt    *element_list;
2211f52046fSDave May   const PetscInt    *element;
222a5f152d1SDave May   PetscScalar        xi_p[2], Ni[4];
223a5f152d1SDave May   const PetscScalar *_coor;
2241f52046fSDave May 
2251f52046fSDave May   PetscFunctionBegin;
2269566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(v_field));
2271f52046fSDave May 
2289566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &v_field_l));
2299566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &denom));
2309566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &denom_l));
2319566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(v_field_l));
2329566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(denom));
2339566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(denom_l));
2341f52046fSDave May 
2359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v_field_l, &_field_l));
2369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(denom_l, &_denom_l));
2371f52046fSDave May 
2389566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coor_l));
2399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coor_l, &_coor));
2401f52046fSDave May 
2419566063dSJacob Faibussowitsch   PetscCall(DMDAGetElements(dm, &nel, &npe, &element_list));
2429566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(swarm, &npoints));
2439566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor));
2449566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell));
2451f52046fSDave May 
2461f52046fSDave May   for (p = 0; p < npoints; p++) {
2471f52046fSDave May     PetscReal         *coor_p;
248a5f152d1SDave May     const PetscScalar *x0;
249a5f152d1SDave May     const PetscScalar *x2;
250a5f152d1SDave May     PetscScalar        dx[2];
2511f52046fSDave May 
2521f52046fSDave May     e       = mpfield_cell[p];
2531f52046fSDave May     coor_p  = &mpfield_coor[2 * p];
2541f52046fSDave May     element = &element_list[npe * e];
2551f52046fSDave May 
2561f52046fSDave May     /* compute local coordinates: (xp-x0)/dx = (xip+1)/2 */
2571f52046fSDave May     x0 = &_coor[2 * element[0]];
2581f52046fSDave May     x2 = &_coor[2 * element[2]];
2591f52046fSDave May 
2601f52046fSDave May     dx[0] = x2[0] - x0[0];
2611f52046fSDave May     dx[1] = x2[1] - x0[1];
2621f52046fSDave May 
2631f52046fSDave May     xi_p[0] = 2.0 * (coor_p[0] - x0[0]) / dx[0] - 1.0;
2641f52046fSDave May     xi_p[1] = 2.0 * (coor_p[1] - x0[1]) / dx[1] - 1.0;
2651f52046fSDave May 
2661f52046fSDave May     /* evaluate basis functions */
2671f52046fSDave May     Ni[0] = 0.25 * (1.0 - xi_p[0]) * (1.0 - xi_p[1]);
2681f52046fSDave May     Ni[1] = 0.25 * (1.0 + xi_p[0]) * (1.0 - xi_p[1]);
2691f52046fSDave May     Ni[2] = 0.25 * (1.0 + xi_p[0]) * (1.0 + xi_p[1]);
2701f52046fSDave May     Ni[3] = 0.25 * (1.0 - xi_p[0]) * (1.0 + xi_p[1]);
2711f52046fSDave May 
2721f52046fSDave May     for (k = 0; k < npe; k++) {
2731f52046fSDave May       _field_l[element[k]] += Ni[k] * swarm_field[p];
2741f52046fSDave May       _denom_l[element[k]] += Ni[k];
2751f52046fSDave May     }
2761f52046fSDave May   }
2771f52046fSDave May 
2789566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell));
2799566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor));
2809566063dSJacob Faibussowitsch   PetscCall(DMDARestoreElements(dm, &nel, &npe, &element_list));
2819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coor_l, &_coor));
2829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v_field_l, &_field_l));
2839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(denom_l, &_denom_l));
2841f52046fSDave May 
2859566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, v_field_l, ADD_VALUES, v_field));
2869566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, v_field_l, ADD_VALUES, v_field));
2879566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(dm, denom_l, ADD_VALUES, denom));
2889566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(dm, denom_l, ADD_VALUES, denom));
2891f52046fSDave May 
2909566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(v_field, v_field, denom));
2911f52046fSDave May 
2929566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &v_field_l));
2939566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &denom_l));
2949566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &denom));
2951f52046fSDave May   PetscFunctionReturn(0);
2961f52046fSDave May }
2971f52046fSDave May 
298*d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[])
299*d71ae5a4SJacob Faibussowitsch {
3001f52046fSDave May   PetscInt        f, dim;
3011f52046fSDave May   DMDAElementType etype;
3021f52046fSDave May 
3031f52046fSDave May   PetscFunctionBegin;
3049566063dSJacob Faibussowitsch   PetscCall(DMDAGetElementType(celldm, &etype));
30508401ef6SPierre Jolivet   PetscCheck(etype != DMDA_ELEMENT_P1, PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "Only Q1 DMDA supported");
3061f52046fSDave May 
3079566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(swarm, &dim));
3081f52046fSDave May   switch (dim) {
3091f52046fSDave May   case 2:
3101f52046fSDave May     for (f = 0; f < nfields; f++) {
3111f52046fSDave May       PetscReal *swarm_field;
3121f52046fSDave May 
3139566063dSJacob Faibussowitsch       PetscCall(DMSwarmDataFieldGetEntries(dfield[f], (void **)&swarm_field));
3149566063dSJacob Faibussowitsch       PetscCall(DMSwarmProjectField_ApproxQ1_DA_2D(swarm, swarm_field, celldm, vecs[f]));
3151f52046fSDave May     }
3161f52046fSDave May     break;
317*d71ae5a4SJacob Faibussowitsch   case 3:
318*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "No support for 3D");
319*d71ae5a4SJacob Faibussowitsch   default:
320*d71ae5a4SJacob Faibussowitsch     break;
3211f52046fSDave May   }
3221f52046fSDave May   PetscFunctionReturn(0);
3231f52046fSDave May }
324