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