xref: /petsc/src/dm/impls/swarm/swarmpic.c (revision 0e2ec84fac7d346010fd148bcbaea49eadbc0158)
1*0e2ec84fSDave May 
2*0e2ec84fSDave May #define PETSCDM_DLL
3*0e2ec84fSDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
4*0e2ec84fSDave May #include <petscsf.h>
5*0e2ec84fSDave May 
6*0e2ec84fSDave May /*
7*0e2ec84fSDave May  Error chceking macto to ensure the swarm type is correct and that a cell DM has been set
8*0e2ec84fSDave May */
9*0e2ec84fSDave May #define DMSWARMPICVALID(dm) \
10*0e2ec84fSDave May { \
11*0e2ec84fSDave May   DM_Swarm *_swarm = (DM_Swarm*)(dm)->data; \
12*0e2ec84fSDave May   if (_swarm->swarm_type != DMSWARM_PIC) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Only valid for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \
13*0e2ec84fSDave May   else \
14*0e2ec84fSDave May     if (!_swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Only valid for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM(dm,celldm)"); \
15*0e2ec84fSDave May }
16*0e2ec84fSDave May 
17*0e2ec84fSDave May /* Coordinate insertition/addition API */
18*0e2ec84fSDave May /*@C
19*0e2ec84fSDave May    DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid
20*0e2ec84fSDave May 
21*0e2ec84fSDave May    Collective on DM
22*0e2ec84fSDave May 
23*0e2ec84fSDave May    Input parameters:
24*0e2ec84fSDave May +  dm - the DMSwarm
25*0e2ec84fSDave May .  min - minimum coordinate values in the x, y, z directions (array of length dim)
26*0e2ec84fSDave May .  max - maximum coordinate values in the x, y, z directions (array of length dim)
27*0e2ec84fSDave May .  npoints - number of points in each spatial direction (array of length dim)
28*0e2ec84fSDave May -  mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES)
29*0e2ec84fSDave May 
30*0e2ec84fSDave May    Level: beginner
31*0e2ec84fSDave May 
32*0e2ec84fSDave May    Notes:
33*0e2ec84fSDave May    When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm
34*0e2ec84fSDave May    to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES,
35*0e2ec84fSDave May    new points will be appended to any already existing in the DMSwarm
36*0e2ec84fSDave May 
37*0e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType
38*0e2ec84fSDave May @*/
39*0e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm,PetscReal min[],PetscReal max[],PetscInt npoints[],InsertMode mode)
40*0e2ec84fSDave May {
41*0e2ec84fSDave May   PetscErrorCode ierr;
42*0e2ec84fSDave May   PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL};
43*0e2ec84fSDave May   PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
44*0e2ec84fSDave May   PetscInt i,j,k,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found;
45*0e2ec84fSDave May   Vec coorlocal;
46*0e2ec84fSDave May   const PetscScalar *_coor;
47*0e2ec84fSDave May   DM celldm;
48*0e2ec84fSDave May   PetscReal dx[3];
49*0e2ec84fSDave May   Vec pos;
50*0e2ec84fSDave May   PetscScalar *_pos;
51*0e2ec84fSDave May   PetscReal *swarm_coor;
52*0e2ec84fSDave May   PetscInt *swarm_cellid;
53*0e2ec84fSDave May   PetscSF sfcell = NULL;
54*0e2ec84fSDave May   const PetscSFNode *LA_sfcell;
55*0e2ec84fSDave May 
56*0e2ec84fSDave May   PetscFunctionBegin;
57*0e2ec84fSDave May   DMSWARMPICVALID(dm);
58*0e2ec84fSDave May   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
59*0e2ec84fSDave May   ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr);
60*0e2ec84fSDave May   ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr);
61*0e2ec84fSDave May   ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr);
62*0e2ec84fSDave May   N = N / bs;
63*0e2ec84fSDave May   ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
64*0e2ec84fSDave May   for (i=0; i<N; i++) {
65*0e2ec84fSDave May     for (b=0; b<bs; b++) {
66*0e2ec84fSDave May       gmin[b] = PetscMin(gmin[b],_coor[bs*i+b]);
67*0e2ec84fSDave May       gmax[b] = PetscMax(gmax[b],_coor[bs*i+b]);
68*0e2ec84fSDave May     }
69*0e2ec84fSDave May   }
70*0e2ec84fSDave May   ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
71*0e2ec84fSDave May 
72*0e2ec84fSDave May   for (b=0; b<bs; b++) {
73*0e2ec84fSDave May     dx[b] = (max[b] - min[b])/((PetscReal)(npoints[b]-1));
74*0e2ec84fSDave May   }
75*0e2ec84fSDave May 
76*0e2ec84fSDave May   /* determine number of points living in the bounding box */
77*0e2ec84fSDave May   n_estimate = 0;
78*0e2ec84fSDave May   if (bs == 2) { npoints[2] = 1; }
79*0e2ec84fSDave May   for (k=0; k<npoints[2]; k++) {
80*0e2ec84fSDave May     for (j=0; j<npoints[1]; j++) {
81*0e2ec84fSDave May       for (i=0; i<npoints[0]; i++) {
82*0e2ec84fSDave May         PetscReal xp[] = {0.0,0.0,0.0};
83*0e2ec84fSDave May         PetscInt ijk[3];
84*0e2ec84fSDave May         PetscBool point_inside = PETSC_TRUE;
85*0e2ec84fSDave May 
86*0e2ec84fSDave May         ijk[0] = i;
87*0e2ec84fSDave May         ijk[1] = j;
88*0e2ec84fSDave May         ijk[2] = k;
89*0e2ec84fSDave May         for (b=0; b<bs; b++) {
90*0e2ec84fSDave May           xp[b] = min[b] + ijk[b] * dx[b];
91*0e2ec84fSDave May         }
92*0e2ec84fSDave May         for (b=0; b<bs; b++) {
93*0e2ec84fSDave May           if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; }
94*0e2ec84fSDave May           if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; }
95*0e2ec84fSDave May         }
96*0e2ec84fSDave May         if (point_inside) { n_estimate++; }
97*0e2ec84fSDave May       }
98*0e2ec84fSDave May     }
99*0e2ec84fSDave May   }
100*0e2ec84fSDave May 
101*0e2ec84fSDave May   /* create candidate list */
102*0e2ec84fSDave May   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&pos);CHKERRQ(ierr);
103*0e2ec84fSDave May   ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr);
104*0e2ec84fSDave May   ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr);
105*0e2ec84fSDave May   ierr = VecSetFromOptions(pos);CHKERRQ(ierr);
106*0e2ec84fSDave May   ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr);
107*0e2ec84fSDave May 
108*0e2ec84fSDave May   n_estimate = 0;
109*0e2ec84fSDave May   for (k=0; k<npoints[2]; k++) {
110*0e2ec84fSDave May     for (j=0; j<npoints[1]; j++) {
111*0e2ec84fSDave May       for (i=0; i<npoints[0]; i++) {
112*0e2ec84fSDave May         PetscReal xp[] = {0.0,0.0,0.0};
113*0e2ec84fSDave May         PetscInt ijk[3];
114*0e2ec84fSDave May         PetscBool point_inside = PETSC_TRUE;
115*0e2ec84fSDave May 
116*0e2ec84fSDave May         ijk[0] = i;
117*0e2ec84fSDave May         ijk[1] = j;
118*0e2ec84fSDave May         ijk[2] = k;
119*0e2ec84fSDave May         for (b=0; b<bs; b++) {
120*0e2ec84fSDave May           xp[b] = min[b] + ijk[b] * dx[b];
121*0e2ec84fSDave May         }
122*0e2ec84fSDave May         for (b=0; b<bs; b++) {
123*0e2ec84fSDave May           if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; }
124*0e2ec84fSDave May           if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; }
125*0e2ec84fSDave May         }
126*0e2ec84fSDave May         if (point_inside) {
127*0e2ec84fSDave May           for (b=0; b<bs; b++) {
128*0e2ec84fSDave May             _pos[bs*n_estimate+b] = xp[b];
129*0e2ec84fSDave May           }
130*0e2ec84fSDave May           n_estimate++;
131*0e2ec84fSDave May         }
132*0e2ec84fSDave May       }
133*0e2ec84fSDave May     }
134*0e2ec84fSDave May   }
135*0e2ec84fSDave May   ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr);
136*0e2ec84fSDave May 
137*0e2ec84fSDave May   /* locate points */
138*0e2ec84fSDave May   ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
139*0e2ec84fSDave May 
140*0e2ec84fSDave May   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
141*0e2ec84fSDave May   n_found = 0;
142*0e2ec84fSDave May   for (p=0; p<n_estimate; p++) {
143*0e2ec84fSDave May     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
144*0e2ec84fSDave May       n_found++;
145*0e2ec84fSDave May     }
146*0e2ec84fSDave May   }
147*0e2ec84fSDave May 
148*0e2ec84fSDave May   /* adjust size */
149*0e2ec84fSDave May   if (mode == ADD_VALUES) {
150*0e2ec84fSDave May     ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr);
151*0e2ec84fSDave May     n_new_est = n_curr + n_found;
152*0e2ec84fSDave May     ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr);
153*0e2ec84fSDave May   }
154*0e2ec84fSDave May   if (mode == INSERT_VALUES) {
155*0e2ec84fSDave May     n_curr = 0;
156*0e2ec84fSDave May     n_new_est = n_found;
157*0e2ec84fSDave May     ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr);
158*0e2ec84fSDave May   }
159*0e2ec84fSDave May 
160*0e2ec84fSDave May   /* initialize new coords, cell owners, pid */
161*0e2ec84fSDave May   ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr);
162*0e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
163*0e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
164*0e2ec84fSDave May   n_found = 0;
165*0e2ec84fSDave May   for (p=0; p<n_estimate; p++) {
166*0e2ec84fSDave May     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
167*0e2ec84fSDave May       for (b=0; b<bs; b++) {
168*0e2ec84fSDave May         swarm_coor[bs*(n_curr + n_found) + b] = _coor[bs*p+b];
169*0e2ec84fSDave May       }
170*0e2ec84fSDave May       swarm_cellid[n_curr + n_found] = LA_sfcell[p].index;
171*0e2ec84fSDave May       n_found++;
172*0e2ec84fSDave May     }
173*0e2ec84fSDave May   }
174*0e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
175*0e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
176*0e2ec84fSDave May   ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr);
177*0e2ec84fSDave May 
178*0e2ec84fSDave May   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
179*0e2ec84fSDave May   ierr = VecDestroy(&pos);CHKERRQ(ierr);
180*0e2ec84fSDave May 
181*0e2ec84fSDave May   PetscFunctionReturn(0);
182*0e2ec84fSDave May }
183*0e2ec84fSDave May 
184*0e2ec84fSDave May /*@C
185*0e2ec84fSDave May    DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list
186*0e2ec84fSDave May 
187*0e2ec84fSDave May    Collective on DM
188*0e2ec84fSDave May 
189*0e2ec84fSDave May    Input parameters:
190*0e2ec84fSDave May +  dm - the DMSwarm
191*0e2ec84fSDave May .  npoints - the number of points to insert
192*0e2ec84fSDave May .  coor - the coordinate values
193*0e2ec84fSDave May .  redundant - if set to PETSC_TRUE, it is assumed that npoints and coor[] are only valid on rank 0 and should be broadcast to other ranks
194*0e2ec84fSDave May -  mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES)
195*0e2ec84fSDave May 
196*0e2ec84fSDave May    Level: beginner
197*0e2ec84fSDave May 
198*0e2ec84fSDave May    Notes:
199*0e2ec84fSDave May    If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within
200*0e2ec84fSDave May    its sub-domain. If they any values within coor[] are not located in the sub-domain, they will be ignored and will not get
201*0e2ec84fSDave May    added to the DMSwarm.
202*0e2ec84fSDave May 
203*0e2ec84fSDave May 
204*0e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType, DMSwarmSetPointsUniformCoordinates()
205*0e2ec84fSDave May @*/
206*0e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm,PetscInt npoints,PetscReal coor[],PetscBool redundant,InsertMode mode)
207*0e2ec84fSDave May {
208*0e2ec84fSDave May   PetscErrorCode ierr;
209*0e2ec84fSDave May   PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL};
210*0e2ec84fSDave May   PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
211*0e2ec84fSDave May   PetscInt i,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found;
212*0e2ec84fSDave May   Vec coorlocal;
213*0e2ec84fSDave May   const PetscScalar *_coor;
214*0e2ec84fSDave May   DM celldm;
215*0e2ec84fSDave May   Vec pos;
216*0e2ec84fSDave May   PetscScalar *_pos;
217*0e2ec84fSDave May   PetscReal *swarm_coor;
218*0e2ec84fSDave May   PetscInt *swarm_cellid;
219*0e2ec84fSDave May   PetscSF sfcell = NULL;
220*0e2ec84fSDave May   const PetscSFNode *LA_sfcell;
221*0e2ec84fSDave May   PetscReal *my_coor;
222*0e2ec84fSDave May   PetscInt my_npoints;
223*0e2ec84fSDave May   PetscMPIInt rank;
224*0e2ec84fSDave May   MPI_Comm comm;
225*0e2ec84fSDave May 
226*0e2ec84fSDave May   PetscFunctionBegin;
227*0e2ec84fSDave May   DMSWARMPICVALID(dm);
228*0e2ec84fSDave May   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
229*0e2ec84fSDave May   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
230*0e2ec84fSDave May 
231*0e2ec84fSDave May   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
232*0e2ec84fSDave May   ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr);
233*0e2ec84fSDave May   ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr);
234*0e2ec84fSDave May   ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr);
235*0e2ec84fSDave May   N = N / bs;
236*0e2ec84fSDave May   ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
237*0e2ec84fSDave May   for (i=0; i<N; i++) {
238*0e2ec84fSDave May     for (b=0; b<bs; b++) {
239*0e2ec84fSDave May       gmin[b] = PetscMin(gmin[b],_coor[bs*i+b]);
240*0e2ec84fSDave May       gmax[b] = PetscMax(gmax[b],_coor[bs*i+b]);
241*0e2ec84fSDave May     }
242*0e2ec84fSDave May   }
243*0e2ec84fSDave May   ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr);
244*0e2ec84fSDave May 
245*0e2ec84fSDave May   /* broadcast points from rank 0 if requested */
246*0e2ec84fSDave May   if (redundant) {
247*0e2ec84fSDave May     my_npoints = npoints;
248*0e2ec84fSDave May     ierr = MPI_Bcast(&my_npoints,1,MPIU_INT,0,comm);CHKERRQ(ierr);
249*0e2ec84fSDave May 
250*0e2ec84fSDave May     if (rank > 0) { /* allocate space */
251*0e2ec84fSDave May       ierr = PetscMalloc1(my_npoints,&my_coor);CHKERRQ(ierr);
252*0e2ec84fSDave May     } else {
253*0e2ec84fSDave May       my_coor = coor;
254*0e2ec84fSDave May     }
255*0e2ec84fSDave May     ierr = MPI_Bcast(my_coor,bs*my_npoints,MPIU_REAL,0,comm);CHKERRQ(ierr);
256*0e2ec84fSDave May   } else {
257*0e2ec84fSDave May     my_npoints = npoints;
258*0e2ec84fSDave May     my_coor = coor;
259*0e2ec84fSDave May   }
260*0e2ec84fSDave May 
261*0e2ec84fSDave May   /* determine the number of points living in the bounding box */
262*0e2ec84fSDave May   n_estimate = 0;
263*0e2ec84fSDave May   for (i=0; i<my_npoints; i++) {
264*0e2ec84fSDave May     PetscBool point_inside = PETSC_TRUE;
265*0e2ec84fSDave May 
266*0e2ec84fSDave May     for (b=0; b<bs; b++) {
267*0e2ec84fSDave May       if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; }
268*0e2ec84fSDave May       if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; }
269*0e2ec84fSDave May     }
270*0e2ec84fSDave May     if (point_inside) { n_estimate++; }
271*0e2ec84fSDave May   }
272*0e2ec84fSDave May 
273*0e2ec84fSDave May   /* create candidate list */
274*0e2ec84fSDave May   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&pos);CHKERRQ(ierr);
275*0e2ec84fSDave May   ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr);
276*0e2ec84fSDave May   ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr);
277*0e2ec84fSDave May   ierr = VecSetFromOptions(pos);CHKERRQ(ierr);
278*0e2ec84fSDave May   ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr);
279*0e2ec84fSDave May 
280*0e2ec84fSDave May   n_estimate = 0;
281*0e2ec84fSDave May   for (i=0; i<my_npoints; i++) {
282*0e2ec84fSDave May     PetscBool point_inside = PETSC_TRUE;
283*0e2ec84fSDave May 
284*0e2ec84fSDave May     for (b=0; b<bs; b++) {
285*0e2ec84fSDave May       if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; }
286*0e2ec84fSDave May       if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; }
287*0e2ec84fSDave May     }
288*0e2ec84fSDave May     if (point_inside) {
289*0e2ec84fSDave May       for (b=0; b<bs; b++) {
290*0e2ec84fSDave May         _pos[bs*n_estimate+b] = my_coor[bs*i+b];
291*0e2ec84fSDave May       }
292*0e2ec84fSDave May       n_estimate++;
293*0e2ec84fSDave May     }
294*0e2ec84fSDave May   }
295*0e2ec84fSDave May   ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr);
296*0e2ec84fSDave May 
297*0e2ec84fSDave May   /* locate points */
298*0e2ec84fSDave May   ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr);
299*0e2ec84fSDave May 
300*0e2ec84fSDave May   ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr);
301*0e2ec84fSDave May   n_found = 0;
302*0e2ec84fSDave May   for (p=0; p<n_estimate; p++) {
303*0e2ec84fSDave May     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
304*0e2ec84fSDave May       n_found++;
305*0e2ec84fSDave May     }
306*0e2ec84fSDave May   }
307*0e2ec84fSDave May 
308*0e2ec84fSDave May   /* adjust size */
309*0e2ec84fSDave May   if (mode == ADD_VALUES) {
310*0e2ec84fSDave May     ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr);
311*0e2ec84fSDave May     n_new_est = n_curr + n_found;
312*0e2ec84fSDave May     ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr);
313*0e2ec84fSDave May   }
314*0e2ec84fSDave May   if (mode == INSERT_VALUES) {
315*0e2ec84fSDave May     n_curr = 0;
316*0e2ec84fSDave May     n_new_est = n_found;
317*0e2ec84fSDave May     ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr);
318*0e2ec84fSDave May   }
319*0e2ec84fSDave May 
320*0e2ec84fSDave May   /* initialize new coords, cell owners, pid */
321*0e2ec84fSDave May   ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr);
322*0e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
323*0e2ec84fSDave May   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
324*0e2ec84fSDave May   n_found = 0;
325*0e2ec84fSDave May   for (p=0; p<n_estimate; p++) {
326*0e2ec84fSDave May     if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) {
327*0e2ec84fSDave May       for (b=0; b<bs; b++) {
328*0e2ec84fSDave May         swarm_coor[bs*(n_curr + n_found) + b] = _coor[bs*p+b];
329*0e2ec84fSDave May       }
330*0e2ec84fSDave May       swarm_cellid[n_curr + n_found] = LA_sfcell[p].index;
331*0e2ec84fSDave May       n_found++;
332*0e2ec84fSDave May     }
333*0e2ec84fSDave May   }
334*0e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
335*0e2ec84fSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr);
336*0e2ec84fSDave May   ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr);
337*0e2ec84fSDave May 
338*0e2ec84fSDave May   if (redundant) {
339*0e2ec84fSDave May     if (rank > 0) {
340*0e2ec84fSDave May       ierr = PetscFree(my_coor);CHKERRQ(ierr);
341*0e2ec84fSDave May     }
342*0e2ec84fSDave May   }
343*0e2ec84fSDave May   ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr);
344*0e2ec84fSDave May   ierr = VecDestroy(&pos);CHKERRQ(ierr);
345*0e2ec84fSDave May 
346*0e2ec84fSDave May   PetscFunctionReturn(0);
347*0e2ec84fSDave May }
348*0e2ec84fSDave May 
349*0e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM,DM,DMSwarmPICLayoutType,PetscInt);
350*0e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM,DM,DMSwarmPICLayoutType,PetscInt);
351*0e2ec84fSDave May 
352*0e2ec84fSDave May /*@C
353*0e2ec84fSDave May    DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell
354*0e2ec84fSDave May 
355*0e2ec84fSDave May    Not collective
356*0e2ec84fSDave May 
357*0e2ec84fSDave May    Input parameters:
358*0e2ec84fSDave May +  dm - the DMSwarm
359*0e2ec84fSDave May .  layout_type - method used to fill each cell with the cell DM
360*0e2ec84fSDave May -  fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type)
361*0e2ec84fSDave May 
362*0e2ec84fSDave May  Level: beginner
363*0e2ec84fSDave May 
364*0e2ec84fSDave May  Notes:
365*0e2ec84fSDave May  The insert method will reset any previous defined points within the DMSwarm
366*0e2ec84fSDave May 
367*0e2ec84fSDave May .seealso: DMSwarmPICLayoutType, DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType
368*0e2ec84fSDave May @*/
369*0e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm,DMSwarmPICLayoutType layout_type,PetscInt fill_param)
370*0e2ec84fSDave May {
371*0e2ec84fSDave May   PetscErrorCode ierr;
372*0e2ec84fSDave May   DM celldm;
373*0e2ec84fSDave May   PetscBool isDA,isPLEX;
374*0e2ec84fSDave May 
375*0e2ec84fSDave May   PetscFunctionBegin;
376*0e2ec84fSDave May   DMSWARMPICVALID(dm);
377*0e2ec84fSDave May   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
378*0e2ec84fSDave May   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr);
379*0e2ec84fSDave May   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr);
380*0e2ec84fSDave May   if (isDA) {
381*0e2ec84fSDave May     ierr = private_DMSwarmInsertPointsUsingCellDM_DA(dm,celldm,layout_type,fill_param);CHKERRQ(ierr);
382*0e2ec84fSDave May   } else if (isPLEX) {
383*0e2ec84fSDave May     ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX(dm,celldm,layout_type,fill_param);CHKERRQ(ierr);
384*0e2ec84fSDave May   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX");
385*0e2ec84fSDave May 
386*0e2ec84fSDave May   PetscFunctionReturn(0);
387*0e2ec84fSDave May }
388*0e2ec84fSDave May 
389*0e2ec84fSDave May /*
390*0e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmAddPointCoordinatesCellWise(DM dm,PetscInt cell,PetscInt npoints,PetscReal xi[],PetscBool proximity_initialization)
391*0e2ec84fSDave May {
392*0e2ec84fSDave May   PetscFunctionBegin;
393*0e2ec84fSDave May   PetscFunctionReturn(0);
394*0e2ec84fSDave May }
395*0e2ec84fSDave May */
396*0e2ec84fSDave May 
397*0e2ec84fSDave May /* Field projection API */
398*0e2ec84fSDave May /*
399*0e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm,PetscInt project_type,PetscInt nfields,const char *fieldnames[],Vec *fields)
400*0e2ec84fSDave May {
401*0e2ec84fSDave May   PetscFunctionBegin;
402*0e2ec84fSDave May   PetscFunctionReturn(0);
403*0e2ec84fSDave May }
404*0e2ec84fSDave May */
405*0e2ec84fSDave May 
406*0e2ec84fSDave May /*@C
407*0e2ec84fSDave May    DMSwarmCreatePointPerCellCount - Count the number of points (particles) per cell
408*0e2ec84fSDave May 
409*0e2ec84fSDave May    Not collective
410*0e2ec84fSDave May 
411*0e2ec84fSDave May    Input parameter:
412*0e2ec84fSDave May .  dm - the DMSwarm
413*0e2ec84fSDave May 
414*0e2ec84fSDave May    Output parameters:
415*0e2ec84fSDave May +  ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore)
416*0e2ec84fSDave May -  count - array of length ncells containing the number of particles per cell
417*0e2ec84fSDave May 
418*0e2ec84fSDave May    Level: beginner
419*0e2ec84fSDave May 
420*0e2ec84fSDave May    Notes:
421*0e2ec84fSDave May    The array count is allocated internally and must be free'd by the user.
422*0e2ec84fSDave May 
423*0e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType
424*0e2ec84fSDave May @*/
425*0e2ec84fSDave May /*
426*0e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm,PetscInt *ncells,PetscInt **count)
427*0e2ec84fSDave May {
428*0e2ec84fSDave May   PetscFunctionBegin;
429*0e2ec84fSDave May   PetscFunctionReturn(0);
430*0e2ec84fSDave May }
431*0e2ec84fSDave May */
432