xref: /petsc/src/dm/impls/plex/pointqueue.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2 
3 PetscErrorCode DMPlexPointQueueCreate(PetscInt size, DMPlexPointQueue *queue) {
4   DMPlexPointQueue q;
5 
6   PetscFunctionBegin;
7   PetscCheck(size >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Queue size %" PetscInt_FMT " must be non-negative", size);
8   PetscCall(PetscCalloc1(1, &q));
9   q->size = size;
10   PetscCall(PetscMalloc1(q->size, &q->points));
11   q->num   = 0;
12   q->front = 0;
13   q->back  = q->size - 1;
14   *queue   = q;
15   PetscFunctionReturn(0);
16 }
17 
18 PetscErrorCode DMPlexPointQueueDestroy(DMPlexPointQueue *queue) {
19   DMPlexPointQueue q = *queue;
20 
21   PetscFunctionBegin;
22   PetscCall(PetscFree(q->points));
23   PetscCall(PetscFree(q));
24   *queue = NULL;
25   PetscFunctionReturn(0);
26 }
27 
28 PetscErrorCode DMPlexPointQueueEnsureSize(DMPlexPointQueue queue) {
29   PetscFunctionBegin;
30   if (queue->num < queue->size) PetscFunctionReturn(0);
31   queue->size *= 2;
32   PetscCall(PetscRealloc(queue->size * sizeof(PetscInt), &queue->points));
33   PetscFunctionReturn(0);
34 }
35 
36 PetscErrorCode DMPlexPointQueueEnqueue(DMPlexPointQueue queue, PetscInt p) {
37   PetscFunctionBegin;
38   PetscCall(DMPlexPointQueueEnsureSize(queue));
39   queue->back                = (queue->back + 1) % queue->size;
40   queue->points[queue->back] = p;
41   ++queue->num;
42   PetscFunctionReturn(0);
43 }
44 
45 PetscErrorCode DMPlexPointQueueDequeue(DMPlexPointQueue queue, PetscInt *p) {
46   PetscFunctionBegin;
47   PetscCheck(queue->num, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot dequeue from an empty queue");
48   *p           = queue->points[queue->front];
49   queue->front = (queue->front + 1) % queue->size;
50   --queue->num;
51   PetscFunctionReturn(0);
52 }
53 
54 PetscErrorCode DMPlexPointQueueFront(DMPlexPointQueue queue, PetscInt *p) {
55   PetscFunctionBegin;
56   PetscCheck(queue->num, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot get the front of an empty queue");
57   *p = queue->points[queue->front];
58   PetscFunctionReturn(0);
59 }
60 
61 PetscErrorCode DMPlexPointQueueBack(DMPlexPointQueue queue, PetscInt *p) {
62   PetscFunctionBegin;
63   PetscCheck(queue->num, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot get the back of an empty queue");
64   *p = queue->points[queue->back];
65   PetscFunctionReturn(0);
66 }
67 
68 PetscBool DMPlexPointQueueEmpty(DMPlexPointQueue queue) {
69   if (!queue->num) return PETSC_TRUE;
70   return PETSC_FALSE;
71 }
72 
73 PetscErrorCode DMPlexPointQueueEmptyCollective(PetscObject obj, DMPlexPointQueue queue, PetscBool *empty) {
74   PetscFunctionBeginHot;
75   *empty = DMPlexPointQueueEmpty(queue);
76   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, empty, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm(obj)));
77   PetscFunctionReturn(0);
78 }
79