xref: /petsc/include/petscbt.h (revision b665b14e20d08dc58a3f47e0addbfcd5129cdb60)
1 #ifndef PETSCBT_H
2 #define PETSCBT_H
3 
4 #include <petscsystypes.h>
5 #include <petscviewertypes.h>
6 #include <petscstring.h>
7 
8 /* SUBMANSEC = Sys */
9 
10 /* convert an index i to an index suitable for indexing a PetscBT, such that
11  * bt[PetscBTIndex(i)] returns the i'th value of the bt */
12 static inline size_t PetscBTIndex_Internal(PetscInt index)
13 {
14   return (size_t)index / PETSC_BITS_PER_BYTE;
15 }
16 
17 static inline char PetscBTMask_Internal(PetscInt index)
18 {
19   return (char)(1 << index % PETSC_BITS_PER_BYTE);
20 }
21 
22 static inline size_t PetscBTLength(PetscInt m)
23 {
24   return (size_t)m / PETSC_BITS_PER_BYTE + 1;
25 }
26 
27 static inline PetscErrorCode PetscBTMemzero(PetscInt m, PetscBT array)
28 {
29   return PetscArrayzero(array, PetscBTLength(m));
30 }
31 
32 static inline PetscErrorCode PetscBTDestroy(PetscBT *array)
33 {
34   return (*array) ? PetscFree(*array) : PETSC_SUCCESS;
35 }
36 
37 static inline PetscErrorCode PetscBTCreate(PetscInt m, PetscBT *array)
38 {
39   return PetscCalloc1(PetscBTLength(m), array);
40 }
41 
42 static inline char PetscBTLookup(PetscBT array, PetscInt index)
43 {
44   return array[PetscBTIndex_Internal(index)] & PetscBTMask_Internal(index);
45 }
46 
47 static inline PetscErrorCode PetscBTSet(PetscBT array, PetscInt index)
48 {
49   PetscFunctionBegin;
50   array[PetscBTIndex_Internal(index)] |= PetscBTMask_Internal(index);
51   PetscFunctionReturn(PETSC_SUCCESS);
52 }
53 
54 static inline PetscErrorCode PetscBTNegate(PetscBT array, PetscInt index)
55 {
56   PetscFunctionBegin;
57   array[PetscBTIndex_Internal(index)] ^= PetscBTMask_Internal(index);
58   PetscFunctionReturn(PETSC_SUCCESS);
59 }
60 
61 static inline PetscErrorCode PetscBTClear(PetscBT array, PetscInt index)
62 {
63   PetscFunctionBegin;
64   array[PetscBTIndex_Internal(index)] &= (char)~PetscBTMask_Internal(index);
65   PetscFunctionReturn(PETSC_SUCCESS);
66 }
67 
68 static inline char PetscBTLookupSet(PetscBT array, PetscInt index)
69 {
70   const char ret = PetscBTLookup(array, index);
71   PetscCallContinue(PetscBTSet(array, index));
72   return ret;
73 }
74 
75 static inline char PetscBTLookupClear(PetscBT array, PetscInt index)
76 {
77   const char ret = PetscBTLookup(array, index);
78   PetscCallContinue(PetscBTClear(array, index));
79   return ret;
80 }
81 
82 PETSC_EXTERN PetscErrorCode PetscBTView(PetscInt, const PetscBT, PetscViewer);
83 
84 #endif /* PETSCBT_H */
85