1 /*$Id: sbaijov.c,v 1.65 2001/08/06 21:15:42 bsmith Exp $*/ 2 3 /* 4 Routines to compute overlapping regions of a parallel MPI matrix. 5 Used for finding submatrices that were shared across processors. 6 */ 7 #include "src/mat/impls/sbaij/mpi/mpisbaij.h" 8 #include "petscbt.h" 9 10 static int MatIncreaseOverlap_MPISBAIJ_Once(Mat,int,IS *); 11 static int MatIncreaseOverlap_MPISBAIJ_Local(Mat,int *,int,int **,PetscBT*); 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ" 15 int MatIncreaseOverlap_MPISBAIJ(Mat C,int is_max,IS is[],int ov) 16 { 17 Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 18 int i,ierr,N=C->N, bs=c->bs; 19 IS *is_new; 20 21 PetscFunctionBegin; 22 ierr = PetscMalloc(is_max*sizeof(IS),&is_new);CHKERRQ(ierr); 23 /* Convert the indices into block format */ 24 ierr = ISCompressIndicesGeneral(N,bs,is_max,is,is_new);CHKERRQ(ierr); 25 if (ov < 0){ SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n");} 26 for (i=0; i<ov; ++i) { 27 ierr = MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr); 28 } 29 for (i=0; i<is_max; i++) {ierr = ISDestroy(is[i]);CHKERRQ(ierr);} 30 ierr = ISExpandIndicesGeneral(N,bs,is_max,is_new,is);CHKERRQ(ierr); 31 for (i=0; i<is_max; i++) {ierr = ISDestroy(is_new[i]);CHKERRQ(ierr);} 32 ierr = PetscFree(is_new);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 typedef enum {MINE,OTHER} WhoseOwner; 37 /* data1, odata1 and odata2 are packed in the format (for communication): 38 data[0] = is_max, no of is 39 data[1] = size of is[0] 40 ... 41 data[is_max] = size of is[is_max-1] 42 data[is_max + 1] = data(is[0]) 43 ... 44 data[is_max+1+sum(size of is[k]), k=0,...,i-1] = data(is[i]) 45 ... 46 data2 is packed in the format (for creating output is[]): 47 data[0] = is_max, no of is 48 data[1] = size of is[0] 49 ... 50 data[is_max] = size of is[is_max-1] 51 data[is_max + 1] = data(is[0]) 52 ... 53 data[is_max + 1 + Mbs*i) = data(is[i]) 54 ... 55 */ 56 #undef __FUNCT__ 57 #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Once" 58 static int MatIncreaseOverlap_MPISBAIJ_Once(Mat C,int is_max,IS is[]) 59 { 60 Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 61 int len,idx,*idx_i,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i, 62 size,rank,Mbs,i,j,k,ierr,nrqs,nrqr,*odata1,*odata2, 63 tag1,tag2,flag,proc_id,**odata2_ptr,*ctable=0,*btable,len_b; 64 int *id_r1,*len_r1,proc_end=0,*iwork,*len_s; 65 char *t_p; 66 MPI_Comm comm; 67 MPI_Request *s_waits1,*s_waits2,r_req; 68 MPI_Status *s_status,r_status; 69 PetscBT *table=0; /* mark indices of this processor's is[] */ 70 PetscBT table_i; 71 PetscBT otable; /* mark indices of other processors' is[] */ 72 int bs=c->bs,Bn = c->B->n,Bnbs = Bn/bs,*Bowners; 73 IS is_local,is_gl; 74 75 PetscFunctionBegin; 76 77 comm = C->comm; 78 size = c->size; 79 rank = c->rank; 80 Mbs = c->Mbs; 81 82 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 83 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 84 85 /* 1. Send this processor's is[] to other processors */ 86 /*---------------------------------------------------*/ 87 /* allocate spaces */ 88 ierr = PetscMalloc(is_max*sizeof(int),&n);CHKERRQ(ierr); 89 len = 0; 90 for (i=0; i<is_max; i++) { 91 ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 92 len += n[i]; 93 } 94 if (len == 0) { 95 is_max = 0; 96 } else { 97 len += 1 + is_max; /* max length of data1 for one processor */ 98 } 99 100 ierr = PetscMalloc((size*len+1)*sizeof(int),&data1);CHKERRQ(ierr); 101 ierr = PetscMalloc(size*sizeof(int*),&data1_start);CHKERRQ(ierr); 102 for (i=0; i<size; i++) data1_start[i] = data1 + i*len; 103 104 ierr = PetscMalloc(size*3*sizeof(int),&len_s);CHKERRQ(ierr); 105 btable = len_s + size; 106 iwork = btable + size; 107 108 if (is_max){ 109 /* create hash table ctable which maps c->row to proc_id) */ 110 ierr = PetscMalloc(Mbs*sizeof(int),&ctable);CHKERRQ(ierr); 111 for (proc_id=0,j=0; proc_id<size; proc_id++) { 112 for (; j<c->rowners[proc_id+1]; j++) { 113 ctable[j] = proc_id; 114 } 115 } 116 117 /* create tables for marking indices */ 118 len_b = is_max*sizeof(PetscBT) + (Mbs/PETSC_BITS_PER_BYTE+1)*is_max*sizeof(char) + 1; 119 ierr = PetscMalloc(len_b,&table);CHKERRQ(ierr); 120 t_p = (char *)(table + is_max); 121 for (i=0; i<is_max; i++) { 122 table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 123 } 124 125 ierr = ISCreateGeneral(comm,Bnbs,c->garray,&is_local);CHKERRQ(ierr); 126 ierr = ISAllGather(is_local, &is_gl);CHKERRQ(ierr); 127 ierr = ISDestroy(is_local);CHKERRQ(ierr); 128 129 ierr = PetscMalloc((size+1)*sizeof(int),&Bowners);CHKERRQ(ierr); 130 ierr = MPI_Allgather(&Bnbs,1,MPI_INT,Bowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 131 Bowners[0] = 0; 132 for (i=0; i<size; i++) Bowners[i+1] += Bowners[i]; 133 134 /* hash table table_i[idx] = 1 if idx is on is[] array 135 = 0 otherwise */ 136 table_i = table[0]; 137 ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 138 for (i=0; i<is_max; i++){ 139 ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr); 140 for (j=0; j<n[i]; j++){ 141 idx = idx_i[j]; 142 ierr = PetscBTSet(table_i,idx);CHKERRQ(ierr); 143 } 144 ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr); 145 } 146 147 /* hash table btable[id_proc] = 1 if a col index of B-matrix in [id_proc] is on is[] array 148 = 0 otherwise */ 149 ierr = ISGetIndices(is_gl,&idx_i); 150 for (i=0; i<size; i++){ 151 btable[i] = 0; 152 for (j = Bowners[i]; j<Bowners[i+1]; j++){ /* go through B cols */ 153 idx = idx_i[j]; 154 if(PetscBTLookup(table_i,idx)){ 155 btable[i] = 1; 156 break; 157 } 158 } 159 } 160 ierr = ISRestoreIndices(is_gl,&idx_i);CHKERRQ(ierr); 161 ierr = PetscFree(Bowners);CHKERRQ(ierr); 162 ierr = ISDestroy(is_gl);CHKERRQ(ierr); 163 } /* if (is_max) */ 164 165 /* evaluate communication - mesg to who, length, and buffer space */ 166 for (i=0; i<size; i++) len_s[i] = 0; 167 168 /* header of data1 */ 169 for (proc_id=0; proc_id<size; proc_id++){ 170 iwork[proc_id] = 0; 171 *data1_start[proc_id] = is_max; 172 data1_start[proc_id]++; 173 for (j=0; j<is_max; j++) { 174 if (proc_id == rank){ 175 *data1_start[proc_id] = n[j]; 176 } else { 177 *data1_start[proc_id] = 0; 178 } 179 data1_start[proc_id]++; 180 } 181 } 182 183 for (i=0; i<is_max; i++) { 184 ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr); 185 ierr = PetscSortInt(n[i],idx_i);CHKERRQ(ierr); 186 for (j=0; j<n[i]; j++){ 187 idx = idx_i[j]; 188 *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */ 189 proc_end = ctable[idx]; 190 for (proc_id=0; proc_id<=proc_end; proc_id++){ /* for others to process */ 191 if (proc_id == rank ) continue; /* done before this loop */ 192 if (proc_id < proc_end && !btable[proc_id]) continue; /* no need for sending idx to [proc_id] */ 193 *data1_start[proc_id] = idx; data1_start[proc_id]++; 194 len_s[proc_id]++; 195 } 196 } 197 /* update header data */ 198 for (proc_id=0; proc_id<=proc_end; proc_id++){ 199 if (proc_id== rank) continue; 200 *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id]; 201 iwork[proc_id] = len_s[proc_id] ; 202 } 203 ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr); 204 } 205 206 nrqs = 0; nrqr = 0; 207 for (i=0; i<size; i++){ 208 data1_start[i] = data1 + i*len; 209 if (len_s[i]){ 210 nrqs++; 211 len_s[i] += 1 + is_max; /* add no. of header msg */ 212 } 213 } 214 215 for (i=0; i<is_max; i++) { 216 ierr = ISDestroy(is[i]);CHKERRQ(ierr); 217 } 218 ierr = PetscFree(n);CHKERRQ(ierr); 219 if (ctable){ierr = PetscFree(ctable);CHKERRQ(ierr);} 220 221 /* Determine the number of messages to expect, their lengths, from from-ids */ 222 ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&nrqr);CHKERRQ(ierr); 223 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr); 224 /* ierr = PetscPrintf(PETSC_COMM_SELF, "[%d] nrqs: %d, nrqr: %d\n",rank,nrqs,nrqr); */ 225 226 /* Now post the sends */ 227 ierr = PetscMalloc(2*size*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 228 s_waits2 = s_waits1 + size; 229 k = 0; 230 for (proc_id=0; proc_id<size; proc_id++){ /* send data1 to processor [proc_id] */ 231 if (len_s[proc_id]){ 232 ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPI_INT,proc_id,tag1,comm,s_waits1+k);CHKERRQ(ierr); 233 k++; 234 } 235 } 236 237 /* 2. Do local work on this processor's is[] */ 238 /*-------------------------------------------*/ 239 ierr = PetscBTCreate(Mbs,otable);CHKERRQ(ierr); 240 ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,&data,table);CHKERRQ(ierr); 241 ierr = PetscFree(data1_start);CHKERRQ(ierr); 242 243 /* 3. Receive other's is[] and process. Then send back */ 244 /*-----------------------------------------------------*/ 245 len = 0; 246 for (i=0; i<nrqr; i++){ 247 if (len_r1[i] > len)len = len_r1[i]; 248 /* ierr = PetscPrintf(PETSC_COMM_SELF, "[%d] expect to recv len=%d from [%d]\n",rank,len_r1[i],id_r1[i]); */ 249 } 250 ierr = PetscMalloc((len+1)*sizeof(int),&odata1);CHKERRQ(ierr); 251 ierr = PetscMalloc(size*sizeof(int**),&odata2_ptr);CHKERRQ(ierr); 252 k = 0; 253 while (k < nrqr){ 254 /* Receive messages */ 255 ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRQ(ierr); 256 if (flag){ 257 ierr = MPI_Get_count(&r_status,MPI_INT,&len);CHKERRQ(ierr); 258 proc_id = r_status.MPI_SOURCE; 259 ierr = MPI_Irecv(odata1,len,MPI_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); 260 ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); 261 /* ierr = PetscPrintf(PETSC_COMM_SELF, " [%d] recv %d from [%d]\n",rank,len,proc_id); */ 262 263 /* Process messages */ 264 ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,&odata2_ptr[k],&otable);CHKERRQ(ierr); 265 odata2 = odata2_ptr[k]; 266 len = 1 + odata2[0]; 267 for (i=0; i<odata2[0]; i++){ 268 len += odata2[1 + i]; 269 } 270 271 /* Send messages back */ 272 ierr = MPI_Isend(odata2,len,MPI_INT,proc_id,tag2,comm,s_waits2+k);CHKERRQ(ierr); 273 /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] send %d back to [%d] \n",rank,len,proc_id); */ 274 k++; 275 } 276 } 277 ierr = PetscFree(odata1);CHKERRQ(ierr); 278 279 /* 4. Receive work done on other processors, then merge */ 280 /*--------------------------------------------------------*/ 281 /* Allocate memory for incoming data */ 282 len = (1+is_max*(Mbs+1)); 283 ierr = PetscMalloc(len*sizeof(int),&data2);CHKERRQ(ierr); 284 285 k = 0; 286 while (k < nrqs){ 287 /* Receive messages */ 288 ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status); 289 if (flag){ 290 ierr = MPI_Get_count(&r_status,MPI_INT,&len);CHKERRQ(ierr); 291 proc_id = r_status.MPI_SOURCE; 292 ierr = MPI_Irecv(data2,len,MPI_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); 293 ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); 294 /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] recv %d from [%d], data2:\n",rank,len,proc_id); */ 295 if (len > 1+is_max){ /* Add data2 into data */ 296 data2_i = data2 + 1 + is_max; 297 for (i=0; i<is_max; i++){ 298 table_i = table[i]; 299 data_i = data + 1 + is_max + Mbs*i; 300 isz = data[1+i]; 301 for (j=0; j<data2[1+i]; j++){ 302 col = data2_i[j]; 303 if (!PetscBTLookupSet(table_i,col)) {data_i[isz++] = col;} 304 } 305 data[1+i] = isz; 306 if (i < is_max - 1) data2_i += data2[1+i]; 307 } 308 } 309 k++; 310 } 311 } 312 313 /* phase 1 sends are complete */ 314 ierr = PetscMalloc(size*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 315 if (nrqs){ 316 ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr); 317 } 318 ierr = PetscFree(data1);CHKERRQ(ierr); 319 320 /* phase 3 sends are complete */ 321 if (nrqr){ 322 ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRQ(ierr); 323 } 324 for (k=0; k<nrqr; k++){ 325 ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr); 326 } 327 ierr = PetscFree(odata2_ptr);CHKERRQ(ierr); 328 329 /* 5. Create new is[] */ 330 /*--------------------*/ 331 for (i=0; i<is_max; i++) { 332 data_i = data + 1 + is_max + Mbs*i; 333 ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,is+i);CHKERRQ(ierr); 334 } 335 336 ierr = PetscFree(data2);CHKERRQ(ierr); 337 ierr = PetscFree(data);CHKERRQ(ierr); 338 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 339 ierr = PetscFree(s_status);CHKERRQ(ierr); 340 if (table) {ierr = PetscFree(table);CHKERRQ(ierr);} 341 ierr = PetscBTDestroy(otable);CHKERRQ(ierr); 342 343 ierr = PetscFree(len_s);CHKERRQ(ierr); 344 ierr = PetscFree(id_r1);CHKERRQ(ierr); 345 ierr = PetscFree(len_r1);CHKERRQ(ierr); 346 347 PetscFunctionReturn(0); 348 } 349 350 #undef __FUNCT__ 351 #define __FUNCT__ "MatIncreaseOverlap_MPISBAIJ_Local" 352 /* 353 MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do 354 the work on the local processor. 355 356 Inputs: 357 C - MAT_MPISBAIJ; 358 data - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format. 359 whose - whose is[] to be processed, 360 MINE: this processor's is[] 361 OTHER: other processor's is[] 362 Output: 363 data_new - whose = MINE: 364 holds input and newly found indices in the same format as data 365 whose = OTHER: 366 only holds the newly found indices 367 table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'. 368 */ 369 static int MatIncreaseOverlap_MPISBAIJ_Local(Mat C,int *data,int whose,int **data_new,PetscBT *table) 370 { 371 Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 372 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data; 373 Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(c->B)->data; 374 int ierr,row,mbs,Mbs,*nidx,*nidx_i,col,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l; 375 int a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n; 376 PetscBT table0; /* mark the indices of input is[] for look up */ 377 PetscBT table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */ 378 379 PetscFunctionBegin; 380 Mbs = c->Mbs; mbs = a->mbs; 381 ai = a->i; aj = a->j; 382 bi = b->i; bj = b->j; 383 garray = c->garray; 384 rstart = c->rstart; 385 is_max = data[0]; 386 387 ierr = PetscBTCreate(Mbs,table0);CHKERRQ(ierr); 388 389 ierr = PetscMalloc((1+is_max*(Mbs+1))*sizeof(int),&nidx);CHKERRQ(ierr); 390 nidx[0] = is_max; 391 392 idx_i = data + is_max + 1; /* ptr to input is[0] array */ 393 nidx_i = nidx + is_max + 1; /* ptr to output is[0] array */ 394 for (i=0; i<is_max; i++) { /* for each is */ 395 isz = 0; 396 n = data[1+i]; /* size of input is[i] */ 397 398 /* initialize table_i, set table0 */ 399 if (whose == MINE){ /* process this processor's is[] */ 400 table_i = table[i]; 401 nidx_i = nidx + 1+ is_max + Mbs*i; 402 } else { /* process other processor's is[] - only use one temp table */ 403 table_i = *table; 404 } 405 ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 406 ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr); 407 if (n > 0) { 408 isz0 = 0; 409 for (j=0; j<n; j++){ 410 col = idx_i[j]; 411 if (col >= Mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"index col %d >= Mbs %d",col,Mbs); 412 if(!PetscBTLookupSet(table_i,col)) { 413 ierr = PetscBTSet(table0,col);CHKERRQ(ierr); 414 if (whose == MINE) {nidx_i[isz0] = col;} 415 isz0++; 416 } 417 } 418 419 if (whose == MINE) {isz = isz0;} 420 k = 0; /* no. of indices from input is[i] that have been examined */ 421 for (row=0; row<mbs; row++){ 422 a_start = ai[row]; a_end = ai[row+1]; 423 b_start = bi[row]; b_end = bi[row+1]; 424 if (PetscBTLookup(table0,row+rstart)){ /* row is on input is[i]: 425 do row search: collect all col in this row */ 426 for (l = a_start; l<a_end ; l++){ /* Amat */ 427 col = aj[l] + rstart; 428 if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;} 429 } 430 for (l = b_start; l<b_end ; l++){ /* Bmat */ 431 col = garray[bj[l]]; 432 if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;} 433 } 434 k++; 435 if (k >= isz0) break; /* for (row=0; row<mbs; row++) */ 436 } else { /* row is not on input is[i]: 437 do col serach: add row onto nidx_i if there is a col in nidx_i */ 438 for (l = a_start; l<a_end ; l++){ /* Amat */ 439 col = aj[l] + rstart; 440 if (PetscBTLookup(table0,col)){ 441 if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;} 442 break; /* for l = start; l<end ; l++) */ 443 } 444 } 445 for (l = b_start; l<b_end ; l++){ /* Bmat */ 446 col = garray[bj[l]]; 447 if (PetscBTLookup(table0,col)){ 448 if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;} 449 break; /* for l = start; l<end ; l++) */ 450 } 451 } 452 } 453 } 454 } /* if (n > 0) */ 455 456 if (i < is_max - 1){ 457 idx_i += n; /* ptr to input is[i+1] array */ 458 nidx_i += isz; /* ptr to output is[i+1] array */ 459 } 460 nidx[1+i] = isz; /* size of new is[i] */ 461 } /* for each is */ 462 *data_new = nidx; 463 ierr = PetscBTDestroy(table0);CHKERRQ(ierr); 464 465 PetscFunctionReturn(0); 466 } 467 468 469