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