1cc4c1da9SBarry Smith #include <petscdmda.h> /*I "petscdmda.h" I*/
2cc4c1da9SBarry Smith #include <petscdmplex.h> /*I "petscdmplex.h" I*/
3cc4c1da9SBarry Smith #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/
4dc5f5ce9SDave May
sort_CompareSwarmPoint(const void * dataA,const void * dataB)566976f2fSJacob Faibussowitsch static int sort_CompareSwarmPoint(const void *dataA, const void *dataB)
6d71ae5a4SJacob Faibussowitsch {
75627991aSBarry Smith SwarmPoint *pointA = (SwarmPoint *)dataA;
85627991aSBarry Smith SwarmPoint *pointB = (SwarmPoint *)dataB;
9dc5f5ce9SDave May
10dc5f5ce9SDave May if (pointA->cell_index < pointB->cell_index) {
11dc5f5ce9SDave May return -1;
12dc5f5ce9SDave May } else if (pointA->cell_index > pointB->cell_index) {
13dc5f5ce9SDave May return 1;
14dc5f5ce9SDave May } else {
15dc5f5ce9SDave May return 0;
16dc5f5ce9SDave May }
17dc5f5ce9SDave May }
18dc5f5ce9SDave May
DMSwarmSortApplyCellIndexSort(DMSwarmSort ctx)1966976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSortApplyCellIndexSort(DMSwarmSort ctx)
20d71ae5a4SJacob Faibussowitsch {
21dc5f5ce9SDave May PetscFunctionBegin;
22810441c8SPierre Jolivet if (ctx->list) qsort(ctx->list, ctx->npoints, sizeof(SwarmPoint), sort_CompareSwarmPoint);
233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
24dc5f5ce9SDave May }
25dc5f5ce9SDave May
DMSwarmSortCreate(DMSwarmSort * _ctx)2666976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSortCreate(DMSwarmSort *_ctx)
27d71ae5a4SJacob Faibussowitsch {
28dc5f5ce9SDave May DMSwarmSort ctx;
29dc5f5ce9SDave May
30dc5f5ce9SDave May PetscFunctionBegin;
319566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx));
32dc5f5ce9SDave May ctx->isvalid = PETSC_FALSE;
33dc5f5ce9SDave May ctx->ncells = 0;
34dc5f5ce9SDave May ctx->npoints = 0;
359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &ctx->pcell_offsets));
369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &ctx->list));
37dc5f5ce9SDave May *_ctx = ctx;
383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
39dc5f5ce9SDave May }
40dc5f5ce9SDave May
DMSwarmSortSetup(DMSwarmSort ctx,DM dm,PetscInt ncells)4166976f2fSJacob Faibussowitsch static PetscErrorCode DMSwarmSortSetup(DMSwarmSort ctx, DM dm, PetscInt ncells)
42d71ae5a4SJacob Faibussowitsch {
4319307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
44dc5f5ce9SDave May PetscInt *swarm_cellid;
45dc5f5ce9SDave May PetscInt p, npoints;
46dc5f5ce9SDave May PetscInt tmp, c, count;
4719307e5cSMatthew G. Knepley const char *cellid;
48dc5f5ce9SDave May
49dc5f5ce9SDave May PetscFunctionBegin;
503ba16761SJacob Faibussowitsch if (!ctx) PetscFunctionReturn(PETSC_SUCCESS);
513ba16761SJacob Faibussowitsch if (ctx->isvalid) PetscFunctionReturn(PETSC_SUCCESS);
52dc5f5ce9SDave May
539566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMSWARM_Sort, 0, 0, 0, 0));
54dc5f5ce9SDave May /* check the number of cells */
55dc5f5ce9SDave May if (ncells != ctx->ncells) {
569566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(PetscInt) * (ncells + 1), &ctx->pcell_offsets));
57dc5f5ce9SDave May ctx->ncells = ncells;
58dc5f5ce9SDave May }
599566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->pcell_offsets, ctx->ncells + 1));
60dc5f5ce9SDave May
61dc5f5ce9SDave May /* get the number of points */
629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &npoints));
63dc5f5ce9SDave May if (npoints != ctx->npoints) {
649566063dSJacob Faibussowitsch PetscCall(PetscRealloc(sizeof(SwarmPoint) * npoints, &ctx->list));
65dc5f5ce9SDave May ctx->npoints = npoints;
66dc5f5ce9SDave May }
679566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->list, npoints));
68dc5f5ce9SDave May
6919307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
7019307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
7119307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
72dc5f5ce9SDave May for (p = 0; p < ctx->npoints; p++) {
73dc5f5ce9SDave May ctx->list[p].point_index = p;
74dc5f5ce9SDave May ctx->list[p].cell_index = swarm_cellid[p];
75dc5f5ce9SDave May }
7619307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
779566063dSJacob Faibussowitsch PetscCall(DMSwarmSortApplyCellIndexSort(ctx));
78dc5f5ce9SDave May
79dc5f5ce9SDave May /* sum points per cell */
80ad540459SPierre Jolivet for (p = 0; p < ctx->npoints; p++) ctx->pcell_offsets[ctx->list[p].cell_index]++;
81dc5f5ce9SDave May
82dc5f5ce9SDave May /* create offset list */
83dc5f5ce9SDave May count = 0;
84dc5f5ce9SDave May for (c = 0; c < ctx->ncells; c++) {
85dc5f5ce9SDave May tmp = ctx->pcell_offsets[c];
86dc5f5ce9SDave May ctx->pcell_offsets[c] = count;
87dc5f5ce9SDave May count = count + tmp;
88dc5f5ce9SDave May }
89dc5f5ce9SDave May ctx->pcell_offsets[c] = count;
90dc5f5ce9SDave May
91dc5f5ce9SDave May ctx->isvalid = PETSC_TRUE;
929566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMSWARM_Sort, 0, 0, 0, 0));
933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
94dc5f5ce9SDave May }
95dc5f5ce9SDave May
DMSwarmSortDestroy(DMSwarmSort * ctx)96*81ee147aSBarry Smith PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *ctx)
97d71ae5a4SJacob Faibussowitsch {
98*81ee147aSBarry Smith DMSwarmSort ictx;
99dc5f5ce9SDave May
100dc5f5ce9SDave May PetscFunctionBegin;
101*81ee147aSBarry Smith if (!ctx) PetscFunctionReturn(PETSC_SUCCESS);
102*81ee147aSBarry Smith if (!*ctx) PetscFunctionReturn(PETSC_SUCCESS);
103*81ee147aSBarry Smith ictx = *ctx;
104*81ee147aSBarry Smith if (ictx->list) PetscCall(PetscFree(ictx->list));
105*81ee147aSBarry Smith if (ictx->pcell_offsets) PetscCall(PetscFree(ictx->pcell_offsets));
106*81ee147aSBarry Smith PetscCall(PetscFree(ictx));
107*81ee147aSBarry Smith *ctx = NULL;
1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
109dc5f5ce9SDave May }
110dc5f5ce9SDave May
111cc4c1da9SBarry Smith /*@
112dc5f5ce9SDave May DMSwarmSortGetNumberOfPointsPerCell - Returns the number of points in a cell
113dc5f5ce9SDave May
11420f4b53cSBarry Smith Not Collective
115dc5f5ce9SDave May
11660225df5SJacob Faibussowitsch Input Parameters:
11719307e5cSMatthew G. Knepley + sw - a `DMSWARM` objects
11819307e5cSMatthew G. Knepley - cell - the cell number in the cell `DM`
119a3b724e8SBarry Smith
120a3b724e8SBarry Smith Output Parameter:
121a3b724e8SBarry Smith . npoints - the number of points in the cell
122dc5f5ce9SDave May
123dc5f5ce9SDave May Level: advanced
124dc5f5ce9SDave May
125dc5f5ce9SDave May Notes:
12620f4b53cSBarry Smith You must call `DMSwarmSortGetAccess()` before you can call `DMSwarmSortGetNumberOfPointsPerCell()`
127dc5f5ce9SDave May
12820f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSortGetAccess()`, `DMSwarmSortGetPointsPerCell()`
129dc5f5ce9SDave May @*/
DMSwarmSortGetNumberOfPointsPerCell(DM sw,PetscInt cell,PetscInt * npoints)13019307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSortGetNumberOfPointsPerCell(DM sw, PetscInt cell, PetscInt *npoints)
131d71ae5a4SJacob Faibussowitsch {
13219307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
133dc5f5ce9SDave May DMSwarmSort ctx;
134dc5f5ce9SDave May
135dc5f5ce9SDave May PetscFunctionBegin;
13619307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
13719307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetSort(celldm, &ctx));
13819307e5cSMatthew G. Knepley PetscCheck(ctx, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "The DMSwarmSort context has not been created. Must call DMSwarmSortGetAccess() first");
13928b400f6SJacob Faibussowitsch PetscCheck(ctx->isvalid, PETSC_COMM_SELF, PETSC_ERR_USER, "SwarmPointSort container is not valid. Must call DMSwarmSortGetAccess() first");
14019307e5cSMatthew G. Knepley PetscCheck(cell < ctx->ncells, PETSC_COMM_SELF, PETSC_ERR_USER, "Cell index (%" PetscInt_FMT ") is greater than max number of local cells (%" PetscInt_FMT ")", cell, ctx->ncells);
14119307e5cSMatthew G. Knepley PetscCheck(cell >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Cell index (%" PetscInt_FMT ") cannot be negative", cell);
14219307e5cSMatthew G. Knepley *npoints = ctx->pcell_offsets[cell + 1] - ctx->pcell_offsets[cell];
1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
144dc5f5ce9SDave May }
145dc5f5ce9SDave May
146dc5f5ce9SDave May /*@C
147dc5f5ce9SDave May DMSwarmSortGetPointsPerCell - Creates an array of point indices for all points in a cell
148dc5f5ce9SDave May
14920f4b53cSBarry Smith Not Collective
150dc5f5ce9SDave May
15160225df5SJacob Faibussowitsch Input Parameters:
15219307e5cSMatthew G. Knepley + sw - a `DMSWARM` object
15319307e5cSMatthew G. Knepley . cell - the cell number in the cell `DM`
154dc5f5ce9SDave May . npoints - the number of points in the cell
155da81f932SPierre Jolivet - pidlist - array of the indices identifying all points in cell e
156dc5f5ce9SDave May
157dc5f5ce9SDave May Level: advanced
158dc5f5ce9SDave May
159fe11354eSMatthew G. Knepley Note:
160fe11354eSMatthew G. Knepley You must call `DMSwarmSortGetAccess()` before you can call `DMSwarmSortGetPointsPerCell()`, and call `DMSwarmRestorePointsPerCell()` afterwards
1615627991aSBarry Smith
162fe11354eSMatthew G. Knepley .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmRestorePointsPerCell()`, `DMSwarmSortGetAccess()`, `DMSwarmSortGetNumberOfPointsPerCell()`
163dc5f5ce9SDave May @*/
DMSwarmSortGetPointsPerCell(DM sw,PetscInt cell,PetscInt * npoints,PetscInt ** pidlist)16419307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSortGetPointsPerCell(DM sw, PetscInt cell, PetscInt *npoints, PetscInt **pidlist)
165d71ae5a4SJacob Faibussowitsch {
16619307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
16719307e5cSMatthew G. Knepley PetscInt pid, pid_unsorted;
168dc5f5ce9SDave May DMSwarmSort ctx;
169dc5f5ce9SDave May
170dc5f5ce9SDave May PetscFunctionBegin;
17119307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
17219307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetSort(celldm, &ctx));
17319307e5cSMatthew G. Knepley PetscCheck(ctx, PetscObjectComm((PetscObject)sw), PETSC_ERR_USER, "The DMSwarmSort context has not been created. Must call DMSwarmSortGetAccess() first");
17419307e5cSMatthew G. Knepley PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cell, npoints));
17519307e5cSMatthew G. Knepley PetscCall(DMGetWorkArray(sw, *npoints, MPIU_SCALAR, pidlist));
17619307e5cSMatthew G. Knepley for (PetscInt p = 0; p < *npoints; ++p) {
17719307e5cSMatthew G. Knepley pid = ctx->pcell_offsets[cell] + p;
178dc5f5ce9SDave May pid_unsorted = ctx->list[pid].point_index;
179fe11354eSMatthew G. Knepley (*pidlist)[p] = pid_unsorted;
180dc5f5ce9SDave May }
181fe11354eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
182fe11354eSMatthew G. Knepley }
183fe11354eSMatthew G. Knepley
184fe11354eSMatthew G. Knepley /*@C
185fe11354eSMatthew G. Knepley DMSwarmSortRestorePointsPerCell - Restores an array of point indices for all points in a cell
186fe11354eSMatthew G. Knepley
187fe11354eSMatthew G. Knepley Not Collective
188fe11354eSMatthew G. Knepley
189fe11354eSMatthew G. Knepley Input Parameters:
190fe11354eSMatthew G. Knepley + dm - a `DMSWARM` object
191fe11354eSMatthew G. Knepley . e - the index of the cell
192fe11354eSMatthew G. Knepley . npoints - the number of points in the cell
193fe11354eSMatthew G. Knepley - pidlist - array of the indices identifying all points in cell e
194fe11354eSMatthew G. Knepley
195fe11354eSMatthew G. Knepley Level: advanced
196fe11354eSMatthew G. Knepley
197fe11354eSMatthew G. Knepley Note:
198fe11354eSMatthew G. Knepley You must call `DMSwarmSortGetAccess()` and `DMSwarmSortGetPointsPerCell()` before you can call `DMSwarmSortRestorePointsPerCell()`
199fe11354eSMatthew G. Knepley
200fe11354eSMatthew G. Knepley .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSortGetPointsPerCell()`, `DMSwarmSortGetAccess()`, `DMSwarmSortGetNumberOfPointsPerCell()`
201fe11354eSMatthew G. Knepley @*/
DMSwarmSortRestorePointsPerCell(DM dm,PetscInt e,PetscInt * npoints,PetscInt ** pidlist)202fe11354eSMatthew G. Knepley PetscErrorCode DMSwarmSortRestorePointsPerCell(DM dm, PetscInt e, PetscInt *npoints, PetscInt **pidlist)
203fe11354eSMatthew G. Knepley {
204fe11354eSMatthew G. Knepley PetscFunctionBegin;
205fe11354eSMatthew G. Knepley PetscCall(DMRestoreWorkArray(dm, *npoints, MPIU_SCALAR, pidlist));
2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
207dc5f5ce9SDave May }
208dc5f5ce9SDave May
209cc4c1da9SBarry Smith /*@
21020f4b53cSBarry Smith DMSwarmSortGetAccess - Setups up a `DMSWARM` point sort context for efficient traversal of points within a cell
211dc5f5ce9SDave May
21220f4b53cSBarry Smith Not Collective
213dc5f5ce9SDave May
21460225df5SJacob Faibussowitsch Input Parameter:
21519307e5cSMatthew G. Knepley . sw - a `DMSWARM` object
216dc5f5ce9SDave May
217dc5f5ce9SDave May Level: advanced
218dc5f5ce9SDave May
21920f4b53cSBarry Smith Notes:
22020f4b53cSBarry Smith Calling `DMSwarmSortGetAccess()` creates a list which enables easy identification of all points contained in a
22120f4b53cSBarry Smith given cell. This method does not explicitly sort the data within the `DMSWARM` based on the cell index associated
22220f4b53cSBarry Smith with a `DMSWARM` point.
22320f4b53cSBarry Smith
22420f4b53cSBarry Smith The sort context is valid only for the `DMSWARM` points defined at the time when `DMSwarmSortGetAccess()` was called.
22520f4b53cSBarry Smith For example, suppose the swarm contained NP points when `DMSwarmSortGetAccess()` was called. If the user subsequently
22620f4b53cSBarry Smith adds 10 additional points to the swarm, the sort context is still valid, but only for the first NP points.
22720f4b53cSBarry Smith The indices associated with the 10 new additional points will not be contained within the sort context.
22820f4b53cSBarry Smith This means that the user can still safely perform queries via `DMSwarmSortGetPointsPerCell()` and
22920f4b53cSBarry Smith `DMSwarmSortGetPointsPerCell()`, however the results return will be based on the first NP points.
23020f4b53cSBarry Smith
23120f4b53cSBarry Smith If any` DMSWARM` re-sizing method is called after `DMSwarmSortGetAccess()` which modifies any of the first NP entries
23220f4b53cSBarry Smith in the `DMSWARM`, the sort context will become invalid. Currently there are no guards to prevent the user from
23320f4b53cSBarry Smith invalidating the sort context. For this reason, we highly recommend you do not use `DMSwarmRemovePointAtIndex()` in
23420f4b53cSBarry Smith between calls to `DMSwarmSortGetAccess()` and `DMSwarmSortRestoreAccess()`.
23520f4b53cSBarry Smith
23620f4b53cSBarry Smith To facilitate safe removal of points using the sort context, we suggest a "two pass" strategy in which the
23720f4b53cSBarry Smith first pass "marks" points for removal, and the second pass actually removes the points from the `DMSWARM`
23820f4b53cSBarry Smith
23920f4b53cSBarry Smith You must call `DMSwarmSortGetAccess()` before you can call `DMSwarmSortGetPointsPerCell()` or `DMSwarmSortGetNumberOfPointsPerCell()`
24020f4b53cSBarry Smith
24120f4b53cSBarry Smith The sort context may become invalid if any re-sizing methods are applied which alter the first NP points
24220f4b53cSBarry Smith within swarm at the time `DMSwarmSortGetAccess()` was called.
24320f4b53cSBarry Smith
24420f4b53cSBarry Smith You must call `DMSwarmSortRestoreAccess()` when you no longer need access to the sort context
24520f4b53cSBarry Smith
24620f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSortRestoreAccess()`
247dc5f5ce9SDave May @*/
DMSwarmSortGetAccess(DM sw)24819307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSortGetAccess(DM sw)
249d71ae5a4SJacob Faibussowitsch {
25019307e5cSMatthew G. Knepley DM dm;
25119307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
25219307e5cSMatthew G. Knepley DMSwarmSort ctx;
25319307e5cSMatthew G. Knepley PetscInt ncells = 0;
254dc5f5ce9SDave May PetscBool isda, isplex, isshell;
255dc5f5ce9SDave May
256dc5f5ce9SDave May PetscFunctionBegin;
25719307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
25819307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetSort(celldm, &ctx));
25919307e5cSMatthew G. Knepley if (!ctx) {
26019307e5cSMatthew G. Knepley PetscCall(DMSwarmSortCreate(&ctx));
26119307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMSetSort(celldm, ctx));
26219307e5cSMatthew G. Knepley }
263dc5f5ce9SDave May
264dc5f5ce9SDave May /* get the number of cells */
26519307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDM(sw, &dm));
26619307e5cSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
26719307e5cSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isplex));
26819307e5cSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSHELL, &isshell));
269dc5f5ce9SDave May if (isda) {
270dc5f5ce9SDave May const PetscInt *element;
27119307e5cSMatthew G. Knepley PetscInt nel, npe;
272dc5f5ce9SDave May
27319307e5cSMatthew G. Knepley PetscCall(DMDAGetElements(dm, &nel, &npe, &element));
274dc5f5ce9SDave May ncells = nel;
27519307e5cSMatthew G. Knepley PetscCall(DMDARestoreElements(dm, &nel, &npe, &element));
276dc5f5ce9SDave May } else if (isplex) {
277dc5f5ce9SDave May PetscInt ps, pe;
278dc5f5ce9SDave May
27919307e5cSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &ps, &pe));
280dc5f5ce9SDave May ncells = pe - ps;
281dc5f5ce9SDave May } else if (isshell) {
282dc5f5ce9SDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *);
283dc5f5ce9SDave May
28419307e5cSMatthew G. Knepley PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells));
285dc5f5ce9SDave May if (method_DMShellGetNumberOfCells) {
28619307e5cSMatthew G. Knepley PetscCall(method_DMShellGetNumberOfCells(dm, &ncells));
2879371c9d4SSatish Balay } else
28819307e5cSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)sw), 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);");
28919307e5cSMatthew G. Knepley } else SETERRQ(PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL");
290dc5f5ce9SDave May
291dc5f5ce9SDave May /* setup */
29219307e5cSMatthew G. Knepley PetscCall(DMSwarmSortSetup(ctx, sw, ncells));
2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
294dc5f5ce9SDave May }
295dc5f5ce9SDave May
296cc4c1da9SBarry Smith /*@
297cc4c1da9SBarry Smith DMSwarmSortRestoreAccess - Invalidates the `DMSWARM` point sorting context previously computed with `DMSwarmSortGetAccess()`
298dc5f5ce9SDave May
29920f4b53cSBarry Smith Not Collective
300dc5f5ce9SDave May
30160225df5SJacob Faibussowitsch Input Parameter:
30219307e5cSMatthew G. Knepley . sw - a `DMSWARM` object
303dc5f5ce9SDave May
304dc5f5ce9SDave May Level: advanced
305dc5f5ce9SDave May
306dc5f5ce9SDave May Note:
30720f4b53cSBarry Smith You must call `DMSwarmSortGetAccess()` before calling `DMSwarmSortRestoreAccess()`
308dc5f5ce9SDave May
30920f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSortGetAccess()`
310dc5f5ce9SDave May @*/
DMSwarmSortRestoreAccess(DM sw)31119307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSortRestoreAccess(DM sw)
312d71ae5a4SJacob Faibussowitsch {
31319307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
31419307e5cSMatthew G. Knepley DMSwarmSort ctx;
315dc5f5ce9SDave May
316dc5f5ce9SDave May PetscFunctionBegin;
31719307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
31819307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetSort(celldm, &ctx));
31919307e5cSMatthew G. Knepley if (!ctx) PetscFunctionReturn(PETSC_SUCCESS);
32019307e5cSMatthew G. Knepley PetscCheck(ctx->isvalid, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "You must call DMSwarmSortGetAccess() before calling DMSwarmSortRestoreAccess()");
32119307e5cSMatthew G. Knepley ctx->isvalid = PETSC_FALSE;
3223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
323dc5f5ce9SDave May }
324b799feefSDave May
325cc4c1da9SBarry Smith /*@
32620f4b53cSBarry Smith DMSwarmSortGetIsValid - Gets the isvalid flag associated with a `DMSWARM` point sorting context
327b799feefSDave May
32820f4b53cSBarry Smith Not Collective
329b799feefSDave May
33060225df5SJacob Faibussowitsch Input Parameter:
33119307e5cSMatthew G. Knepley . sw - a `DMSWARM` object
332b799feefSDave May
33360225df5SJacob Faibussowitsch Output Parameter:
334b799feefSDave May . isvalid - flag indicating whether the sort context is up-to-date
335b799feefSDave May
336b799feefSDave May Level: advanced
337b799feefSDave May
33820f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSortGetAccess()`
339b799feefSDave May @*/
DMSwarmSortGetIsValid(DM sw,PetscBool * isvalid)34019307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSortGetIsValid(DM sw, PetscBool *isvalid)
341d71ae5a4SJacob Faibussowitsch {
34219307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
34319307e5cSMatthew G. Knepley DMSwarmSort ctx;
344b799feefSDave May
345b799feefSDave May PetscFunctionBegin;
34619307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
34719307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetSort(celldm, &ctx));
34819307e5cSMatthew G. Knepley if (!ctx) {
349b799feefSDave May *isvalid = PETSC_FALSE;
3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
351b799feefSDave May }
35219307e5cSMatthew G. Knepley *isvalid = ctx->isvalid;
3533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
354b799feefSDave May }
355b799feefSDave May
356cc4c1da9SBarry Smith /*@
35720f4b53cSBarry Smith DMSwarmSortGetSizes - Gets the sizes associated with a `DMSWARM` point sorting context
358b799feefSDave May
35920f4b53cSBarry Smith Not Collective
360b799feefSDave May
36160225df5SJacob Faibussowitsch Input Parameter:
36219307e5cSMatthew G. Knepley . sw - a `DMSWARM` object
363b799feefSDave May
36460225df5SJacob Faibussowitsch Output Parameters:
36520f4b53cSBarry Smith + ncells - number of cells within the sort context (pass `NULL` to ignore)
36620f4b53cSBarry Smith - npoints - number of points used to create the sort context (pass `NULL` to ignore)
367b799feefSDave May
368b799feefSDave May Level: advanced
369b799feefSDave May
37020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSortGetAccess()`
371b799feefSDave May @*/
DMSwarmSortGetSizes(DM sw,PetscInt * ncells,PetscInt * npoints)37219307e5cSMatthew G. Knepley PetscErrorCode DMSwarmSortGetSizes(DM sw, PetscInt *ncells, PetscInt *npoints)
373d71ae5a4SJacob Faibussowitsch {
37419307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
37519307e5cSMatthew G. Knepley DMSwarmSort ctx;
376b799feefSDave May
377b799feefSDave May PetscFunctionBegin;
37819307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(sw, &celldm));
37919307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetSort(celldm, &ctx));
38019307e5cSMatthew G. Knepley if (!ctx) {
3815627991aSBarry Smith if (ncells) *ncells = 0;
3825627991aSBarry Smith if (npoints) *npoints = 0;
3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
384b799feefSDave May }
38519307e5cSMatthew G. Knepley if (ncells) *ncells = ctx->ncells;
38619307e5cSMatthew G. Knepley if (npoints) *npoints = ctx->npoints;
3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
388b799feefSDave May }
389