xref: /petsc/include/petscbt.h (revision fe998a80077c9ee0917a39496df43fc256e1b478)
1 
2 #if !defined(__PETSCBT_H)
3 #define __PETSCBT_H
4 
5 #include <petscconf.h>
6 #include <petscviewer.h>
7 
8 /*S
9      PetscBT - PETSc bitarrays
10 
11      Level: advanced
12 
13      PetscBTCreate(m,&bt)         - creates a bit array with enough room to hold m values
14      PetscBTDestroy(&bt)          - destroys the bit array
15      PetscBTMemzero(m,bt)         - zeros the entire bit array (sets all values to false)
16      PetscBTSet(bt,index)         - sets a particular entry as true
17      PetscBTClear(bt,index)       - sets a particular entry as false
18      PetscBTLookup(bt,index)      - returns the value
19      PetscBTLookupSet(bt,index)   - returns the value and then sets it true
20      PetscBTLookupClear(bt,index) - returns the value and then sets it false
21      PetscBTLength(m)             - returns number of bytes in array with m bits
22      PetscBTView(m,bt,viewer)     - prints all the entries in a bit array
23 
24     We do not currently check error flags on PetscBTSet(), PetscBTClear(), PetscBTLookup(),
25     PetcBTLookupSet(), PetscBTLength() cause error checking would cost hundreds more cycles then
26     the operation.
27 
28 S*/
29 typedef char* PetscBT;
30 
31 
32 PETSC_STATIC_INLINE PetscInt  PetscBTLength(PetscInt m)
33 {
34   return  ((m)/PETSC_BITS_PER_BYTE+1);
35 }
36 
37 PETSC_STATIC_INLINE PetscErrorCode PetscBTMemzero(PetscInt m,PetscBT array)
38 {
39   return PetscMemzero(array,sizeof(char)*((m)/PETSC_BITS_PER_BYTE+1));
40 }
41 
42 PETSC_STATIC_INLINE PetscErrorCode PetscBTDestroy(PetscBT *array)
43 {
44   return PetscFree(*array);
45 }
46 
47 PETSC_STATIC_INLINE char PetscBTLookup(PetscBT array,PetscInt index)
48 {
49   char      BT_mask,BT_c;
50   PetscInt  BT_idx;
51 
52  return  (BT_idx        = (index)/PETSC_BITS_PER_BYTE,
53           BT_c          = array[BT_idx],
54           BT_mask       = (char)1 << ((index)%PETSC_BITS_PER_BYTE),
55           BT_c & BT_mask);
56 }
57 
58 PETSC_STATIC_INLINE PetscErrorCode PetscBTView(PetscInt m,const PetscBT bt,PetscViewer viewer)
59 {
60   PetscInt       i;
61   PetscErrorCode ierr;
62 
63   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PETSC_COMM_SELF,&viewer);CHKERRQ(ierr);}
64   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
65   for (i=0; i<m; i++) {
66     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D %d\n",i,PetscBTLookup(bt,i));CHKERRQ(ierr);
67   }
68   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
69   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
70   return 0;
71 }
72 
73 PETSC_STATIC_INLINE PetscErrorCode PetscBTCreate(PetscInt m,PetscBT *array)
74 {
75   return (PetscMalloc1(m/PETSC_BITS_PER_BYTE+1,array) || PetscBTMemzero(m,*array));
76 }
77 
78 PETSC_STATIC_INLINE char PetscBTLookupSet(PetscBT array,PetscInt index)
79 {
80   char      BT_mask,BT_c;
81   PetscInt  BT_idx;
82 
83   return (BT_idx        = (index)/PETSC_BITS_PER_BYTE,
84           BT_c          = array[BT_idx],
85           BT_mask       = (char)1 << ((index)%PETSC_BITS_PER_BYTE),
86           array[BT_idx] = BT_c | BT_mask,
87           BT_c & BT_mask);
88 }
89 
90 PETSC_STATIC_INLINE PetscErrorCode PetscBTSet(PetscBT array,PetscInt index)
91 {
92   char      BT_mask,BT_c;
93   PetscInt  BT_idx;
94 
95   BT_idx        = (index)/PETSC_BITS_PER_BYTE;
96   BT_c          = array[BT_idx];
97   BT_mask       = (char)1 << ((index)%PETSC_BITS_PER_BYTE);
98   array[BT_idx] = BT_c | BT_mask;
99   return 0;
100 }
101 
102 PETSC_STATIC_INLINE PetscErrorCode PetscBTNegate(PetscBT array,PetscInt index)
103 {
104   const PetscInt BT_idx  = (index)/PETSC_BITS_PER_BYTE;
105   const char     BT_mask = (char)1 << ((index)%PETSC_BITS_PER_BYTE);
106 
107   array[BT_idx] ^= BT_mask;
108   return 0;
109 }
110 
111 PETSC_STATIC_INLINE char PetscBTLookupClear(PetscBT array,PetscInt index)
112 {
113   char      BT_mask,BT_c;
114   PetscInt  BT_idx;
115 
116   return (BT_idx        = (index)/PETSC_BITS_PER_BYTE,
117           BT_c          = array[BT_idx],
118           BT_mask       = (char)1 << ((index)%PETSC_BITS_PER_BYTE),
119           array[BT_idx] = BT_c & (~BT_mask),
120           BT_c & BT_mask);
121 }
122 
123 PETSC_STATIC_INLINE PetscErrorCode PetscBTClear(PetscBT array,PetscInt index)
124 {
125   char      BT_mask,BT_c;
126   PetscInt  BT_idx;
127 
128   BT_idx        = (index)/PETSC_BITS_PER_BYTE;
129   BT_c          = array[BT_idx];
130   BT_mask       = (char)1 << ((index)%PETSC_BITS_PER_BYTE);
131   array[BT_idx] = BT_c & (~BT_mask);
132  return 0;
133 }
134 
135 
136 #endif
137