xref: /petsc/include/petscbt.h (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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*9371c9d4SSatish Balay static inline size_t PetscBTIndex_Internal(PetscInt index) {
38c133393eSBarry Smith   return (size_t)index / PETSC_BITS_PER_BYTE;
39652f91beSJacob Faibussowitsch }
40652f91beSJacob Faibussowitsch 
41*9371c9d4SSatish Balay static inline char PetscBTMask_Internal(PetscInt index) {
42652f91beSJacob Faibussowitsch   return 1 << index % PETSC_BITS_PER_BYTE;
43652f91beSJacob Faibussowitsch }
44652f91beSJacob Faibussowitsch 
45*9371c9d4SSatish Balay static inline size_t PetscBTLength(PetscInt m) {
46c133393eSBarry Smith   return (size_t)m / PETSC_BITS_PER_BYTE + 1;
4753b8de81SBarry Smith }
4882502324SSatish Balay 
49*9371c9d4SSatish Balay static inline PetscErrorCode PetscBTMemzero(PetscInt m, PetscBT array) {
50652f91beSJacob Faibussowitsch   return PetscArrayzero(array, PetscBTLength(m));
5153b8de81SBarry Smith }
523b71518fSBarry Smith 
53*9371c9d4SSatish Balay static inline PetscErrorCode PetscBTDestroy(PetscBT *array) {
54652f91beSJacob Faibussowitsch   return (*array) ? PetscFree(*array) : 0;
5553b8de81SBarry Smith }
56ca71c51bSBarry Smith 
57*9371c9d4SSatish Balay static inline PetscErrorCode PetscBTCreate(PetscInt m, PetscBT *array) {
58652f91beSJacob Faibussowitsch   return PetscCalloc1(PetscBTLength(m), array);
5953b8de81SBarry Smith }
6053b8de81SBarry Smith 
61*9371c9d4SSatish Balay static inline char PetscBTLookup(PetscBT array, PetscInt index) {
62652f91beSJacob Faibussowitsch   return array[PetscBTIndex_Internal(index)] & PetscBTMask_Internal(index);
6353b8de81SBarry Smith }
6453b8de81SBarry Smith 
65*9371c9d4SSatish Balay static inline PetscErrorCode PetscBTSet(PetscBT array, PetscInt index) {
66652f91beSJacob Faibussowitsch   PetscFunctionBegin;
67652f91beSJacob Faibussowitsch   array[PetscBTIndex_Internal(index)] |= PetscBTMask_Internal(index);
68652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
6953b8de81SBarry Smith }
7053b8de81SBarry Smith 
71*9371c9d4SSatish Balay static inline PetscErrorCode PetscBTNegate(PetscBT array, PetscInt index) {
72652f91beSJacob Faibussowitsch   PetscFunctionBegin;
73652f91beSJacob Faibussowitsch   array[PetscBTIndex_Internal(index)] ^= PetscBTMask_Internal(index);
74652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
7531b932fcSLisandro Dalcin }
7631b932fcSLisandro Dalcin 
77*9371c9d4SSatish Balay static inline PetscErrorCode PetscBTClear(PetscBT array, PetscInt index) {
78652f91beSJacob Faibussowitsch   PetscFunctionBegin;
79090438c4SLisandro Dalcin   array[PetscBTIndex_Internal(index)] &= (char)~PetscBTMask_Internal(index);
80652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
8153b8de81SBarry Smith }
8253b8de81SBarry Smith 
83*9371c9d4SSatish Balay static inline char PetscBTLookupSet(PetscBT array, PetscInt index) {
84652f91beSJacob Faibussowitsch   const char ret = PetscBTLookup(array, index);
859566063dSJacob Faibussowitsch   PetscCallContinue(PetscBTSet(array, index));
86652f91beSJacob Faibussowitsch   return ret;
87652f91beSJacob Faibussowitsch }
88652f91beSJacob Faibussowitsch 
89*9371c9d4SSatish Balay static inline char PetscBTLookupClear(PetscBT array, PetscInt index) {
90652f91beSJacob Faibussowitsch   const char ret = PetscBTLookup(array, index);
919566063dSJacob Faibussowitsch   PetscCallContinue(PetscBTClear(array, index));
92652f91beSJacob Faibussowitsch   return ret;
93652f91beSJacob Faibussowitsch }
94652f91beSJacob Faibussowitsch 
95*9371c9d4SSatish Balay static inline PetscErrorCode PetscBTView(PetscInt m, const PetscBT bt, PetscViewer viewer) {
96652f91beSJacob Faibussowitsch   PetscFunctionBegin;
97652f91beSJacob Faibussowitsch   if (m < 1) PetscFunctionReturn(0);
989566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &viewer));
999566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
100*9371c9d4SSatish Balay   for (PetscInt i = 0; i < m; ++i) { PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT " %d\n", i, (int)PetscBTLookup(bt, i))); }
1019566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
1029566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
103652f91beSJacob Faibussowitsch   PetscFunctionReturn(0);
104652f91beSJacob Faibussowitsch }
105652f91beSJacob Faibussowitsch #endif /* PETSCBT_H */
106