xref: /petsc/src/dm/impls/swarm/swarmpic_plex.c (revision b41ce5d507ea9a58bfa83cf403107a702e77a67d)
1 #include <petscdm.h>
2 #include <petscdmplex.h>
3 #include <petscdmswarm.h>
4 #include "data_bucket.h"
5 
6 PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[3],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np)
7 {
8   PetscReal v12[2],v23[2],v31[2];
9   PetscInt i;
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   if (depth == max) {
14     PetscReal cx[2];
15 
16     cx[0] = (v1[0] + v2[0] + v3[0])/3.0;
17     cx[1] = (v1[1] + v2[1] + v3[1])/3.0;
18 
19     xi[2*(*np)+0] = cx[0];
20     xi[2*(*np)+1] = cx[1];
21     (*np)++;
22     PetscFunctionReturn(0);
23   }
24 
25   /* calculate midpoints of each side */
26   for (i = 0; i < 2; i++) {
27     v12[i] = (v1[i]+v2[i])/2.0;
28     v23[i] = (v2[i]+v3[i])/2.0;
29     v31[i] = (v3[i]+v1[i])/2.0;
30   }
31 
32   /* recursively subdivide new triangles */
33   ierr = subdivide_triangle(v1,v12,v31,depth+1,max,xi,np);CHKERRQ(ierr);
34   ierr = subdivide_triangle(v2,v23,v12,depth+1,max,xi,np);CHKERRQ(ierr);
35   ierr = subdivide_triangle(v3,v31,v23,depth+1,max,xi,np);CHKERRQ(ierr);
36   ierr = subdivide_triangle(v12,v23,v31,depth+1,max,xi,np);CHKERRQ(ierr);
37   PetscFunctionReturn(0);
38 }
39 
40 
41 PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub)
42 {
43   PetscErrorCode ierr;
44   const PetscInt dim = 2;
45   PetscInt q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth;
46   PetscReal *xi;
47   PetscReal **basis;
48   Vec coorlocal;
49   PetscSection coordSection;
50   PetscScalar *elcoor = NULL;
51   PetscReal *swarm_coor;
52   PetscInt *swarm_cellid;
53   PetscReal v1[2],v2[2],v3[2];
54 
55   PetscFunctionBegin;
56 
57   npoints_q = 1;
58   for (d=0; d<nsub; d++) { npoints_q *= 4; }
59   ierr = PetscMalloc1(dim*npoints_q,&xi);CHKERRQ(ierr);
60 
61   v1[0] = 0.0;  v1[1] = 0.0;
62   v2[0] = 1.0;  v2[1] = 0.0;
63   v3[0] = 0.0;  v3[1] = 1.0;
64   depth = 0;
65   pcnt = 0;
66   ierr = subdivide_triangle(v1,v2,v3,depth,nsub,xi,&pcnt);CHKERRQ(ierr);
67 
68   npe = 3; /* nodes per element (triangle) */
69   ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr);
70   for (q=0; q<npoints_q; q++) {
71     ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr);
72 
73     basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1];
74     basis[q][1] = xi[dim*q+0];
75     basis[q][2] = xi[dim*q+1];
76   }
77 
78   /* 0->cell, 1->edge, 2->vert */
79   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
80   nel = pe - ps;
81 
82   ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr);
83   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
84   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
85 
86   ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr);
87   ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr);
88 
89   pcnt = 0;
90   for (e=0; e<nel; e++) {
91     ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
92 
93     for (q=0; q<npoints_q; q++) {
94       for (d=0; d<dim; d++) {
95         swarm_coor[dim*pcnt+d] = 0.0;
96         for (k=0; k<npe; k++) {
97           swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]);
98         }
99       }
100       swarm_cellid[pcnt] = e;
101       pcnt++;
102     }
103     ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
104   }
105   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
106   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
107 
108   ierr = PetscFree(xi);CHKERRQ(ierr);
109   for (q=0; q<npoints_q; q++) {
110     ierr = PetscFree(basis[q]);CHKERRQ(ierr);
111   }
112   ierr = PetscFree(basis);CHKERRQ(ierr);
113 
114   PetscFunctionReturn(0);
115 }
116 
117 PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints)
118 {
119   PetscErrorCode ierr;
120   const PetscInt dim = 2;
121   PetscInt ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k;
122   PetscReal *xi,ds,ds2;
123   PetscReal **basis;
124   Vec coorlocal;
125   PetscSection coordSection;
126   PetscScalar *elcoor = NULL;
127   PetscReal *swarm_coor;
128   PetscInt *swarm_cellid;
129 
130   PetscFunctionBegin;
131   ierr = PetscMalloc1(dim*npoints*npoints,&xi);CHKERRQ(ierr);
132   pcnt = 0;
133   ds = 1.0/((PetscReal)(npoints-1));
134   ds2 = 1.0/((PetscReal)(npoints));
135   for (jj = 0; jj<npoints; jj++) {
136     for (ii=0; ii<npoints-jj; ii++) {
137       xi[dim*pcnt+0] = ii * ds;
138       xi[dim*pcnt+1] = jj * ds;
139 
140       xi[dim*pcnt+0] *= (1.0 - 1.2*ds2);
141       xi[dim*pcnt+1] *= (1.0 - 1.2*ds2);
142 
143       xi[dim*pcnt+0] += 0.35*ds2;
144       xi[dim*pcnt+1] += 0.35*ds2;
145 
146       pcnt++;
147     }
148   }
149   npoints_q = pcnt;
150 
151   npe = 3; /* nodes per element (triangle) */
152   ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr);
153   for (q=0; q<npoints_q; q++) {
154     ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr);
155 
156     basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1];
157     basis[q][1] = xi[dim*q+0];
158     basis[q][2] = xi[dim*q+1];
159   }
160 
161   /* 0->cell, 1->edge, 2->vert */
162   ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr);
163   nel = pe - ps;
164 
165   ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr);
166   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
167   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
168 
169   ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr);
170   ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr);
171 
172   pcnt = 0;
173   for (e=0; e<nel; e++) {
174     ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
175 
176     for (q=0; q<npoints_q; q++) {
177       for (d=0; d<dim; d++) {
178         swarm_coor[dim*pcnt+d] = 0.0;
179         for (k=0; k<npe; k++) {
180           swarm_coor[dim*pcnt+d] += basis[q][k] * PetscRealPart(elcoor[dim*k+d]);
181         }
182       }
183       swarm_cellid[pcnt] = e;
184       pcnt++;
185     }
186     ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr);
187   }
188   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
189   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
190 
191   ierr = PetscFree(xi);CHKERRQ(ierr);
192   for (q=0; q<npoints_q; q++) {
193     ierr = PetscFree(basis[q]);CHKERRQ(ierr);
194   }
195   ierr = PetscFree(basis);CHKERRQ(ierr);
196 
197   PetscFunctionReturn(0);
198 }
199 
200 PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param)
201 {
202   PetscErrorCode ierr;
203   PetscInt dim;
204 
205   PetscFunctionBegin;
206   ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr);
207   if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No 3D support for PLEX");
208   switch (layout) {
209     case DMSWARMPIC_LAYOUT_REGULAR:
210       ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm,celldm,layout_param);CHKERRQ(ierr);
211       break;
212     case DMSWARMPIC_LAYOUT_GAUSS:
213       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Gauss layout not supported for PLEX");
214       break;
215     case DMSWARMPIC_LAYOUT_SUBDIVISION:
216       ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(dm,celldm,layout_param);CHKERRQ(ierr);
217       break;
218   }
219   PetscFunctionReturn(0);
220 }
221 
222 /*
223 typedef struct {
224   PetscReal x,y;
225 } Point2d;
226 
227 static PetscErrorCode signp2d(Point2d p1, Point2d p2, Point2d p3,PetscReal *s)
228 {
229   PetscFunctionBegin;
230   *s = (p1.x - p3.x) * (p2.y - p3.y) - (p2.x - p3.x) * (p1.y - p3.y);
231   PetscFunctionReturn(0);
232 }
233 */
234 /*
235 static PetscErrorCode PointInTriangle(Point2d pt, Point2d v1, Point2d v2, Point2d v3,PetscBool *v)
236 {
237   PetscReal s1,s2,s3;
238   PetscBool b1, b2, b3;
239 
240   PetscFunctionBegin;
241   signp2d(pt, v1, v2,&s1); b1 = s1 < 0.0f;
242   signp2d(pt, v2, v3,&s2); b2 = s2 < 0.0f;
243   signp2d(pt, v3, v1,&s3); b3 = s3 < 0.0f;
244 
245   *v = ((b1 == b2) && (b2 == b3));
246   PetscFunctionReturn(0);
247 }
248 */
249 /*
250 static PetscErrorCode _ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscReal coords[],PetscReal xip[],PetscReal *dJ)
251 {
252   PetscReal x1,y1,x2,y2,x3,y3;
253   PetscReal c,b[2],A[2][2],inv[2][2],detJ,od;
254 
255   PetscFunctionBegin;
256   x1 = coords[2*0+0];
257   x2 = coords[2*1+0];
258   x3 = coords[2*2+0];
259 
260   y1 = coords[2*0+1];
261   y2 = coords[2*1+1];
262   y3 = coords[2*2+1];
263 
264   c = x1 - 0.5*x1 - 0.5*x1 + 0.5*x2 + 0.5*x3;
265   b[0] = xp[0] - c;
266   c = y1 - 0.5*y1 - 0.5*y1 + 0.5*y2 + 0.5*y3;
267   b[1] = xp[1] - c;
268 
269   A[0][0] = -0.5*x1 + 0.5*x2;   A[0][1] = -0.5*x1 + 0.5*x3;
270   A[1][0] = -0.5*y1 + 0.5*y2;   A[1][1] = -0.5*y1 + 0.5*y3;
271 
272   detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0];
273   *dJ = PetscAbsReal(detJ);
274   od = 1.0/detJ;
275 
276   inv[0][0] =  A[1][1] * od;
277   inv[0][1] = -A[0][1] * od;
278   inv[1][0] = -A[1][0] * od;
279   inv[1][1] =  A[0][0] * od;
280 
281   xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1];
282   xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1];
283   PetscFunctionReturn(0);
284 }
285 */
286 
287 static PetscErrorCode ComputeLocalCoordinateAffine2d(PetscReal xp[],PetscScalar coords[],PetscReal xip[],PetscReal *dJ)
288 {
289   PetscReal x1,y1,x2,y2,x3,y3;
290   PetscReal b[2],A[2][2],inv[2][2],detJ,od;
291 
292   PetscFunctionBegin;
293   x1 = PetscRealPart(coords[2*0+0]);
294   x2 = PetscRealPart(coords[2*1+0]);
295   x3 = PetscRealPart(coords[2*2+0]);
296 
297   y1 = PetscRealPart(coords[2*0+1]);
298   y2 = PetscRealPart(coords[2*1+1]);
299   y3 = PetscRealPart(coords[2*2+1]);
300 
301   b[0] = xp[0] - x1;
302   b[1] = xp[1] - y1;
303 
304   A[0][0] = x2-x1;   A[0][1] = x3-x1;
305   A[1][0] = y2-y1;   A[1][1] = y3-y1;
306 
307   detJ = A[0][0]*A[1][1] - A[0][1]*A[1][0];
308   *dJ = PetscAbsReal(detJ);
309   od = 1.0/detJ;
310 
311   inv[0][0] =  A[1][1] * od;
312   inv[0][1] = -A[0][1] * od;
313   inv[1][0] = -A[1][0] * od;
314   inv[1][1] =  A[0][0] * od;
315 
316   xip[0] = inv[0][0]*b[0] + inv[0][1]*b[1];
317   xip[1] = inv[1][0]*b[0] + inv[1][1]*b[1];
318   PetscFunctionReturn(0);
319 }
320 
321 PetscErrorCode DMSwarmProjectField_ApproxP1_PLEX_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field)
322 {
323   PetscErrorCode ierr;
324   const PetscReal PLEX_C_EPS = 1.0e-8;
325   Vec v_field_l,denom_l,coor_l,denom;
326   PetscInt k,p,e,npoints;
327   PetscInt *mpfield_cell;
328   PetscReal *mpfield_coor;
329   PetscReal xi_p[2];
330   PetscScalar Ni[3];
331   PetscSection coordSection;
332   PetscScalar *elcoor = NULL;
333 
334   PetscFunctionBegin;
335   ierr = VecZeroEntries(v_field);CHKERRQ(ierr);
336 
337   ierr = DMGetLocalVector(dm,&v_field_l);CHKERRQ(ierr);
338   ierr = DMGetGlobalVector(dm,&denom);CHKERRQ(ierr);
339   ierr = DMGetLocalVector(dm,&denom_l);CHKERRQ(ierr);
340   ierr = VecZeroEntries(v_field_l);CHKERRQ(ierr);
341   ierr = VecZeroEntries(denom);CHKERRQ(ierr);
342   ierr = VecZeroEntries(denom_l);CHKERRQ(ierr);
343 
344   ierr = DMGetCoordinatesLocal(dm,&coor_l);CHKERRQ(ierr);
345   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
346 
347   ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr);
348   ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr);
349   ierr = DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr);
350 
351   for (p=0; p<npoints; p++) {
352     PetscReal   *coor_p,dJ;
353     PetscScalar elfield[3];
354     PetscBool   point_located;
355 
356     e       = mpfield_cell[p];
357     coor_p  = &mpfield_coor[2*p];
358 
359     ierr = DMPlexVecGetClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr);
360 
361 /*
362     while (!point_located && (failed_counter < 25)) {
363       ierr = PointInTriangle(point, coords[0], coords[1], coords[2], &point_located);CHKERRQ(ierr);
364       point.x = coor_p[0];
365       point.y = coor_p[1];
366       point.x += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0);
367       point.y += 1.0e-10 * (2.0 * rand()/((double)RAND_MAX)-1.0);
368       failed_counter++;
369     }
370 
371     if (!point_located){
372         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);
373     }
374 
375     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);
376     else {
377       ierr = _ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr);
378       xi_p[0] = 0.5*(xi_p[0] + 1.0);
379       xi_p[1] = 0.5*(xi_p[1] + 1.0);
380 
381       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]);
382 
383     }
384 */
385 
386     ierr = ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ);CHKERRQ(ierr);
387     /*
388     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]);
389     */
390     /*
391      point_located = PETSC_TRUE;
392     if (xi_p[0] < 0.0) {
393       if (xi_p[0] > -PLEX_C_EPS) {
394         xi_p[0] = 0.0;
395       } else {
396         point_located = PETSC_FALSE;
397       }
398     }
399     if (xi_p[1] < 0.0) {
400       if (xi_p[1] > -PLEX_C_EPS) {
401         xi_p[1] = 0.0;
402       } else {
403         point_located = PETSC_FALSE;
404       }
405     }
406     if (xi_p[1] > (1.0-xi_p[0])) {
407       if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) {
408         xi_p[1] = 1.0 - xi_p[0];
409       } else {
410         point_located = PETSC_FALSE;
411       }
412     }
413     if (!point_located){
414       PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]);
415       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]);
416     }
417     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);
418     */
419 
420     Ni[0] = 1.0 - xi_p[0] - xi_p[1];
421     Ni[1] = xi_p[0];
422     Ni[2] = xi_p[1];
423 
424     point_located = PETSC_TRUE;
425     for (k=0; k<3; k++) {
426       if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE;
427       if (PetscRealPart(Ni[k]) > (1.0+PLEX_C_EPS)) point_located = PETSC_FALSE;
428     }
429     if (!point_located){
430       PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",(double)xi_p[0],(double)xi_p[1]);
431       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]));
432     }
433     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);
434 
435     for (k=0; k<3; k++) {
436       Ni[k] = Ni[k] * dJ;
437       elfield[k] = Ni[k] * swarm_field[p];
438     }
439 
440     ierr = DMPlexVecRestoreClosure(dm,coordSection,coor_l,e,NULL,&elcoor);CHKERRQ(ierr);
441 
442     ierr = DMPlexVecSetClosure(dm, NULL,v_field_l, e, elfield, ADD_VALUES);CHKERRQ(ierr);
443     ierr = DMPlexVecSetClosure(dm, NULL,denom_l, e, Ni, ADD_VALUES);CHKERRQ(ierr);
444   }
445 
446   ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr);
447   ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr);
448 
449   ierr = DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr);
450   ierr = DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr);
451   ierr = DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr);
452   ierr = DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr);
453 
454   ierr = VecPointwiseDivide(v_field,v_field,denom);CHKERRQ(ierr);
455 
456   ierr = DMRestoreLocalVector(dm,&v_field_l);CHKERRQ(ierr);
457   ierr = DMRestoreLocalVector(dm,&denom_l);CHKERRQ(ierr);
458   ierr = DMRestoreGlobalVector(dm,&denom);CHKERRQ(ierr);
459 
460   PetscFunctionReturn(0);
461 }
462 
463 PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DataField dfield[],Vec vecs[])
464 {
465   PetscErrorCode ierr;
466   PetscInt f,dim;
467 
468   PetscFunctionBegin;
469   ierr = DMGetDimension(swarm,&dim);CHKERRQ(ierr);
470   switch (dim) {
471     case 2:
472       for (f=0; f<nfields; f++) {
473         PetscReal *swarm_field;
474 
475         ierr = DataFieldGetEntries(dfield[f],(void**)&swarm_field);CHKERRQ(ierr);
476         ierr = DMSwarmProjectField_ApproxP1_PLEX_2D(swarm,swarm_field,celldm,vecs[f]);CHKERRQ(ierr);
477       }
478       break;
479     case 3:
480       SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D");
481       break;
482     default:
483       break;
484   }
485 
486   PetscFunctionReturn(0);
487 }
488