1 2 /* 3 Routines to compute overlapping regions of a parallel MPI matrix. 4 Used for finding submatrices that were shared across processors. 5 */ 6 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 7 #include <petscbt.h> 8 9 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat,PetscInt,IS*); 10 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat,PetscInt*,PetscInt,PetscInt*,PetscBT*); 11 12 PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat C,PetscInt is_max,IS is[],PetscInt ov) 13 { 14 PetscErrorCode ierr; 15 PetscInt i,N=C->cmap->N, bs=C->rmap->bs,M=C->rmap->N,Mbs=M/bs,*nidx,isz,iov; 16 IS *is_new,*is_row; 17 Mat *submats; 18 Mat_MPISBAIJ *c=(Mat_MPISBAIJ*)C->data; 19 Mat_SeqSBAIJ *asub_i; 20 PetscBT table; 21 PetscInt *ai,brow,nz,nis,l,nmax,nstages_local,nstages,max_no,pos; 22 const PetscInt *idx; 23 PetscBool flg,*allcolumns,*allrows; 24 25 PetscFunctionBegin; 26 ierr = PetscMalloc1(is_max,&is_new);CHKERRQ(ierr); 27 /* Convert the indices into block format */ 28 ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,is_max,is,is_new);CHKERRQ(ierr); 29 if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n"); 30 31 /* ----- previous non-scalable implementation ----- */ 32 flg = PETSC_FALSE; 33 ierr = PetscOptionsHasName(NULL,NULL, "-IncreaseOverlap_old", &flg);CHKERRQ(ierr); 34 if (flg) { /* previous non-scalable implementation */ 35 printf("use previous non-scalable implementation...\n"); 36 for (i=0; i<ov; ++i) { 37 ierr = MatIncreaseOverlap_MPISBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr); 38 } 39 } else { /* implementation using modified BAIJ routines */ 40 41 ierr = PetscMalloc1(Mbs+1,&nidx);CHKERRQ(ierr); 42 ierr = PetscBTCreate(Mbs,&table);CHKERRQ(ierr); /* for column search */ 43 ierr = PetscMalloc2(is_max+1,&allcolumns,is_max+1,&allrows);CHKERRQ(ierr); 44 45 /* Create is_row */ 46 ierr = PetscMalloc1(is_max,&is_row);CHKERRQ(ierr); 47 ierr = ISCreateStride(PETSC_COMM_SELF,Mbs,0,1,&is_row[0]);CHKERRQ(ierr); 48 49 allrows[0] = PETSC_TRUE; 50 for (i=1; i<is_max; i++) { 51 is_row[i] = is_row[0]; /* reuse is_row[0] */ 52 allrows[i] = PETSC_TRUE; 53 } 54 55 /* Allocate memory to hold all the submatrices - Modified from MatGetSubMatrices_MPIBAIJ() */ 56 ierr = PetscMalloc1(is_max+1,&submats);CHKERRQ(ierr); 57 58 /* Check for special case: each processor gets entire matrix columns */ 59 for (i=0; i<is_max; i++) { 60 ierr = ISIdentity(is_new[i],&flg);CHKERRQ(ierr); 61 ierr = ISGetLocalSize(is_new[i],&isz);CHKERRQ(ierr); 62 if (flg && isz == Mbs) { 63 allcolumns[i] = PETSC_TRUE; 64 } else { 65 allcolumns[i] = PETSC_FALSE; 66 } 67 } 68 69 /* Determine the number of stages through which submatrices are done */ 70 nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 71 if (!nmax) nmax = 1; 72 nstages_local = is_max/nmax + ((is_max % nmax) ? 1 : 0); 73 74 /* Make sure every processor loops through the nstages */ 75 ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 76 77 for (iov=0; iov<ov; ++iov) { 78 /* 1) Get submats for column search */ 79 for (i=0,pos=0; i<nstages; i++) { 80 if (pos+nmax <= is_max) max_no = nmax; 81 else if (pos == is_max) max_no = 0; 82 else max_no = is_max-pos; 83 c->ijonly = PETSC_TRUE; 84 ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,is_row+pos,is_new+pos,MAT_INITIAL_MATRIX,allrows+pos,allcolumns+pos,submats+pos);CHKERRQ(ierr); 85 pos += max_no; 86 } 87 88 /* 2) Row search */ 89 ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,is_max,is_new);CHKERRQ(ierr); 90 91 /* 3) Column search */ 92 for (i=0; i<is_max; i++) { 93 asub_i = (Mat_SeqSBAIJ*)submats[i]->data; 94 ai = asub_i->i;; 95 96 /* put is_new obtained from MatIncreaseOverlap_MPIBAIJ() to table */ 97 ierr = PetscBTMemzero(Mbs,table);CHKERRQ(ierr); 98 99 ierr = ISGetIndices(is_new[i],&idx);CHKERRQ(ierr); 100 ierr = ISGetLocalSize(is_new[i],&nis);CHKERRQ(ierr); 101 for (l=0; l<nis; l++) { 102 ierr = PetscBTSet(table,idx[l]);CHKERRQ(ierr); 103 nidx[l] = idx[l]; 104 } 105 isz = nis; 106 107 /* add column entries to table */ 108 for (brow=0; brow<Mbs; brow++) { 109 nz = ai[brow+1] - ai[brow]; 110 if (nz) { 111 if (!PetscBTLookupSet(table,brow)) nidx[isz++] = brow; 112 } 113 } 114 ierr = ISRestoreIndices(is_new[i],&idx);CHKERRQ(ierr); 115 ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr); 116 117 /* create updated is_new */ 118 ierr = ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,is_new+i);CHKERRQ(ierr); 119 } 120 121 /* Free tmp spaces */ 122 for (i=0; i<is_max; i++) { 123 ierr = MatDestroy(&submats[i]);CHKERRQ(ierr); 124 } 125 } 126 ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 127 ierr = PetscBTDestroy(&table);CHKERRQ(ierr); 128 ierr = PetscFree(submats);CHKERRQ(ierr); 129 ierr = ISDestroy(&is_row[0]);CHKERRQ(ierr); 130 ierr = PetscFree(is_row);CHKERRQ(ierr); 131 ierr = PetscFree(nidx);CHKERRQ(ierr); 132 133 } 134 135 for (i=0; i<is_max; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 136 ierr = ISExpandIndicesGeneral(N,N,bs,is_max,is_new,is);CHKERRQ(ierr); 137 138 for (i=0; i<is_max; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 139 ierr = PetscFree(is_new);CHKERRQ(ierr); 140 PetscFunctionReturn(0); 141 } 142 143 typedef enum {MINE,OTHER} WhoseOwner; 144 /* data1, odata1 and odata2 are packed in the format (for communication): 145 data[0] = is_max, no of is 146 data[1] = size of is[0] 147 ... 148 data[is_max] = size of is[is_max-1] 149 data[is_max + 1] = data(is[0]) 150 ... 151 data[is_max+1+sum(size of is[k]), k=0,...,i-1] = data(is[i]) 152 ... 153 data2 is packed in the format (for creating output is[]): 154 data[0] = is_max, no of is 155 data[1] = size of is[0] 156 ... 157 data[is_max] = size of is[is_max-1] 158 data[is_max + 1] = data(is[0]) 159 ... 160 data[is_max + 1 + Mbs*i) = data(is[i]) 161 ... 162 */ 163 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,PetscInt is_max,IS is[]) 164 { 165 Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 166 PetscErrorCode ierr; 167 PetscMPIInt size,rank,tag1,tag2,*len_s,nrqr,nrqs,*id_r1,*len_r1,flag,len,*iwork; 168 const PetscInt *idx_i; 169 PetscInt idx,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i; 170 PetscInt Mbs,i,j,k,*odata1,*odata2; 171 PetscInt proc_id,**odata2_ptr,*ctable=0,*btable,len_max,len_est; 172 PetscInt proc_end=0,len_unused,nodata2; 173 PetscInt ois_max; /* max no of is[] in each of processor */ 174 char *t_p; 175 MPI_Comm comm; 176 MPI_Request *s_waits1,*s_waits2,r_req; 177 MPI_Status *s_status,r_status; 178 PetscBT *table; /* mark indices of this processor's is[] */ 179 PetscBT table_i; 180 PetscBT otable; /* mark indices of other processors' is[] */ 181 PetscInt bs=C->rmap->bs,Bn = c->B->cmap->n,Bnbs = Bn/bs,*Bowners; 182 IS garray_local,garray_gl; 183 184 PetscFunctionBegin; 185 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 186 size = c->size; 187 rank = c->rank; 188 Mbs = c->Mbs; 189 190 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 191 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 192 193 /* create tables used in 194 step 1: table[i] - mark c->garray of proc [i] 195 step 3: table[i] - mark indices of is[i] when whose=MINE 196 table[0] - mark incideces of is[] when whose=OTHER */ 197 len = PetscMax(is_max, size);CHKERRQ(ierr); 198 ierr = PetscMalloc2(len,&table,(Mbs/PETSC_BITS_PER_BYTE+1)*len,&t_p);CHKERRQ(ierr); 199 for (i=0; i<len; i++) { 200 table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 201 } 202 203 ierr = MPIU_Allreduce(&is_max,&ois_max,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 204 205 /* 1. Send this processor's is[] to other processors */ 206 /*---------------------------------------------------*/ 207 /* allocate spaces */ 208 ierr = PetscMalloc1(is_max,&n);CHKERRQ(ierr); 209 len = 0; 210 for (i=0; i<is_max; i++) { 211 ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 212 len += n[i]; 213 } 214 if (!len) { 215 is_max = 0; 216 } else { 217 len += 1 + is_max; /* max length of data1 for one processor */ 218 } 219 220 221 ierr = PetscMalloc1(size*len+1,&data1);CHKERRQ(ierr); 222 ierr = PetscMalloc1(size,&data1_start);CHKERRQ(ierr); 223 for (i=0; i<size; i++) data1_start[i] = data1 + i*len; 224 225 ierr = PetscMalloc4(size,&len_s,size,&btable,size,&iwork,size+1,&Bowners);CHKERRQ(ierr); 226 227 /* gather c->garray from all processors */ 228 ierr = ISCreateGeneral(comm,Bnbs,c->garray,PETSC_COPY_VALUES,&garray_local);CHKERRQ(ierr); 229 ierr = ISAllGather(garray_local, &garray_gl);CHKERRQ(ierr); 230 ierr = ISDestroy(&garray_local);CHKERRQ(ierr); 231 ierr = MPI_Allgather(&Bnbs,1,MPIU_INT,Bowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 232 233 Bowners[0] = 0; 234 for (i=0; i<size; i++) Bowners[i+1] += Bowners[i]; 235 236 if (is_max) { 237 /* hash table ctable which maps c->row to proc_id) */ 238 ierr = PetscMalloc1(Mbs,&ctable);CHKERRQ(ierr); 239 for (proc_id=0,j=0; proc_id<size; proc_id++) { 240 for (; j<C->rmap->range[proc_id+1]/bs; j++) ctable[j] = proc_id; 241 } 242 243 /* hash tables marking c->garray */ 244 ierr = ISGetIndices(garray_gl,&idx_i);CHKERRQ(ierr); 245 for (i=0; i<size; i++) { 246 table_i = table[i]; 247 ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 248 for (j = Bowners[i]; j<Bowners[i+1]; j++) { /* go through B cols of proc[i]*/ 249 ierr = PetscBTSet(table_i,idx_i[j]);CHKERRQ(ierr); 250 } 251 } 252 ierr = ISRestoreIndices(garray_gl,&idx_i);CHKERRQ(ierr); 253 } /* if (is_max) */ 254 ierr = ISDestroy(&garray_gl);CHKERRQ(ierr); 255 256 /* evaluate communication - mesg to who, length, and buffer space */ 257 for (i=0; i<size; i++) len_s[i] = 0; 258 259 /* header of data1 */ 260 for (proc_id=0; proc_id<size; proc_id++) { 261 iwork[proc_id] = 0; 262 *data1_start[proc_id] = is_max; 263 data1_start[proc_id]++; 264 for (j=0; j<is_max; j++) { 265 if (proc_id == rank) { 266 *data1_start[proc_id] = n[j]; 267 } else { 268 *data1_start[proc_id] = 0; 269 } 270 data1_start[proc_id]++; 271 } 272 } 273 274 for (i=0; i<is_max; i++) { 275 ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr); 276 for (j=0; j<n[i]; j++) { 277 idx = idx_i[j]; 278 *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */ 279 proc_end = ctable[idx]; 280 for (proc_id=0; proc_id<=proc_end; proc_id++) { /* for others to process */ 281 if (proc_id == rank) continue; /* done before this loop */ 282 if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx)) continue; /* no need for sending idx to [proc_id] */ 283 *data1_start[proc_id] = idx; data1_start[proc_id]++; 284 len_s[proc_id]++; 285 } 286 } 287 /* update header data */ 288 for (proc_id=0; proc_id<size; proc_id++) { 289 if (proc_id== rank) continue; 290 *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id]; 291 iwork[proc_id] = len_s[proc_id]; 292 } 293 ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr); 294 } 295 296 nrqs = 0; nrqr = 0; 297 for (i=0; i<size; i++) { 298 data1_start[i] = data1 + i*len; 299 if (len_s[i]) { 300 nrqs++; 301 len_s[i] += 1 + is_max; /* add no. of header msg */ 302 } 303 } 304 305 for (i=0; i<is_max; i++) { 306 ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 307 } 308 ierr = PetscFree(n);CHKERRQ(ierr); 309 ierr = PetscFree(ctable);CHKERRQ(ierr); 310 311 /* Determine the number of messages to expect, their lengths, from from-ids */ 312 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrqr);CHKERRQ(ierr); 313 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr); 314 315 /* Now post the sends */ 316 ierr = PetscMalloc2(size,&s_waits1,size,&s_waits2);CHKERRQ(ierr); 317 k = 0; 318 for (proc_id=0; proc_id<size; proc_id++) { /* send data1 to processor [proc_id] */ 319 if (len_s[proc_id]) { 320 ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPIU_INT,proc_id,tag1,comm,s_waits1+k);CHKERRQ(ierr); 321 k++; 322 } 323 } 324 325 /* 2. Receive other's is[] and process. Then send back */ 326 /*-----------------------------------------------------*/ 327 len = 0; 328 for (i=0; i<nrqr; i++) { 329 if (len_r1[i] > len) len = len_r1[i]; 330 } 331 ierr = PetscFree(len_r1);CHKERRQ(ierr); 332 ierr = PetscFree(id_r1);CHKERRQ(ierr); 333 334 for (proc_id=0; proc_id<size; proc_id++) len_s[proc_id] = iwork[proc_id] = 0; 335 336 ierr = PetscMalloc1(len+1,&odata1);CHKERRQ(ierr); 337 ierr = PetscMalloc1(size,&odata2_ptr);CHKERRQ(ierr); 338 ierr = PetscBTCreate(Mbs,&otable);CHKERRQ(ierr); 339 340 len_max = ois_max*(Mbs+1); /* max space storing all is[] for each receive */ 341 len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */ 342 ierr = PetscMalloc1(len_est+1,&odata2);CHKERRQ(ierr); 343 nodata2 = 0; /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */ 344 345 odata2_ptr[nodata2] = odata2; 346 347 len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max */ 348 349 k = 0; 350 while (k < nrqr) { 351 /* Receive messages */ 352 ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRQ(ierr); 353 if (flag) { 354 ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr); 355 proc_id = r_status.MPI_SOURCE; 356 ierr = MPI_Irecv(odata1,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); 357 ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); 358 359 /* Process messages */ 360 /* make sure there is enough unused space in odata2 array */ 361 if (len_unused < len_max) { /* allocate more space for odata2 */ 362 ierr = PetscMalloc1(len_est+1,&odata2);CHKERRQ(ierr); 363 364 odata2_ptr[++nodata2] = odata2; 365 366 len_unused = len_est; 367 } 368 369 ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);CHKERRQ(ierr); 370 len = 1 + odata2[0]; 371 for (i=0; i<odata2[0]; i++) len += odata2[1 + i]; 372 373 /* Send messages back */ 374 ierr = MPI_Isend(odata2,len,MPIU_INT,proc_id,tag2,comm,s_waits2+k);CHKERRQ(ierr); 375 k++; 376 odata2 += len; 377 len_unused -= len; 378 len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */ 379 } 380 } 381 ierr = PetscFree(odata1);CHKERRQ(ierr); 382 ierr = PetscBTDestroy(&otable);CHKERRQ(ierr); 383 384 /* 3. Do local work on this processor's is[] */ 385 /*-------------------------------------------*/ 386 /* make sure there is enough unused space in odata2(=data) array */ 387 len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */ 388 if (len_unused < len_max) { /* allocate more space for odata2 */ 389 ierr = PetscMalloc1(len_est+1,&odata2);CHKERRQ(ierr); 390 391 odata2_ptr[++nodata2] = odata2; 392 } 393 394 data = odata2; 395 ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);CHKERRQ(ierr); 396 ierr = PetscFree(data1_start);CHKERRQ(ierr); 397 398 /* 4. Receive work done on other processors, then merge */ 399 /*------------------------------------------------------*/ 400 /* get max number of messages that this processor expects to recv */ 401 ierr = MPIU_Allreduce(len_s,iwork,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 402 ierr = PetscMalloc1(iwork[rank]+1,&data2);CHKERRQ(ierr); 403 ierr = PetscFree4(len_s,btable,iwork,Bowners);CHKERRQ(ierr); 404 405 k = 0; 406 while (k < nrqs) { 407 /* Receive messages */ 408 ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status);CHKERRQ(ierr); 409 if (flag) { 410 ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr); 411 412 proc_id = r_status.MPI_SOURCE; 413 414 ierr = MPI_Irecv(data2,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); 415 ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); 416 if (len > 1+is_max) { /* Add data2 into data */ 417 data2_i = data2 + 1 + is_max; 418 for (i=0; i<is_max; i++) { 419 table_i = table[i]; 420 data_i = data + 1 + is_max + Mbs*i; 421 isz = data[1+i]; 422 for (j=0; j<data2[1+i]; j++) { 423 col = data2_i[j]; 424 if (!PetscBTLookupSet(table_i,col)) data_i[isz++] = col; 425 } 426 data[1+i] = isz; 427 if (i < is_max - 1) data2_i += data2[1+i]; 428 } 429 } 430 k++; 431 } 432 } 433 ierr = PetscFree(data2);CHKERRQ(ierr); 434 ierr = PetscFree2(table,t_p);CHKERRQ(ierr); 435 436 /* phase 1 sends are complete */ 437 ierr = PetscMalloc1(size,&s_status);CHKERRQ(ierr); 438 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 439 ierr = PetscFree(data1);CHKERRQ(ierr); 440 441 /* phase 2 sends are complete */ 442 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRQ(ierr);} 443 ierr = PetscFree2(s_waits1,s_waits2);CHKERRQ(ierr); 444 ierr = PetscFree(s_status);CHKERRQ(ierr); 445 446 /* 5. Create new is[] */ 447 /*--------------------*/ 448 for (i=0; i<is_max; i++) { 449 data_i = data + 1 + is_max + Mbs*i; 450 ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 451 } 452 for (k=0; k<=nodata2; k++) { 453 ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr); 454 } 455 ierr = PetscFree(odata2_ptr);CHKERRQ(ierr); 456 PetscFunctionReturn(0); 457 } 458 459 /* 460 MatIncreaseOverlap_MPISBAIJ_Local - Called by MatIncreaseOverlap, to do 461 the work on the local processor. 462 463 Inputs: 464 C - MAT_MPISBAIJ; 465 data - holds is[]. See MatIncreaseOverlap_MPISBAIJ_Once() for the format. 466 whose - whose is[] to be processed, 467 MINE: this processor's is[] 468 OTHER: other processor's is[] 469 Output: 470 nidx - whose = MINE: 471 holds input and newly found indices in the same format as data 472 whose = OTHER: 473 only holds the newly found indices 474 table - table[i]: mark the indices of is[i], i=0,...,is_max. Used only in the case 'whose=MINE'. 475 */ 476 /* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */ 477 static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,PetscInt *data,PetscInt whose,PetscInt *nidx,PetscBT *table) 478 { 479 Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; 480 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data; 481 Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(c->B)->data; 482 PetscErrorCode ierr; 483 PetscInt row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l; 484 PetscInt a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n; 485 PetscBT table0; /* mark the indices of input is[] for look up */ 486 PetscBT table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */ 487 488 PetscFunctionBegin; 489 Mbs = c->Mbs; mbs = a->mbs; 490 ai = a->i; aj = a->j; 491 bi = b->i; bj = b->j; 492 garray = c->garray; 493 rstart = c->rstartbs; 494 is_max = data[0]; 495 496 ierr = PetscBTCreate(Mbs,&table0);CHKERRQ(ierr); 497 498 nidx[0] = is_max; 499 idx_i = data + is_max + 1; /* ptr to input is[0] array */ 500 nidx_i = nidx + is_max + 1; /* ptr to output is[0] array */ 501 for (i=0; i<is_max; i++) { /* for each is */ 502 isz = 0; 503 n = data[1+i]; /* size of input is[i] */ 504 505 /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */ 506 if (whose == MINE) { /* process this processor's is[] */ 507 table_i = table[i]; 508 nidx_i = nidx + 1+ is_max + Mbs*i; 509 } else { /* process other processor's is[] - only use one temp table */ 510 table_i = table[0]; 511 } 512 ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); 513 ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr); 514 if (n==0) { 515 nidx[1+i] = 0; /* size of new is[i] */ 516 continue; 517 } 518 519 isz0 = 0; col_max = 0; 520 for (j=0; j<n; j++) { 521 col = idx_i[j]; 522 if (col >= Mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index col %D >= Mbs %D",col,Mbs); 523 if (!PetscBTLookupSet(table_i,col)) { 524 ierr = PetscBTSet(table0,col);CHKERRQ(ierr); 525 if (whose == MINE) nidx_i[isz0] = col; 526 if (col_max < col) col_max = col; 527 isz0++; 528 } 529 } 530 531 if (whose == MINE) isz = isz0; 532 k = 0; /* no. of indices from input is[i] that have been examined */ 533 for (row=0; row<mbs; row++) { 534 a_start = ai[row]; a_end = ai[row+1]; 535 b_start = bi[row]; b_end = bi[row+1]; 536 if (PetscBTLookup(table0,row+rstart)) { /* row is on input is[i]: 537 do row search: collect all col in this row */ 538 for (l = a_start; l<a_end ; l++) { /* Amat */ 539 col = aj[l] + rstart; 540 if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col; 541 } 542 for (l = b_start; l<b_end ; l++) { /* Bmat */ 543 col = garray[bj[l]]; 544 if (!PetscBTLookupSet(table_i,col)) nidx_i[isz++] = col; 545 } 546 k++; 547 if (k >= isz0) break; /* for (row=0; row<mbs; row++) */ 548 } else { /* row is not on input is[i]: 549 do col serach: add row onto nidx_i if there is a col in nidx_i */ 550 for (l = a_start; l<a_end; l++) { /* Amat */ 551 col = aj[l] + rstart; 552 if (col > col_max) break; 553 if (PetscBTLookup(table0,col)) { 554 if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart; 555 break; /* for l = start; l<end ; l++) */ 556 } 557 } 558 for (l = b_start; l<b_end; l++) { /* Bmat */ 559 col = garray[bj[l]]; 560 if (col > col_max) break; 561 if (PetscBTLookup(table0,col)) { 562 if (!PetscBTLookupSet(table_i,row+rstart)) nidx_i[isz++] = row+rstart; 563 break; /* for l = start; l<end ; l++) */ 564 } 565 } 566 } 567 } 568 569 if (i < is_max - 1) { 570 idx_i += n; /* ptr to input is[i+1] array */ 571 nidx_i += isz; /* ptr to output is[i+1] array */ 572 } 573 nidx[1+i] = isz; /* size of new is[i] */ 574 } /* for each is */ 575 ierr = PetscBTDestroy(&table0);CHKERRQ(ierr); 576 PetscFunctionReturn(0); 577 } 578 579 580