1a1db5e05SBarry Smith #define PETSCKSP_DLL 2a1db5e05SBarry Smith 3a1db5e05SBarry Smith /********************************bit_mask.c************************************ 4a1db5e05SBarry Smith 5a1db5e05SBarry Smith Author: Henry M. Tufo III 6a1db5e05SBarry Smith 7a1db5e05SBarry Smith e-mail: hmt@cs.brown.edu 8a1db5e05SBarry Smith 9a1db5e05SBarry Smith snail-mail: 10a1db5e05SBarry Smith Division of Applied Mathematics 11a1db5e05SBarry Smith Brown University 12a1db5e05SBarry Smith Providence, RI 02912 13a1db5e05SBarry Smith 14a1db5e05SBarry Smith Last Modification: 15a1db5e05SBarry Smith 11.21.97 16a1db5e05SBarry Smith *********************************bit_mask.c***********************************/ 17a1db5e05SBarry Smith #include "../src/ksp/pc/impls/tfs/tfs.h" 18a1db5e05SBarry Smith 19a1db5e05SBarry Smith 20a1db5e05SBarry Smith /*********************************bit_mask.c***********************************/ 21a1db5e05SBarry Smith PetscErrorCode bm_to_proc( char *ptr, PetscInt p_mask, PetscInt *msg_list) 22a1db5e05SBarry Smith { 23a1db5e05SBarry Smith PetscInt i, tmp; 24a1db5e05SBarry Smith 25a1db5e05SBarry Smith PetscFunctionBegin; 26a1db5e05SBarry Smith if (msg_list) 27a1db5e05SBarry Smith { 28a1db5e05SBarry Smith /* low to high */ 29a1db5e05SBarry Smith ptr+=(p_mask-1); 30a1db5e05SBarry Smith for (i=p_mask-1;i>=0;i--) 31a1db5e05SBarry Smith { 32a1db5e05SBarry Smith tmp = BYTE*(p_mask-i-1); 33a1db5e05SBarry Smith if (*ptr&BIT_0) 34a1db5e05SBarry Smith {*msg_list = tmp; msg_list++;} 35a1db5e05SBarry Smith if (*ptr&BIT_1) 36a1db5e05SBarry Smith {*msg_list = tmp+1; msg_list++;} 37a1db5e05SBarry Smith if (*ptr&BIT_2) 38a1db5e05SBarry Smith {*msg_list = tmp+2; msg_list++;} 39a1db5e05SBarry Smith if (*ptr&BIT_3) 40a1db5e05SBarry Smith {*msg_list = tmp+3; msg_list++;} 41a1db5e05SBarry Smith if (*ptr&BIT_4) 42a1db5e05SBarry Smith {*msg_list = tmp+4; msg_list++;} 43a1db5e05SBarry Smith if (*ptr&BIT_5) 44a1db5e05SBarry Smith {*msg_list = tmp+5; msg_list++;} 45a1db5e05SBarry Smith if (*ptr&BIT_6) 46a1db5e05SBarry Smith {*msg_list = tmp+6; msg_list++;} 47a1db5e05SBarry Smith if (*ptr&BIT_7) 48a1db5e05SBarry Smith {*msg_list = tmp+7; msg_list++;} 49a1db5e05SBarry Smith ptr --; 50a1db5e05SBarry Smith } 51a1db5e05SBarry Smith } 52a1db5e05SBarry Smith PetscFunctionReturn(0); 53a1db5e05SBarry Smith } 54a1db5e05SBarry Smith 55a1db5e05SBarry Smith /*********************************bit_mask.c***********************************/ 56a1db5e05SBarry Smith PetscInt ct_bits( char *ptr, PetscInt n) 57a1db5e05SBarry Smith { 58a1db5e05SBarry Smith PetscInt i, tmp=0; 59a1db5e05SBarry Smith 60a1db5e05SBarry Smith PetscFunctionBegin; 61a1db5e05SBarry Smith for(i=0;i<n;i++) 62a1db5e05SBarry Smith { 63a1db5e05SBarry Smith if (*ptr&128) {tmp++;} 64a1db5e05SBarry Smith if (*ptr&64) {tmp++;} 65a1db5e05SBarry Smith if (*ptr&32) {tmp++;} 66a1db5e05SBarry Smith if (*ptr&16) {tmp++;} 67a1db5e05SBarry Smith if (*ptr&8) {tmp++;} 68a1db5e05SBarry Smith if (*ptr&4) {tmp++;} 69a1db5e05SBarry Smith if (*ptr&2) {tmp++;} 70a1db5e05SBarry Smith if (*ptr&1) {tmp++;} 71a1db5e05SBarry Smith ptr++; 72a1db5e05SBarry Smith } 73a1db5e05SBarry Smith 74a1db5e05SBarry Smith return(tmp); 75a1db5e05SBarry Smith } 76a1db5e05SBarry Smith 77a1db5e05SBarry Smith /*********************************bit_mask.c***********************************/ 78a1db5e05SBarry Smith PetscInt 79a1db5e05SBarry Smith div_ceil( PetscInt numer, PetscInt denom) 80a1db5e05SBarry Smith { 81a1db5e05SBarry Smith PetscInt rt_val; 82a1db5e05SBarry Smith 83a1db5e05SBarry Smith if ((numer<0)||(denom<=0)) 84e32f2f54SBarry Smith {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"div_ceil() :: numer=%D ! >=0, denom=%D ! >0",numer,denom);} 85a1db5e05SBarry Smith 86a1db5e05SBarry Smith /* if integer division remainder then increment */ 87a1db5e05SBarry Smith rt_val = numer/denom; 88a1db5e05SBarry Smith if (numer%denom) 89a1db5e05SBarry Smith {rt_val++;} 90a1db5e05SBarry Smith 91a1db5e05SBarry Smith return(rt_val); 92a1db5e05SBarry Smith } 93a1db5e05SBarry Smith 94a1db5e05SBarry Smith /*********************************bit_mask.c***********************************/ 95a1db5e05SBarry Smith PetscInt 96a1db5e05SBarry Smith len_bit_mask( PetscInt num_items) 97a1db5e05SBarry Smith { 98a1db5e05SBarry Smith PetscInt rt_val, tmp; 99a1db5e05SBarry Smith 100a1db5e05SBarry Smith if (num_items<0) 101e32f2f54SBarry Smith {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Value Sent To len_bit_mask() Must be >= 0!");} 102a1db5e05SBarry Smith 103a1db5e05SBarry Smith /* mod BYTE ceiling function */ 104a1db5e05SBarry Smith rt_val = num_items/BYTE; 105a1db5e05SBarry Smith if (num_items%BYTE) 106a1db5e05SBarry Smith {rt_val++;} 107a1db5e05SBarry Smith 108a1db5e05SBarry Smith /* make mults of sizeof int */ 109a1db5e05SBarry Smith if ((tmp=rt_val%sizeof(PetscInt))) 110a1db5e05SBarry Smith {rt_val+=(sizeof(PetscInt)-tmp);} 111a1db5e05SBarry Smith 112a1db5e05SBarry Smith return(rt_val); 113a1db5e05SBarry Smith } 114a1db5e05SBarry Smith 115a1db5e05SBarry Smith /*********************************bit_mask.c***********************************/ 116a1db5e05SBarry Smith PetscErrorCode set_bit_mask( PetscInt *bm, PetscInt len, PetscInt val) 117a1db5e05SBarry Smith { 118a1db5e05SBarry Smith PetscInt i, offset; 119a1db5e05SBarry Smith char mask = 1; 120a1db5e05SBarry Smith char *cptr; 121a1db5e05SBarry Smith 122*e7e72b3dSBarry Smith if (len_bit_mask(val)>len) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"The Bit Mask Isn't That Large!"); 123a1db5e05SBarry Smith 124a1db5e05SBarry Smith cptr = (char *) bm; 125a1db5e05SBarry Smith 126a1db5e05SBarry Smith offset = len/sizeof(PetscInt); 127a1db5e05SBarry Smith for (i=0;i<offset;i++) 128a1db5e05SBarry Smith {*bm=0; bm++;} 129a1db5e05SBarry Smith 130a1db5e05SBarry Smith offset = val%BYTE; 131a1db5e05SBarry Smith for (i=0;i<offset;i++) 132a1db5e05SBarry Smith {mask <<= 1;} 133a1db5e05SBarry Smith 134a1db5e05SBarry Smith offset = len - val/BYTE - 1; 135a1db5e05SBarry Smith cptr[offset] = mask; 136a1db5e05SBarry Smith PetscFunctionReturn(0); 137a1db5e05SBarry Smith } 138a1db5e05SBarry Smith 139a1db5e05SBarry Smith /*********************************bit_mask.c***********************************/ 140a1db5e05SBarry Smith PetscInt len_buf(PetscInt item_size, PetscInt num_items) 141a1db5e05SBarry Smith { 142a1db5e05SBarry Smith PetscInt rt_val, tmp; 143a1db5e05SBarry Smith 144a1db5e05SBarry Smith PetscFunctionBegin; 145a1db5e05SBarry Smith rt_val = item_size * num_items; 146a1db5e05SBarry Smith 147a1db5e05SBarry Smith /* double precision align for now ... consider page later */ 148a1db5e05SBarry Smith if ((tmp = (rt_val%(PetscInt)sizeof(double)))) 149a1db5e05SBarry Smith {rt_val += (sizeof(double) - tmp);} 150a1db5e05SBarry Smith 151a1db5e05SBarry Smith return(rt_val); 152a1db5e05SBarry Smith } 153a1db5e05SBarry Smith 154a1db5e05SBarry Smith 155a1db5e05SBarry Smith 156