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