Actual source code: petscbt.h
1: #ifndef PETSCBT_H
2: #define PETSCBT_H
4: #include <petscviewer.h>
6: /* SUBMANSEC = Sys */
8: /*S
9: PetscBT - PETSc bitarrays, efficient storage of arrays of boolean values
11: Level: advanced
13: Notes:
14: The following routines do not have their own manual pages
16: .vb
17: PetscBTCreate(m,&bt) - creates a bit array with enough room to hold m values
18: PetscBTDestroy(&bt) - destroys the bit array
19: PetscBTMemzero(m,bt) - zeros the entire bit array (sets all values to false)
20: PetscBTSet(bt,index) - sets a particular entry as true
21: PetscBTClear(bt,index) - sets a particular entry as false
22: PetscBTLookup(bt,index) - returns the value
23: PetscBTLookupSet(bt,index) - returns the value and then sets it true
24: PetscBTLookupClear(bt,index) - returns the value and then sets it false
25: PetscBTLength(m) - returns number of bytes in array with m bits
26: PetscBTView(m,bt,viewer) - prints all the entries in a bit array
27: .ve
29: PETSc does not check error flags on `PetscBTLookup()`, `PetcBTLookupSet()`, `PetscBTLength()` because error checking
30: would cost hundreds more cycles then the operation.
32: S*/
33: typedef char *PetscBT;
35: /* convert an index i to an index suitable for indexing a PetscBT, such that
36: * bt[PetscBTIndex(i)] returns the i'th value of the bt */
37: static inline size_t PetscBTIndex_Internal(PetscInt index)
38: {
39: return (size_t)index / PETSC_BITS_PER_BYTE;
40: }
42: static inline char PetscBTMask_Internal(PetscInt index)
43: {
44: return (char)(1 << index % PETSC_BITS_PER_BYTE);
45: }
47: static inline size_t PetscBTLength(PetscInt m)
48: {
49: return (size_t)m / PETSC_BITS_PER_BYTE + 1;
50: }
52: static inline PetscErrorCode PetscBTMemzero(PetscInt m, PetscBT array)
53: {
54: return PetscArrayzero(array, PetscBTLength(m));
55: }
57: static inline PetscErrorCode PetscBTDestroy(PetscBT *array)
58: {
59: return (*array) ? PetscFree(*array) : 0;
60: }
62: static inline PetscErrorCode PetscBTCreate(PetscInt m, PetscBT *array)
63: {
64: return PetscCalloc1(PetscBTLength(m), array);
65: }
67: static inline char PetscBTLookup(PetscBT array, PetscInt index)
68: {
69: return array[PetscBTIndex_Internal(index)] & PetscBTMask_Internal(index);
70: }
72: static inline PetscErrorCode PetscBTSet(PetscBT array, PetscInt index)
73: {
74: array[PetscBTIndex_Internal(index)] |= PetscBTMask_Internal(index);
75: return 0;
76: }
78: static inline PetscErrorCode PetscBTNegate(PetscBT array, PetscInt index)
79: {
80: array[PetscBTIndex_Internal(index)] ^= PetscBTMask_Internal(index);
81: return 0;
82: }
84: static inline PetscErrorCode PetscBTClear(PetscBT array, PetscInt index)
85: {
86: array[PetscBTIndex_Internal(index)] &= (char)~PetscBTMask_Internal(index);
87: return 0;
88: }
90: static inline char PetscBTLookupSet(PetscBT array, PetscInt index)
91: {
92: const char ret = PetscBTLookup(array, index);
93: PetscBTSet(array, index);
94: return ret;
95: }
97: static inline char PetscBTLookupClear(PetscBT array, PetscInt index)
98: {
99: const char ret = PetscBTLookup(array, index);
100: PetscBTClear(array, index);
101: return ret;
102: }
104: static inline PetscErrorCode PetscBTView(PetscInt m, const PetscBT bt, PetscViewer viewer)
105: {
106: if (m < 1) return 0;
107: if (!viewer) PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &viewer);
108: PetscViewerASCIIPushSynchronized(viewer);
109: for (PetscInt i = 0; i < m; ++i) PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT " %d\n", i, (int)PetscBTLookup(bt, i));
110: PetscViewerFlush(viewer);
111: PetscViewerASCIIPopSynchronized(viewer);
112: return 0;
113: }
114: #endif /* PETSCBT_H */