xref: /petsc/src/dm/impls/swarm/swarmpic_sort.c (revision 3d96e9964ff330fd2a9eee374bcd4376da7efe60)
1 
2 #include <petsc.h>
3 #include <petscdm.h>
4 #include <petscdmda.h>
5 #include <petscdmplex.h>
6 #include <petscdmswarm.h>
7 #include <petsc/private/dmswarmimpl.h>
8 
9 int sort_CompareSwarmPoint(const void *dataA,const void *dataB)
10 {
11   SwarmPoint *pointA;
12   SwarmPoint *pointB;
13 
14   pointA = (SwarmPoint*)dataA;
15   pointB = (SwarmPoint*)dataB;
16 
17   if (pointA->cell_index < pointB->cell_index) {
18     return -1;
19   } else if (pointA->cell_index > pointB->cell_index) {
20     return 1;
21   } else {
22     return 0;
23   }
24 }
25 
26 PetscErrorCode DMSwarmSortApplyCellIndexSort(DMSwarmSort ctx)
27 {
28   PetscFunctionBegin;
29   qsort(ctx->list,ctx->npoints,sizeof(SwarmPoint),sort_CompareSwarmPoint);
30   PetscFunctionReturn(0);
31 }
32 
33 PetscErrorCode DMSwarmSortCreate(DMSwarmSort *_ctx)
34 {
35   PetscErrorCode ierr;
36   DMSwarmSort    ctx;
37 
38   PetscFunctionBegin;
39   ierr = PetscMalloc1(1,&ctx);CHKERRQ(ierr);
40   ierr = PetscMemzero(ctx,sizeof(struct _p_DMSwarmSort));CHKERRQ(ierr);
41   ctx->isvalid = PETSC_FALSE;
42   ctx->ncells = 0;
43   ctx->npoints = 0;
44   ierr = PetscMalloc1(1,&ctx->pcell_offsets);CHKERRQ(ierr);
45   ierr = PetscMalloc1(1,&ctx->list);CHKERRQ(ierr);
46   *_ctx = ctx;
47   PetscFunctionReturn(0);
48 }
49 
50 PetscErrorCode DMSwarmSortSetup(DMSwarmSort ctx,DM dm,PetscInt ncells)
51 {
52   PetscInt        *swarm_cellid;
53   PetscInt        p,npoints;
54   PetscInt        tmp,c,count;
55   PetscErrorCode  ierr;
56 
57   PetscFunctionBegin;
58   if (!ctx) PetscFunctionReturn(0);
59   if (ctx->isvalid) PetscFunctionReturn(0);
60 
61   ierr = PetscLogEventBegin(DMSWARM_Sort,0,0,0,0);CHKERRQ(ierr);
62   /* check the number of cells */
63   if (ncells != ctx->ncells) {
64     ierr = PetscRealloc(sizeof(PetscInt)*(ncells + 1),&ctx->pcell_offsets);CHKERRQ(ierr);
65     ctx->ncells = ncells;
66   }
67   ierr = PetscMemzero(ctx->pcell_offsets,sizeof(PetscInt)*(ctx->ncells + 1));CHKERRQ(ierr);
68 
69   /* get the number of points */
70   ierr = DMSwarmGetLocalSize(dm,&npoints);CHKERRQ(ierr);
71   if (npoints != ctx->npoints) {
72     ierr = PetscRealloc(sizeof(SwarmPoint)*npoints,&ctx->list);CHKERRQ(ierr);
73     ctx->npoints = npoints;
74   }
75   ierr = PetscMemzero(ctx->list,sizeof(SwarmPoint)*npoints);CHKERRQ(ierr);
76 
77   ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
78   for (p=0; p<ctx->npoints; p++) {
79     ctx->list[p].point_index = p;
80     ctx->list[p].cell_index  = swarm_cellid[p];
81   }
82   ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr);
83 
84   ierr = DMSwarmSortApplyCellIndexSort(ctx);CHKERRQ(ierr);
85 
86   /* sum points per cell */
87   for (p=0; p<ctx->npoints; p++) {
88     ctx->pcell_offsets[ ctx->list[p].cell_index ]++;
89   }
90 
91   /* create offset list */
92   count = 0;
93   for (c=0; c<ctx->ncells; c++) {
94     tmp = ctx->pcell_offsets[c];
95     ctx->pcell_offsets[c] = count;
96     count = count + tmp;
97   }
98   ctx->pcell_offsets[c] = count;
99 
100   ctx->isvalid = PETSC_TRUE;
101   ierr = PetscLogEventEnd(DMSWARM_Sort,0,0,0,0);CHKERRQ(ierr);
102   PetscFunctionReturn(0);
103 }
104 
105 PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx)
106 {
107   DMSwarmSort     ctx;
108   PetscErrorCode  ierr;
109 
110   PetscFunctionBegin;
111   if (!_ctx) PetscFunctionReturn(0);
112   if (!*_ctx) PetscFunctionReturn(0);
113   ctx = *_ctx;
114   if (ctx->list)      {
115     ierr = PetscFree(ctx->list);CHKERRQ(ierr);
116   }
117   if (ctx->pcell_offsets) {
118     ierr = PetscFree(ctx->pcell_offsets);CHKERRQ(ierr);
119   }
120   ierr = PetscFree(ctx);CHKERRQ(ierr);
121   *_ctx = NULL;
122   PetscFunctionReturn(0);
123 }
124 
125 /*@C
126    DMSwarmSortGetNumberOfPointsPerCell - Returns the number of points in a cell
127 
128    Not collective
129 
130    Input parameters:
131 +  dm - a DMSwarm objects
132 .  e - the index of the cell
133 -  npoints - the number of points in the cell
134 
135    Level: advanced
136 
137    Notes:
138    You must call DMSwarmSortGetAccess() before you can call DMSwarmSortGetNumberOfPointsPerCell()
139 
140 .seealso: DMSwarmSetType(), DMSwarmSortGetAccess(), DMSwarmSortGetPointsPerCell()
141 @*/
142 PETSC_EXTERN PetscErrorCode DMSwarmSortGetNumberOfPointsPerCell(DM dm,PetscInt e,PetscInt *npoints)
143 {
144   DM_Swarm     *swarm = (DM_Swarm*)dm->data;
145   PetscInt     points_per_cell;
146   DMSwarmSort  ctx;
147 
148   PetscFunctionBegin;
149   ctx = swarm->sort_context;
150   if (!ctx) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"The DMSwarmSort context has not been created. Must call DMSwarmSortGetAccess() first");
151   if (!ctx->isvalid) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"SwarmPointSort container is not valid. Must call DMSwarmSortGetAccess() first");
152   if (e >= ctx->ncells) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cell index (%D) is greater than max number of local cells (%D)",e,ctx->ncells);
153   if (e < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Cell index (%D) cannot be negative",e);
154   points_per_cell = ctx->pcell_offsets[e+1] - ctx->pcell_offsets[e];
155   *npoints = points_per_cell;
156   PetscFunctionReturn(0);
157 }
158 
159 /*@C
160    DMSwarmSortGetPointsPerCell - Creates an array of point indices for all points in a cell
161 
162    Not collective
163 
164    Input parameters:
165 +  dm - a DMSwarm object
166 .  e - the index of the cell
167 .  npoints - the number of points in the cell
168 -  pidlist - array of the indices indentifying all points in cell e
169 
170    Level: advanced
171 
172    Notes:
173    - You must call DMSwarmSortGetAccess() before you can call DMSwarmSortGetPointsPerCell()
174    - The array pidlist is internally created and must be free'd by the user
175 
176 .seealso: DMSwarmSetType(), DMSwarmSortGetAccess(), DMSwarmSortGetNumberOfPointsPerCell()
177 @*/
178 PETSC_EXTERN PetscErrorCode DMSwarmSortGetPointsPerCell(DM dm,PetscInt e,PetscInt *npoints,PetscInt **pidlist)
179 {
180   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
181   PetscErrorCode ierr;
182   PetscInt       points_per_cell;
183   PetscInt       p,pid,pid_unsorted;
184   PetscInt       *plist;
185   DMSwarmSort    ctx;
186 
187   PetscFunctionBegin;
188   ctx = swarm->sort_context;
189   if (!ctx) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"The DMSwarmSort context has not been created. Must call DMSwarmSortGetAccess() first");
190   ierr = DMSwarmSortGetNumberOfPointsPerCell(dm,e,&points_per_cell);CHKERRQ(ierr);
191   ierr = PetscMalloc1(points_per_cell,&plist);CHKERRQ(ierr);
192   for (p=0; p<points_per_cell; p++) {
193     pid = ctx->pcell_offsets[e] + p;
194     pid_unsorted = ctx->list[pid].point_index;
195     plist[p] = pid_unsorted;
196   }
197   *npoints = points_per_cell;
198   *pidlist = plist;
199 
200   PetscFunctionReturn(0);
201 }
202 
203 /*@C
204    DMSwarmSortGetAccess - Setups up a DMSwarm point sort context for efficient traversal of points within a cell
205 
206    Not collective
207 
208    Input parameter:
209 .  dm - a DMSwarm object
210 
211    Calling DMSwarmSortGetAccess() creates a list which enables easy identification of all points contained in a
212    given cell. This method does not explicitly sort the data within the DMSwarm based on the cell index associated
213    with a DMSwarm point.
214 
215    The sort context is valid only for the DMSwarm points defined at the time when DMSwarmSortGetAccess() was called.
216    For example, suppose the swarm contained NP points when DMSwarmSortGetAccess() was called. If the user subsequently
217    adds 10 additional points to the swarm, the sort context is still valid, but only for the first NP points.
218    The indices associated with the 10 new additional points will not be contained within the sort context.
219    This means that the user can still safely perform queries via DMSwarmSortGetPointsPerCell() and
220    DMSwarmSortGetPointsPerCell(), however the results return will be based on the first NP points.
221 
222    If any DMSwam re-sizing method is called after DMSwarmSortGetAccess() which modifies any of the first NP entries
223    in the DMSwarm, the sort context will become invalid. Currently there are no guards to prevent the user from
224    invalidating the sort context. For this reason, we highly recommend you do not use DMSwarmRemovePointAtIndex() in
225    between calls to DMSwarmSortGetAccess() and DMSwarmSortRestoreAccess().
226 
227    To facilitate safe removal of points using the sort context, we suggest a "two pass" strategy in which the
228    first pass "marks" points for removal, and the second pass actually removes the points from the DMSwarm.
229 
230    Notes:
231    - You must call DMSwarmSortGetAccess() before you can call DMSwarmSortGetPointsPerCell() or DMSwarmSortGetNumberOfPointsPerCell()
232    - The sort context may become invalid if any re-sizing methods are applied which alter the first NP points
233      within swarm at the time DMSwarmSortGetAccess() was called.
234    - You must call DMSwarmSortRestoreAccess() when you no longer need access to the sort context
235 
236    Level: advanced
237 
238 .seealso: DMSwarmSetType(), DMSwarmSortRestoreAccess()
239 @*/
240 PETSC_EXTERN PetscErrorCode DMSwarmSortGetAccess(DM dm)
241 {
242   DM_Swarm        *swarm = (DM_Swarm*)dm->data;
243   PetscErrorCode  ierr;
244   PetscInt        ncells;
245   DM              celldm;
246   PetscBool       isda,isplex,isshell;
247 
248   PetscFunctionBegin;
249   if (!swarm->sort_context) {
250     ierr = DMSwarmSortCreate(&swarm->sort_context);CHKERRQ(ierr);
251   }
252 
253   /* get the number of cells */
254   ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr);
255   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isda);CHKERRQ(ierr);
256   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isplex);CHKERRQ(ierr);
257   ierr = PetscObjectTypeCompare((PetscObject)celldm,DMSHELL,&isshell);CHKERRQ(ierr);
258   ncells = 0;
259   if (isda) {
260     PetscInt nel,npe;
261     const PetscInt *element;
262 
263     ierr = DMDAGetElements(celldm,&nel,&npe,&element);CHKERRQ(ierr);
264     ncells = nel;
265     ierr = DMDARestoreElements(celldm,&nel,&npe,&element);CHKERRQ(ierr);
266   } else if (isplex) {
267     PetscInt ps,pe;
268 
269     ierr = DMPlexGetHeightStratum(celldm,0,&ps,&pe);CHKERRQ(ierr);
270     ncells = pe - ps;
271   } else if (isshell) {
272     PetscErrorCode (*method_DMShellGetNumberOfCells)(DM,PetscInt*);
273 
274     ierr = PetscObjectQueryFunction((PetscObject)celldm,"DMGetNumberOfCells_C",&method_DMShellGetNumberOfCells);CHKERRQ(ierr);
275     if (method_DMShellGetNumberOfCells) {
276       ierr = method_DMShellGetNumberOfCells(celldm,&ncells);CHKERRQ(ierr);
277     } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot determine the number of cells for the DMSHELL object. User must provide a method via PetscObjectComposeFunction( (PetscObject)shelldm, \"DMGetNumberOfCells_C\", your_function_to_compute_number_of_cells );");
278   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL");
279 
280   /* setup */
281   ierr = DMSwarmSortSetup(swarm->sort_context,dm,ncells);CHKERRQ(ierr);
282 
283   PetscFunctionReturn(0);
284 }
285 
286 /*@C
287    DMSwarmSortRestoreAccess - Invalidates the DMSwarm point sorting context
288 
289    Not collective
290 
291    Input parameter:
292 .  dm - a DMSwarm object
293 
294    Level: advanced
295 
296    Note:
297    You must call DMSwarmSortGetAccess() before calling DMSwarmSortRestoreAccess()
298 
299 .seealso: DMSwarmSetType(), DMSwarmSortGetAccess()
300 @*/
301 PETSC_EXTERN PetscErrorCode DMSwarmSortRestoreAccess(DM dm)
302 {
303   DM_Swarm *swarm = (DM_Swarm*)dm->data;
304 
305   PetscFunctionBegin;
306   if (!swarm->sort_context) PetscFunctionReturn(0);
307   if (!swarm->sort_context->isvalid) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"You must call DMSwarmSortGetAccess() before calling DMSwarmSortRestoreAccess()");
308   swarm->sort_context->isvalid = PETSC_FALSE;
309   PetscFunctionReturn(0);
310 }
311 
312 /*@C
313    DMSwarmSortGetIsValid - Gets the isvalid flag associated with a DMSwarm point sorting context
314 
315    Not collective
316 
317    Input parameter:
318 .  dm - a DMSwarm object
319 
320    Output parameter:
321 .  isvalid - flag indicating whether the sort context is up-to-date
322 
323  Level: advanced
324 
325 .seealso: DMSwarmSetType(), DMSwarmSortGetAccess()
326 @*/
327 PETSC_EXTERN PetscErrorCode DMSwarmSortGetIsValid(DM dm,PetscBool *isvalid)
328 {
329   DM_Swarm *swarm = (DM_Swarm*)dm->data;
330 
331   PetscFunctionBegin;
332   if (!swarm->sort_context) {
333     *isvalid = PETSC_FALSE;
334     PetscFunctionReturn(0);
335   }
336   *isvalid = swarm->sort_context->isvalid;
337   PetscFunctionReturn(0);
338 }
339 
340 /*@C
341    DMSwarmSortGetSizes - Gets the sizes associated with a DMSwarm point sorting context
342 
343    Not collective
344 
345    Input parameter:
346 .  dm - a DMSwarm object
347 
348    Output parameters:
349 +  ncells - number of cells within the sort context (pass NULL to ignore)
350 -  npoints - number of points used to create the sort context (pass NULL to ignore)
351 
352    Level: advanced
353 
354 .seealso: DMSwarmSetType(), DMSwarmSortGetAccess()
355 @*/
356 PETSC_EXTERN PetscErrorCode DMSwarmSortGetSizes(DM dm,PetscInt *ncells,PetscInt *npoints)
357 {
358   DM_Swarm *swarm = (DM_Swarm*)dm->data;
359 
360   PetscFunctionBegin;
361   if (!swarm->sort_context) {
362     if (ncells) { *ncells = 0; }
363     if (npoints) { *npoints = 0; }
364     PetscFunctionReturn(0);
365   }
366   if (ncells) { *ncells = swarm->sort_context->ncells; }
367   if (npoints) { *npoints = swarm->sort_context->npoints; }
368   PetscFunctionReturn(0);
369 }
370