xref: /petsc/include/petscbt.h (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1652f91beSJacob Faibussowitsch #ifndef PETSCBT_H
226bd1501SBarry Smith #define PETSCBT_H
319fee000SSatish Balay 
4665c2dedSJed Brown #include <petscviewer.h>
56420c192SJed Brown 
6ac09b921SBarry Smith /* SUBMANSEC = Sys */
7ac09b921SBarry Smith 
8f621e05eSBarry Smith /*S
987497f52SBarry Smith      PetscBT - PETSc bitarrays, efficient storage of arrays of boolean values
10f621e05eSBarry Smith 
11f621e05eSBarry Smith      Level: advanced
1219fee000SSatish Balay 
1387497f52SBarry Smith      Notes:
1487497f52SBarry Smith      The following routines do not have their own manual pages
1587497f52SBarry Smith 
1687497f52SBarry Smith .vb
1794bacf5dSBarry Smith      PetscBTCreate(m,&bt)         - creates a bit array with enough room to hold m values
1894bacf5dSBarry Smith      PetscBTDestroy(&bt)          - destroys the bit array
196831982aSBarry Smith      PetscBTMemzero(m,bt)         - zeros the entire bit array (sets all values to false)
206831982aSBarry Smith      PetscBTSet(bt,index)         - sets a particular entry as true
216831982aSBarry Smith      PetscBTClear(bt,index)       - sets a particular entry as false
22f1af5d2fSBarry Smith      PetscBTLookup(bt,index)      - returns the value
23f1af5d2fSBarry Smith      PetscBTLookupSet(bt,index)   - returns the value and then sets it true
2431b932fcSLisandro Dalcin      PetscBTLookupClear(bt,index) - returns the value and then sets it false
256831982aSBarry Smith      PetscBTLength(m)             - returns number of bytes in array with m bits
266831982aSBarry Smith      PetscBTView(m,bt,viewer)     - prints all the entries in a bit array
2787497f52SBarry Smith .ve
28eec0b4cfSBarry Smith 
2987497f52SBarry Smith     PETSc does not check error flags on `PetscBTLookup()`, `PetcBTLookupSet()`, `PetscBTLength()` because error checking
30652f91beSJacob Faibussowitsch     would cost hundreds more cycles then the operation.
316831982aSBarry Smith 
32b9617806SBarry Smith S*/
33521d7252SBarry Smith typedef char *PetscBT;
34ca71c51bSBarry Smith 
35652f91beSJacob Faibussowitsch /* convert an index i to an index suitable for indexing a PetscBT, such that
36652f91beSJacob Faibussowitsch  * bt[PetscBTIndex(i)] returns the i'th value of the bt */
37*d71ae5a4SJacob Faibussowitsch static inline size_t PetscBTIndex_Internal(PetscInt index)
38*d71ae5a4SJacob Faibussowitsch {
39c133393eSBarry Smith   return (size_t)index / PETSC_BITS_PER_BYTE;
40652f91beSJacob Faibussowitsch }
41652f91beSJacob Faibussowitsch 
42*d71ae5a4SJacob Faibussowitsch static inline char PetscBTMask_Internal(PetscInt index)
43*d71ae5a4SJacob Faibussowitsch {
44652f91beSJacob Faibussowitsch   return 1 << index % PETSC_BITS_PER_BYTE;
45652f91beSJacob Faibussowitsch }
46652f91beSJacob Faibussowitsch 
47*d71ae5a4SJacob Faibussowitsch static inline size_t PetscBTLength(PetscInt m)
48*d71ae5a4SJacob Faibussowitsch {
49c133393eSBarry Smith   return (size_t)m / PETSC_BITS_PER_BYTE + 1;
5053b8de81SBarry Smith }
5182502324SSatish Balay 
52*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBTMemzero(PetscInt m, PetscBT array)
53*d71ae5a4SJacob Faibussowitsch {
54652f91beSJacob Faibussowitsch   return PetscArrayzero(array, PetscBTLength(m));
5553b8de81SBarry Smith }
563b71518fSBarry Smith 
57*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBTDestroy(PetscBT *array)
58*d71ae5a4SJacob Faibussowitsch {
59652f91beSJacob Faibussowitsch   return (*array) ? PetscFree(*array) : 0;
6053b8de81SBarry Smith }
61ca71c51bSBarry Smith 
62*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBTCreate(PetscInt m, PetscBT *array)
63*d71ae5a4SJacob Faibussowitsch {
64652f91beSJacob Faibussowitsch   return PetscCalloc1(PetscBTLength(m), array);
6553b8de81SBarry Smith }
6653b8de81SBarry Smith 
67*d71ae5a4SJacob Faibussowitsch static inline char PetscBTLookup(PetscBT array, PetscInt index)
68*d71ae5a4SJacob Faibussowitsch {
69652f91beSJacob Faibussowitsch   return array[PetscBTIndex_Internal(index)] & PetscBTMask_Internal(index);
7053b8de81SBarry Smith }
7153b8de81SBarry Smith 
72*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBTSet(PetscBT array, PetscInt index)
73*d71ae5a4SJacob Faibussowitsch {
74652f91beSJacob Faibussowitsch   PetscFunctionBegin;
75652f91beSJacob Faibussowitsch   array[PetscBTIndex_Internal(index)] |= PetscBTMask_Internal(index);
76652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
7753b8de81SBarry Smith }
7853b8de81SBarry Smith 
79*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBTNegate(PetscBT array, PetscInt index)
80*d71ae5a4SJacob Faibussowitsch {
81652f91beSJacob Faibussowitsch   PetscFunctionBegin;
82652f91beSJacob Faibussowitsch   array[PetscBTIndex_Internal(index)] ^= PetscBTMask_Internal(index);
83652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
8431b932fcSLisandro Dalcin }
8531b932fcSLisandro Dalcin 
86*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBTClear(PetscBT array, PetscInt index)
87*d71ae5a4SJacob Faibussowitsch {
88652f91beSJacob Faibussowitsch   PetscFunctionBegin;
89090438c4SLisandro Dalcin   array[PetscBTIndex_Internal(index)] &= (char)~PetscBTMask_Internal(index);
90652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
9153b8de81SBarry Smith }
9253b8de81SBarry Smith 
93*d71ae5a4SJacob Faibussowitsch static inline char PetscBTLookupSet(PetscBT array, PetscInt index)
94*d71ae5a4SJacob Faibussowitsch {
95652f91beSJacob Faibussowitsch   const char ret = PetscBTLookup(array, index);
969566063dSJacob Faibussowitsch   PetscCallContinue(PetscBTSet(array, index));
97652f91beSJacob Faibussowitsch   return ret;
98652f91beSJacob Faibussowitsch }
99652f91beSJacob Faibussowitsch 
100*d71ae5a4SJacob Faibussowitsch static inline char PetscBTLookupClear(PetscBT array, PetscInt index)
101*d71ae5a4SJacob Faibussowitsch {
102652f91beSJacob Faibussowitsch   const char ret = PetscBTLookup(array, index);
1039566063dSJacob Faibussowitsch   PetscCallContinue(PetscBTClear(array, index));
104652f91beSJacob Faibussowitsch   return ret;
105652f91beSJacob Faibussowitsch }
106652f91beSJacob Faibussowitsch 
107*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscBTView(PetscInt m, const PetscBT bt, PetscViewer viewer)
108*d71ae5a4SJacob Faibussowitsch {
109652f91beSJacob Faibussowitsch   PetscFunctionBegin;
110652f91beSJacob Faibussowitsch   if (m < 1) PetscFunctionReturn(0);
1119566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &viewer));
1129566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
11348a46eb9SPierre Jolivet   for (PetscInt i = 0; i < m; ++i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT " %d\n", i, (int)PetscBTLookup(bt, i)));
1149566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
1159566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
116652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
117652f91beSJacob Faibussowitsch }
118652f91beSJacob Faibussowitsch #endif /* PETSCBT_H */
119