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