Actual source code: mpiov.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1: /*
  2:    Routines to compute overlapping regions of a parallel MPI matrix
  3:   and to find submatrices that were shared across processors.
  4: */
  5:  #include <../src/mat/impls/aij/seq/aij.h>
  6:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  7:  #include <petscbt.h>
  8:  #include <petscsf.h>

 10: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Once(Mat,PetscInt,IS*);
 11: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Local(Mat,PetscInt,char**,PetscInt*,PetscInt**,PetscTable*);
 12: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Receive(Mat,PetscInt,PetscInt**,PetscInt**,PetscInt*);
 13: extern PetscErrorCode MatGetRow_MPIAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
 14: extern PetscErrorCode MatRestoreRow_MPIAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);

 16: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Once_Scalable(Mat,PetscInt,IS*);
 17: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Local_Scalable(Mat,PetscInt,IS*);
 18: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Send_Scalable(Mat,PetscInt,PetscMPIInt,PetscMPIInt *,PetscInt *, PetscInt *,PetscInt **,PetscInt **);
 19: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Receive_Scalable(Mat,PetscInt,IS*,PetscInt,PetscInt *);


 22: PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov)
 23: {
 25:   PetscInt       i;

 28:   if (ov < 0) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
 29:   for (i=0; i<ov; ++i) {
 30:     MatIncreaseOverlap_MPIAIJ_Once(C,imax,is);
 31:   }
 32:   return(0);
 33: }

 35: PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat C,PetscInt imax,IS is[],PetscInt ov)
 36: {
 38:   PetscInt       i;

 41:   if (ov < 0) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
 42:   for (i=0; i<ov; ++i) {
 43:     MatIncreaseOverlap_MPIAIJ_Once_Scalable(C,imax,is);
 44:   }
 45:   return(0);
 46: }


 49: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Once_Scalable(Mat mat,PetscInt nidx,IS is[])
 50: {
 52:   MPI_Comm       comm;
 53:   PetscInt       *length,length_i,tlength,*remoterows,nrrows,reducednrrows,*rrow_ranks,*rrow_isids,i,j,owner;
 54:   PetscInt       *tosizes,*tosizes_temp,*toffsets,*fromsizes,*todata,*fromdata;
 55:   PetscInt       nrecvrows,*sbsizes = 0,*sbdata = 0;
 56:   const PetscInt *indices_i,**indices;
 57:   PetscLayout    rmap;
 58:   PetscMPIInt    rank,size,*toranks,*fromranks,nto,nfrom;
 59:   PetscSF        sf;
 60:   PetscSFNode    *remote;

 63:   PetscObjectGetComm((PetscObject)mat,&comm);
 64:   MPI_Comm_rank(comm,&rank);
 65:   MPI_Comm_size(comm,&size);
 66:   /* get row map to determine where rows should be going */
 67:   MatGetLayouts(mat,&rmap,NULL);
 68:   /* retrieve IS data and put all together so that we
 69:    * can optimize communication
 70:    *  */
 71:   PetscMalloc2(nidx,(PetscInt ***)&indices,nidx,&length);
 72:   for (i=0,tlength=0; i<nidx; i++){
 73:     ISGetLocalSize(is[i],&length[i]);
 74:     tlength += length[i];
 75:     ISGetIndices(is[i],&indices[i]);
 76:   }
 77:   /* find these rows on remote processors */
 78:   PetscCalloc3(tlength,&remoterows,tlength,&rrow_ranks,tlength,&rrow_isids);
 79:   PetscCalloc3(size,&toranks,2*size,&tosizes,size,&tosizes_temp);
 80:   nrrows = 0;
 81:   for (i=0; i<nidx; i++) {
 82:     length_i  = length[i];
 83:     indices_i = indices[i];
 84:     for (j=0; j<length_i; j++) {
 85:       owner = -1;
 86:       PetscLayoutFindOwner(rmap,indices_i[j],&owner);
 87:       /* remote processors */
 88:       if (owner != rank) {
 89:         tosizes_temp[owner]++; /* number of rows to owner */
 90:         rrow_ranks[nrrows]  = owner; /* processor */
 91:         rrow_isids[nrrows]   = i; /* is id */
 92:         remoterows[nrrows++] = indices_i[j]; /* row */
 93:       }
 94:     }
 95:     ISRestoreIndices(is[i],&indices[i]);
 96:   }
 97:   PetscFree2(*(PetscInt***)&indices,length);
 98:   /* test if we need to exchange messages
 99:    * generally speaking, we do not need to exchange
100:    * data when overlap is 1
101:    * */
102:   MPIU_Allreduce(&nrrows,&reducednrrows,1,MPIU_INT,MPIU_MAX,comm);
103:   /* we do not have any messages
104:    * It usually corresponds to overlap 1
105:    * */
106:   if (!reducednrrows) {
107:     PetscFree3(toranks,tosizes,tosizes_temp);
108:     PetscFree3(remoterows,rrow_ranks,rrow_isids);
109:     MatIncreaseOverlap_MPIAIJ_Local_Scalable(mat,nidx,is);
110:     return(0);
111:   }
112:   nto = 0;
113:   /* send sizes and ranks for building a two-sided communcation */
114:   for (i=0; i<size; i++) {
115:     if (tosizes_temp[i]) {
116:       tosizes[nto*2]  = tosizes_temp[i]*2; /* size */
117:       tosizes_temp[i] = nto; /* a map from processor to index */
118:       toranks[nto++]  = i; /* processor */
119:     }
120:   }
121:   PetscMalloc1(nto+1,&toffsets);
122:   toffsets[0] = 0;
123:   for (i=0; i<nto; i++) {
124:     toffsets[i+1]  = toffsets[i]+tosizes[2*i]; /* offsets */
125:     tosizes[2*i+1] = toffsets[i]; /* offsets to send */
126:   }
127:   /* send information to other processors */
128:   PetscCommBuildTwoSided(comm,2,MPIU_INT,nto,toranks,tosizes,&nfrom,&fromranks,&fromsizes);
129:   nrecvrows = 0;
130:   for (i=0; i<nfrom; i++) nrecvrows += fromsizes[2*i];
131:   PetscMalloc1(nrecvrows,&remote);
132:   nrecvrows = 0;
133:   for (i=0; i<nfrom; i++) {
134:     for (j=0; j<fromsizes[2*i]; j++) {
135:       remote[nrecvrows].rank    = fromranks[i];
136:       remote[nrecvrows++].index = fromsizes[2*i+1]+j;
137:     }
138:   }
139:   PetscSFCreate(comm,&sf);
140:   PetscSFSetGraph(sf,nrecvrows,nrecvrows,NULL,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);
141:   /* use two-sided communication by default since OPENMPI has some bugs for one-sided one */
142:   PetscSFSetType(sf,PETSCSFBASIC);
143:   PetscSFSetFromOptions(sf);
144:   /* message pair <no of is, row>  */
145:   PetscCalloc2(2*nrrows,&todata,nrecvrows,&fromdata);
146:   for (i=0; i<nrrows; i++) {
147:     owner = rrow_ranks[i]; /* processor */
148:     j     = tosizes_temp[owner]; /* index */
149:     todata[toffsets[j]++] = rrow_isids[i];
150:     todata[toffsets[j]++] = remoterows[i];
151:   }
152:   PetscFree3(toranks,tosizes,tosizes_temp);
153:   PetscFree3(remoterows,rrow_ranks,rrow_isids);
154:   PetscFree(toffsets);
155:   PetscSFBcastBegin(sf,MPIU_INT,todata,fromdata);
156:   PetscSFBcastEnd(sf,MPIU_INT,todata,fromdata);
157:   PetscSFDestroy(&sf);
158:   /* send rows belonging to the remote so that then we could get the overlapping data back */
159:   MatIncreaseOverlap_MPIAIJ_Send_Scalable(mat,nidx,nfrom,fromranks,fromsizes,fromdata,&sbsizes,&sbdata);
160:   PetscFree2(todata,fromdata);
161:   PetscFree(fromsizes);
162:   PetscCommBuildTwoSided(comm,2,MPIU_INT,nfrom,fromranks,sbsizes,&nto,&toranks,&tosizes);
163:   PetscFree(fromranks);
164:   nrecvrows = 0;
165:   for (i=0; i<nto; i++) nrecvrows += tosizes[2*i];
166:   PetscCalloc1(nrecvrows,&todata);
167:   PetscMalloc1(nrecvrows,&remote);
168:   nrecvrows = 0;
169:   for (i=0; i<nto; i++) {
170:     for (j=0; j<tosizes[2*i]; j++) {
171:       remote[nrecvrows].rank    = toranks[i];
172:       remote[nrecvrows++].index = tosizes[2*i+1]+j;
173:     }
174:   }
175:   PetscSFCreate(comm,&sf);
176:   PetscSFSetGraph(sf,nrecvrows,nrecvrows,NULL,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);
177:   /* use two-sided communication by default since OPENMPI has some bugs for one-sided one */
178:   PetscSFSetType(sf,PETSCSFBASIC);
179:   PetscSFSetFromOptions(sf);
180:   /* overlap communication and computation */
181:   PetscSFBcastBegin(sf,MPIU_INT,sbdata,todata);
182:   MatIncreaseOverlap_MPIAIJ_Local_Scalable(mat,nidx,is);
183:   PetscSFBcastEnd(sf,MPIU_INT,sbdata,todata);
184:   PetscSFDestroy(&sf);
185:   PetscFree2(sbdata,sbsizes);
186:   MatIncreaseOverlap_MPIAIJ_Receive_Scalable(mat,nidx,is,nrecvrows,todata);
187:   PetscFree(toranks);
188:   PetscFree(tosizes);
189:   PetscFree(todata);
190:   return(0);
191: }

193: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Receive_Scalable(Mat mat,PetscInt nidx, IS is[], PetscInt nrecvs, PetscInt *recvdata)
194: {
195:   PetscInt         *isz,isz_i,i,j,is_id, data_size;
196:   PetscInt          col,lsize,max_lsize,*indices_temp, *indices_i;
197:   const PetscInt   *indices_i_temp;
198:   MPI_Comm         *iscomms;
199:   PetscErrorCode    ierr;

202:   max_lsize = 0;
203:   PetscMalloc1(nidx,&isz);
204:   for (i=0; i<nidx; i++){
205:     ISGetLocalSize(is[i],&lsize);
206:     max_lsize = lsize>max_lsize ? lsize:max_lsize;
207:     isz[i]    = lsize;
208:   }
209:   PetscMalloc2((max_lsize+nrecvs)*nidx,&indices_temp,nidx,&iscomms);
210:   for (i=0; i<nidx; i++){
211:     PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]),&iscomms[i],NULL);
212:     ISGetIndices(is[i],&indices_i_temp);
213:     PetscArraycpy(indices_temp+i*(max_lsize+nrecvs),indices_i_temp, isz[i]);
214:     ISRestoreIndices(is[i],&indices_i_temp);
215:     ISDestroy(&is[i]);
216:   }
217:   /* retrieve information to get row id and its overlap */
218:   for (i=0; i<nrecvs; ){
219:     is_id     = recvdata[i++];
220:     data_size = recvdata[i++];
221:     indices_i = indices_temp+(max_lsize+nrecvs)*is_id;
222:     isz_i     = isz[is_id];
223:     for (j=0; j< data_size; j++){
224:       col = recvdata[i++];
225:       indices_i[isz_i++] = col;
226:     }
227:     isz[is_id] = isz_i;
228:   }
229:   /* remove duplicate entities */
230:   for (i=0; i<nidx; i++){
231:     indices_i = indices_temp+(max_lsize+nrecvs)*i;
232:     isz_i     = isz[i];
233:     PetscSortRemoveDupsInt(&isz_i,indices_i);
234:     ISCreateGeneral(iscomms[i],isz_i,indices_i,PETSC_COPY_VALUES,&is[i]);
235:     PetscCommDestroy(&iscomms[i]);
236:   }
237:   PetscFree(isz);
238:   PetscFree2(indices_temp,iscomms);
239:   return(0);
240: }

242: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Send_Scalable(Mat mat,PetscInt nidx, PetscMPIInt nfrom,PetscMPIInt *fromranks,PetscInt *fromsizes, PetscInt *fromrows, PetscInt **sbrowsizes, PetscInt **sbrows)
243: {
244:   PetscLayout       rmap,cmap;
245:   PetscInt          i,j,k,l,*rows_i,*rows_data_ptr,**rows_data,max_fszs,rows_pos,*rows_pos_i;
246:   PetscInt          is_id,tnz,an,bn,rstart,cstart,row,start,end,col,totalrows,*sbdata;
247:   PetscInt         *indv_counts,indvc_ij,*sbsizes,*indices_tmp,*offsets;
248:   const PetscInt   *gcols,*ai,*aj,*bi,*bj;
249:   Mat               amat,bmat;
250:   PetscMPIInt       rank;
251:   PetscBool         done;
252:   MPI_Comm          comm;
253:   PetscErrorCode    ierr;

256:   PetscObjectGetComm((PetscObject)mat,&comm);
257:   MPI_Comm_rank(comm,&rank);
258:   MatMPIAIJGetSeqAIJ(mat,&amat,&bmat,&gcols);
259:   /* Even if the mat is symmetric, we still assume it is not symmetric */
260:   MatGetRowIJ(amat,0,PETSC_FALSE,PETSC_FALSE,&an,&ai,&aj,&done);
261:   if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"can not get row IJ \n");
262:   MatGetRowIJ(bmat,0,PETSC_FALSE,PETSC_FALSE,&bn,&bi,&bj,&done);
263:   if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"can not get row IJ \n");
264:   /* total number of nonzero values is used to estimate the memory usage in the next step */
265:   tnz  = ai[an]+bi[bn];
266:   MatGetLayouts(mat,&rmap,&cmap);
267:   PetscLayoutGetRange(rmap,&rstart,NULL);
268:   PetscLayoutGetRange(cmap,&cstart,NULL);
269:   /* to find the longest message */
270:   max_fszs = 0;
271:   for (i=0; i<nfrom; i++) max_fszs = fromsizes[2*i]>max_fszs ? fromsizes[2*i]:max_fszs;
272:   /* better way to estimate number of nonzero in the mat??? */
273:   PetscCalloc5(max_fszs*nidx,&rows_data_ptr,nidx,&rows_data,nidx,&rows_pos_i,nfrom*nidx,&indv_counts,tnz,&indices_tmp);
274:   for (i=0; i<nidx; i++) rows_data[i] = rows_data_ptr+max_fszs*i;
275:   rows_pos  = 0;
276:   totalrows = 0;
277:   for (i=0; i<nfrom; i++){
278:     PetscArrayzero(rows_pos_i,nidx);
279:     /* group data together */
280:     for (j=0; j<fromsizes[2*i]; j+=2){
281:       is_id                       = fromrows[rows_pos++];/* no of is */
282:       rows_i                      = rows_data[is_id];
283:       rows_i[rows_pos_i[is_id]++] = fromrows[rows_pos++];/* row */
284:     }
285:     /* estimate a space to avoid multiple allocations  */
286:     for (j=0; j<nidx; j++){
287:       indvc_ij = 0;
288:       rows_i   = rows_data[j];
289:       for (l=0; l<rows_pos_i[j]; l++){
290:         row    = rows_i[l]-rstart;
291:         start  = ai[row];
292:         end    = ai[row+1];
293:         for (k=start; k<end; k++){ /* Amat */
294:           col = aj[k] + cstart;
295:           indices_tmp[indvc_ij++] = col;/* do not count the rows from the original rank */
296:         }
297:         start = bi[row];
298:         end   = bi[row+1];
299:         for (k=start; k<end; k++) { /* Bmat */
300:           col = gcols[bj[k]];
301:           indices_tmp[indvc_ij++] = col;
302:         }
303:       }
304:       PetscSortRemoveDupsInt(&indvc_ij,indices_tmp);
305:       indv_counts[i*nidx+j] = indvc_ij;
306:       totalrows            += indvc_ij;
307:     }
308:   }
309:   /* message triple <no of is, number of rows, rows> */
310:   PetscCalloc2(totalrows+nidx*nfrom*2,&sbdata,2*nfrom,&sbsizes);
311:   totalrows = 0;
312:   rows_pos  = 0;
313:   /* use this code again */
314:   for (i=0;i<nfrom;i++){
315:     PetscArrayzero(rows_pos_i,nidx);
316:     for (j=0; j<fromsizes[2*i]; j+=2){
317:       is_id                       = fromrows[rows_pos++];
318:       rows_i                      = rows_data[is_id];
319:       rows_i[rows_pos_i[is_id]++] = fromrows[rows_pos++];
320:     }
321:     /* add data  */
322:     for (j=0; j<nidx; j++){
323:       if (!indv_counts[i*nidx+j]) continue;
324:       indvc_ij = 0;
325:       sbdata[totalrows++] = j;
326:       sbdata[totalrows++] = indv_counts[i*nidx+j];
327:       sbsizes[2*i]       += 2;
328:       rows_i              = rows_data[j];
329:       for (l=0; l<rows_pos_i[j]; l++){
330:         row   = rows_i[l]-rstart;
331:         start = ai[row];
332:         end   = ai[row+1];
333:         for (k=start; k<end; k++){ /* Amat */
334:           col = aj[k] + cstart;
335:           indices_tmp[indvc_ij++] = col;
336:         }
337:         start = bi[row];
338:         end   = bi[row+1];
339:         for (k=start; k<end; k++) { /* Bmat */
340:           col = gcols[bj[k]];
341:           indices_tmp[indvc_ij++] = col;
342:         }
343:       }
344:       PetscSortRemoveDupsInt(&indvc_ij,indices_tmp);
345:       sbsizes[2*i]  += indvc_ij;
346:       PetscArraycpy(sbdata+totalrows,indices_tmp,indvc_ij);
347:       totalrows += indvc_ij;
348:     }
349:   }
350:   PetscMalloc1(nfrom+1,&offsets);
351:   offsets[0] = 0;
352:   for (i=0; i<nfrom; i++){
353:     offsets[i+1]   = offsets[i] + sbsizes[2*i];
354:     sbsizes[2*i+1] = offsets[i];
355:   }
356:   PetscFree(offsets);
357:   if (sbrowsizes) *sbrowsizes = sbsizes;
358:   if (sbrows) *sbrows = sbdata;
359:   PetscFree5(rows_data_ptr,rows_data,rows_pos_i,indv_counts,indices_tmp);
360:   MatRestoreRowIJ(amat,0,PETSC_FALSE,PETSC_FALSE,&an,&ai,&aj,&done);
361:   MatRestoreRowIJ(bmat,0,PETSC_FALSE,PETSC_FALSE,&bn,&bi,&bj,&done);
362:   return(0);
363: }

365: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Local_Scalable(Mat mat,PetscInt nidx, IS is[])
366: {
367:   const PetscInt   *gcols,*ai,*aj,*bi,*bj, *indices;
368:   PetscInt          tnz,an,bn,i,j,row,start,end,rstart,cstart,col,k,*indices_temp;
369:   PetscInt          lsize,lsize_tmp,owner;
370:   PetscMPIInt       rank;
371:   Mat               amat,bmat;
372:   PetscBool         done;
373:   PetscLayout       cmap,rmap;
374:   MPI_Comm          comm;
375:   PetscErrorCode    ierr;

378:   PetscObjectGetComm((PetscObject)mat,&comm);
379:   MPI_Comm_rank(comm,&rank);
380:   MatMPIAIJGetSeqAIJ(mat,&amat,&bmat,&gcols);
381:   MatGetRowIJ(amat,0,PETSC_FALSE,PETSC_FALSE,&an,&ai,&aj,&done);
382:   if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"can not get row IJ \n");
383:   MatGetRowIJ(bmat,0,PETSC_FALSE,PETSC_FALSE,&bn,&bi,&bj,&done);
384:   if (!done) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"can not get row IJ \n");
385:   /* is it a safe way to compute number of nonzero values ? */
386:   tnz  = ai[an]+bi[bn];
387:   MatGetLayouts(mat,&rmap,&cmap);
388:   PetscLayoutGetRange(rmap,&rstart,NULL);
389:   PetscLayoutGetRange(cmap,&cstart,NULL);
390:   /* it is a better way to estimate memory than the old implementation
391:    * where global size of matrix is used
392:    * */
393:   PetscMalloc1(tnz,&indices_temp);
394:   for (i=0; i<nidx; i++) {
395:     MPI_Comm iscomm;

397:     ISGetLocalSize(is[i],&lsize);
398:     ISGetIndices(is[i],&indices);
399:     lsize_tmp = 0;
400:     for (j=0; j<lsize; j++) {
401:       owner = -1;
402:       row   = indices[j];
403:       PetscLayoutFindOwner(rmap,row,&owner);
404:       if (owner != rank) continue;
405:       /* local number */
406:       row  -= rstart;
407:       start = ai[row];
408:       end   = ai[row+1];
409:       for (k=start; k<end; k++) { /* Amat */
410:         col = aj[k] + cstart;
411:         indices_temp[lsize_tmp++] = col;
412:       }
413:       start = bi[row];
414:       end   = bi[row+1];
415:       for (k=start; k<end; k++) { /* Bmat */
416:         col = gcols[bj[k]];
417:         indices_temp[lsize_tmp++] = col;
418:       }
419:     }
420:    ISRestoreIndices(is[i],&indices);
421:    PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]),&iscomm,NULL);
422:    ISDestroy(&is[i]);
423:    PetscSortRemoveDupsInt(&lsize_tmp,indices_temp);
424:    ISCreateGeneral(iscomm,lsize_tmp,indices_temp,PETSC_COPY_VALUES,&is[i]);
425:    PetscCommDestroy(&iscomm);
426:   }
427:   PetscFree(indices_temp);
428:   MatRestoreRowIJ(amat,0,PETSC_FALSE,PETSC_FALSE,&an,&ai,&aj,&done);
429:   MatRestoreRowIJ(bmat,0,PETSC_FALSE,PETSC_FALSE,&bn,&bi,&bj,&done);
430:   return(0);
431: }


434: /*
435:   Sample message format:
436:   If a processor A wants processor B to process some elements corresponding
437:   to index sets is[1],is[5]
438:   mesg [0] = 2   (no of index sets in the mesg)
439:   -----------
440:   mesg [1] = 1 => is[1]
441:   mesg [2] = sizeof(is[1]);
442:   -----------
443:   mesg [3] = 5  => is[5]
444:   mesg [4] = sizeof(is[5]);
445:   -----------
446:   mesg [5]
447:   mesg [n]  datas[1]
448:   -----------
449:   mesg[n+1]
450:   mesg[m]  data(is[5])
451:   -----------

453:   Notes:
454:   nrqs - no of requests sent (or to be sent out)
455:   nrqr - no of requests recieved (which have to be or which have been processed
456: */
457: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Once(Mat C,PetscInt imax,IS is[])
458: {
459:   Mat_MPIAIJ     *c = (Mat_MPIAIJ*)C->data;
460:   PetscMPIInt    *w1,*w2,nrqr,*w3,*w4,*onodes1,*olengths1,*onodes2,*olengths2;
461:   const PetscInt **idx,*idx_i;
462:   PetscInt       *n,**data,len;
463: #if defined(PETSC_USE_CTABLE)
464:   PetscTable     *table_data,table_data_i;
465:   PetscInt       *tdata,tcount,tcount_max;
466: #else
467:   PetscInt       *data_i,*d_p;
468: #endif
470:   PetscMPIInt    size,rank,tag1,tag2;
471:   PetscInt       M,i,j,k,**rbuf,row,proc = 0,nrqs,msz,**outdat,**ptr;
472:   PetscInt       *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2;
473:   PetscBT        *table;
474:   MPI_Comm       comm;
475:   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
476:   MPI_Status     *s_status,*recv_status;
477:   MPI_Comm       *iscomms;
478:   char           *t_p;

481:   PetscObjectGetComm((PetscObject)C,&comm);
482:   size = c->size;
483:   rank = c->rank;
484:   M    = C->rmap->N;

486:   PetscObjectGetNewTag((PetscObject)C,&tag1);
487:   PetscObjectGetNewTag((PetscObject)C,&tag2);

489:   PetscMalloc2(imax,(PetscInt***)&idx,imax,&n);

491:   for (i=0; i<imax; i++) {
492:     ISGetIndices(is[i],&idx[i]);
493:     ISGetLocalSize(is[i],&n[i]);
494:   }

496:   /* evaluate communication - mesg to who,length of mesg, and buffer space
497:      required. Based on this, buffers are allocated, and data copied into them  */
498:   PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);
499:   for (i=0; i<imax; i++) {
500:     PetscArrayzero(w4,size); /* initialise work vector*/
501:     idx_i = idx[i];
502:     len   = n[i];
503:     for (j=0; j<len; j++) {
504:       row = idx_i[j];
505:       if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries");
506:       PetscLayoutFindOwner(C->rmap,row,&proc);
507:       w4[proc]++;
508:     }
509:     for (j=0; j<size; j++) {
510:       if (w4[j]) { w1[j] += w4[j]; w3[j]++;}
511:     }
512:   }

514:   nrqs     = 0;              /* no of outgoing messages */
515:   msz      = 0;              /* total mesg length (for all proc */
516:   w1[rank] = 0;              /* no mesg sent to intself */
517:   w3[rank] = 0;
518:   for (i=0; i<size; i++) {
519:     if (w1[i])  {w2[i] = 1; nrqs++;} /* there exists a message to proc i */
520:   }
521:   /* pa - is list of processors to communicate with */
522:   PetscMalloc1(nrqs+1,&pa);
523:   for (i=0,j=0; i<size; i++) {
524:     if (w1[i]) {pa[j] = i; j++;}
525:   }

527:   /* Each message would have a header = 1 + 2*(no of IS) + data */
528:   for (i=0; i<nrqs; i++) {
529:     j      = pa[i];
530:     w1[j] += w2[j] + 2*w3[j];
531:     msz   += w1[j];
532:   }

534:   /* Determine the number of messages to expect, their lengths, from from-ids */
535:   PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
536:   PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);

538:   /* Now post the Irecvs corresponding to these messages */
539:   PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);

541:   /* Allocate Memory for outgoing messages */
542:   PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);
543:   PetscArrayzero(outdat,size);
544:   PetscArrayzero(ptr,size);

546:   {
547:     PetscInt *iptr = tmp,ict  = 0;
548:     for (i=0; i<nrqs; i++) {
549:       j         = pa[i];
550:       iptr     +=  ict;
551:       outdat[j] = iptr;
552:       ict       = w1[j];
553:     }
554:   }

556:   /* Form the outgoing messages */
557:   /* plug in the headers */
558:   for (i=0; i<nrqs; i++) {
559:     j            = pa[i];
560:     outdat[j][0] = 0;
561:     PetscArrayzero(outdat[j]+1,2*w3[j]);
562:     ptr[j]       = outdat[j] + 2*w3[j] + 1;
563:   }

565:   /* Memory for doing local proc's work */
566:   {
567:     PetscInt M_BPB_imax = 0;
568: #if defined(PETSC_USE_CTABLE)
569:     PetscIntMultError((M/PETSC_BITS_PER_BYTE+1),imax, &M_BPB_imax);
570:     PetscMalloc1(imax,&table_data);
571:     for (i=0; i<imax; i++) {
572:       PetscTableCreate(n[i]+1,M+1,&table_data[i]);
573:     }
574:     PetscCalloc4(imax,&table, imax,&data, imax,&isz, M_BPB_imax,&t_p);
575:     for (i=0; i<imax; i++) {
576:       table[i] = t_p + (M/PETSC_BITS_PER_BYTE+1)*i;
577:     }
578: #else
579:     PetscInt Mimax = 0;
580:     PetscIntMultError(M,imax, &Mimax);
581:     PetscIntMultError((M/PETSC_BITS_PER_BYTE+1),imax, &M_BPB_imax);
582:     PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mimax,&d_p, M_BPB_imax,&t_p);
583:     for (i=0; i<imax; i++) {
584:       table[i] = t_p + (M/PETSC_BITS_PER_BYTE+1)*i;
585:       data[i]  = d_p + M*i;
586:     }
587: #endif
588:   }

590:   /* Parse the IS and update local tables and the outgoing buf with the data */
591:   {
592:     PetscInt n_i,isz_i,*outdat_j,ctr_j;
593:     PetscBT  table_i;

595:     for (i=0; i<imax; i++) {
596:       PetscArrayzero(ctr,size);
597:       n_i     = n[i];
598:       table_i = table[i];
599:       idx_i   = idx[i];
600: #if defined(PETSC_USE_CTABLE)
601:       table_data_i = table_data[i];
602: #else
603:       data_i  = data[i];
604: #endif
605:       isz_i   = isz[i];
606:       for (j=0; j<n_i; j++) {   /* parse the indices of each IS */
607:         row  = idx_i[j];
608:         PetscLayoutFindOwner(C->rmap,row,&proc);
609:         if (proc != rank) { /* copy to the outgoing buffer */
610:           ctr[proc]++;
611:           *ptr[proc] = row;
612:           ptr[proc]++;
613:         } else if (!PetscBTLookupSet(table_i,row)) {
614: #if defined(PETSC_USE_CTABLE)
615:           PetscTableAdd(table_data_i,row+1,isz_i+1,INSERT_VALUES);
616: #else
617:           data_i[isz_i] = row; /* Update the local table */
618: #endif
619:           isz_i++;
620:         }
621:       }
622:       /* Update the headers for the current IS */
623:       for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
624:         if ((ctr_j = ctr[j])) {
625:           outdat_j        = outdat[j];
626:           k               = ++outdat_j[0];
627:           outdat_j[2*k]   = ctr_j;
628:           outdat_j[2*k-1] = i;
629:         }
630:       }
631:       isz[i] = isz_i;
632:     }
633:   }

635:   /*  Now  post the sends */
636:   PetscMalloc1(nrqs+1,&s_waits1);
637:   for (i=0; i<nrqs; ++i) {
638:     j    = pa[i];
639:     MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);
640:   }

642:   /* No longer need the original indices */
643:   PetscMalloc1(imax,&iscomms);
644:   for (i=0; i<imax; ++i) {
645:     ISRestoreIndices(is[i],idx+i);
646:     PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]),&iscomms[i],NULL);
647:   }
648:   PetscFree2(*(PetscInt***)&idx,n);

650:   for (i=0; i<imax; ++i) {
651:     ISDestroy(&is[i]);
652:   }

654:   /* Do Local work */
655: #if defined(PETSC_USE_CTABLE)
656:   MatIncreaseOverlap_MPIAIJ_Local(C,imax,table,isz,NULL,table_data);
657: #else
658:   MatIncreaseOverlap_MPIAIJ_Local(C,imax,table,isz,data,NULL);
659: #endif

661:   /* Receive messages */
662:   PetscMalloc1(nrqr+1,&recv_status);
663:   if (nrqr) {MPI_Waitall(nrqr,r_waits1,recv_status);}

665:   PetscMalloc1(nrqs+1,&s_status);
666:   if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status);}

668:   /* Phase 1 sends are complete - deallocate buffers */
669:   PetscFree4(outdat,ptr,tmp,ctr);
670:   PetscFree4(w1,w2,w3,w4);

672:   PetscMalloc1(nrqr+1,&xdata);
673:   PetscMalloc1(nrqr+1,&isz1);
674:   MatIncreaseOverlap_MPIAIJ_Receive(C,nrqr,rbuf,xdata,isz1);
675:   PetscFree(rbuf[0]);
676:   PetscFree(rbuf);


679:   /* Send the data back */
680:   /* Do a global reduction to know the buffer space req for incoming messages */
681:   {
682:     PetscMPIInt *rw1;

684:     PetscCalloc1(size,&rw1);

686:     for (i=0; i<nrqr; ++i) {
687:       proc = recv_status[i].MPI_SOURCE;

689:       if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch");
690:       rw1[proc] = isz1[i];
691:     }
692:     PetscFree(onodes1);
693:     PetscFree(olengths1);

695:     /* Determine the number of messages to expect, their lengths, from from-ids */
696:     PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);
697:     PetscFree(rw1);
698:   }
699:   /* Now post the Irecvs corresponding to these messages */
700:   PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);

702:   /* Now  post the sends */
703:   PetscMalloc1(nrqr+1,&s_waits2);
704:   for (i=0; i<nrqr; ++i) {
705:     j    = recv_status[i].MPI_SOURCE;
706:     MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);
707:   }

709:   /* receive work done on other processors */
710:   {
711:     PetscInt    is_no,ct1,max,*rbuf2_i,isz_i,jmax;
712:     PetscMPIInt idex;
713:     PetscBT     table_i;
714:     MPI_Status  *status2;

716:     PetscMalloc1((PetscMax(nrqr,nrqs)+1),&status2);
717:     for (i=0; i<nrqs; ++i) {
718:       MPI_Waitany(nrqs,r_waits2,&idex,status2+i);
719:       /* Process the message */
720:       rbuf2_i = rbuf2[idex];
721:       ct1     = 2*rbuf2_i[0]+1;
722:       jmax    = rbuf2[idex][0];
723:       for (j=1; j<=jmax; j++) {
724:         max     = rbuf2_i[2*j];
725:         is_no   = rbuf2_i[2*j-1];
726:         isz_i   = isz[is_no];
727:         table_i = table[is_no];
728: #if defined(PETSC_USE_CTABLE)
729:         table_data_i = table_data[is_no];
730: #else
731:         data_i  = data[is_no];
732: #endif
733:         for (k=0; k<max; k++,ct1++) {
734:           row = rbuf2_i[ct1];
735:           if (!PetscBTLookupSet(table_i,row)) {
736: #if defined(PETSC_USE_CTABLE)
737:             PetscTableAdd(table_data_i,row+1,isz_i+1,INSERT_VALUES);
738: #else
739:             data_i[isz_i] = row;
740: #endif
741:             isz_i++;
742:           }
743:         }
744:         isz[is_no] = isz_i;
745:       }
746:     }

748:     if (nrqr) {MPI_Waitall(nrqr,s_waits2,status2);}
749:     PetscFree(status2);
750:   }

752: #if defined(PETSC_USE_CTABLE)
753:   tcount_max = 0;
754:   for (i=0; i<imax; ++i) {
755:     table_data_i = table_data[i];
756:     PetscTableGetCount(table_data_i,&tcount);
757:     if (tcount_max < tcount) tcount_max = tcount;
758:   }
759:   PetscMalloc1(tcount_max+1,&tdata);
760: #endif

762:   for (i=0; i<imax; ++i) {
763: #if defined(PETSC_USE_CTABLE)
764:     PetscTablePosition tpos;
765:     table_data_i = table_data[i];

767:     PetscTableGetHeadPosition(table_data_i,&tpos);
768:     while (tpos) {
769:       PetscTableGetNext(table_data_i,&tpos,&k,&j);
770:       tdata[--j] = --k;
771:     }
772:     ISCreateGeneral(iscomms[i],isz[i],tdata,PETSC_COPY_VALUES,is+i);
773: #else
774:     ISCreateGeneral(iscomms[i],isz[i],data[i],PETSC_COPY_VALUES,is+i);
775: #endif
776:     PetscCommDestroy(&iscomms[i]);
777:   }

779:   PetscFree(iscomms);
780:   PetscFree(onodes2);
781:   PetscFree(olengths2);

783:   PetscFree(pa);
784:   PetscFree(rbuf2[0]);
785:   PetscFree(rbuf2);
786:   PetscFree(s_waits1);
787:   PetscFree(r_waits1);
788:   PetscFree(s_waits2);
789:   PetscFree(r_waits2);
790:   PetscFree(s_status);
791:   PetscFree(recv_status);
792:   PetscFree(xdata[0]);
793:   PetscFree(xdata);
794:   PetscFree(isz1);
795: #if defined(PETSC_USE_CTABLE)
796:   for (i=0; i<imax; i++) {
797:     PetscTableDestroy((PetscTable*)&table_data[i]);
798:   }
799:   PetscFree(table_data);
800:   PetscFree(tdata);
801:   PetscFree4(table,data,isz,t_p);
802: #else
803:   PetscFree5(table,data,isz,d_p,t_p);
804: #endif
805:   return(0);
806: }

808: /*
809:    MatIncreaseOverlap_MPIAIJ_Local - Called by MatincreaseOverlap, to do
810:        the work on the local processor.

812:      Inputs:
813:       C      - MAT_MPIAIJ;
814:       imax - total no of index sets processed at a time;
815:       table  - an array of char - size = m bits.

817:      Output:
818:       isz    - array containing the count of the solution elements corresponding
819:                to each index set;
820:       data or table_data  - pointer to the solutions
821: */
822: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data,PetscTable *table_data)
823: {
824:   Mat_MPIAIJ *c = (Mat_MPIAIJ*)C->data;
825:   Mat        A  = c->A,B = c->B;
826:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data;
827:   PetscInt   start,end,val,max,rstart,cstart,*ai,*aj;
828:   PetscInt   *bi,*bj,*garray,i,j,k,row,isz_i;
829:   PetscBT    table_i;
830: #if defined(PETSC_USE_CTABLE)
831:   PetscTable         table_data_i;
832:   PetscErrorCode     ierr;
833:   PetscTablePosition tpos;
834:   PetscInt           tcount,*tdata;
835: #else
836:   PetscInt           *data_i;
837: #endif

840:   rstart = C->rmap->rstart;
841:   cstart = C->cmap->rstart;
842:   ai     = a->i;
843:   aj     = a->j;
844:   bi     = b->i;
845:   bj     = b->j;
846:   garray = c->garray;

848:   for (i=0; i<imax; i++) {
849: #if defined(PETSC_USE_CTABLE)
850:     /* copy existing entries of table_data_i into tdata[] */
851:     table_data_i = table_data[i];
852:     PetscTableGetCount(table_data_i,&tcount);
853:     if (tcount != isz[i]) SETERRQ3(PETSC_COMM_SELF,0," tcount %d != isz[%d] %d",tcount,i,isz[i]);

855:     PetscMalloc1(tcount,&tdata);
856:     PetscTableGetHeadPosition(table_data_i,&tpos);
857:     while (tpos) {
858:       PetscTableGetNext(table_data_i,&tpos,&row,&j);
859:       tdata[--j] = --row;
860:       if (j > tcount - 1) SETERRQ2(PETSC_COMM_SELF,0," j %d >= tcount %d",j,tcount);
861:     }
862: #else
863:     data_i  = data[i];
864: #endif
865:     table_i = table[i];
866:     isz_i   = isz[i];
867:     max     = isz[i];

869:     for (j=0; j<max; j++) {
870: #if defined(PETSC_USE_CTABLE)
871:       row   = tdata[j] - rstart;
872: #else
873:       row   = data_i[j] - rstart;
874: #endif
875:       start = ai[row];
876:       end   = ai[row+1];
877:       for (k=start; k<end; k++) { /* Amat */
878:         val = aj[k] + cstart;
879:         if (!PetscBTLookupSet(table_i,val)) {
880: #if defined(PETSC_USE_CTABLE)
881:           PetscTableAdd(table_data_i,val+1,isz_i+1,INSERT_VALUES);
882: #else
883:           data_i[isz_i] = val;
884: #endif
885:           isz_i++;
886:         }
887:       }
888:       start = bi[row];
889:       end   = bi[row+1];
890:       for (k=start; k<end; k++) { /* Bmat */
891:         val = garray[bj[k]];
892:         if (!PetscBTLookupSet(table_i,val)) {
893: #if defined(PETSC_USE_CTABLE)
894:           PetscTableAdd(table_data_i,val+1,isz_i+1,INSERT_VALUES);
895: #else
896:           data_i[isz_i] = val;
897: #endif
898:           isz_i++;
899:         }
900:       }
901:     }
902:     isz[i] = isz_i;

904: #if defined(PETSC_USE_CTABLE)
905:     PetscFree(tdata);
906: #endif
907:   }
908:   return(0);
909: }

911: /*
912:       MatIncreaseOverlap_MPIAIJ_Receive - Process the recieved messages,
913:          and return the output

915:          Input:
916:            C    - the matrix
917:            nrqr - no of messages being processed.
918:            rbuf - an array of pointers to the recieved requests

920:          Output:
921:            xdata - array of messages to be sent back
922:            isz1  - size of each message

924:   For better efficiency perhaps we should malloc separately each xdata[i],
925: then if a remalloc is required we need only copy the data for that one row
926: rather then all previous rows as it is now where a single large chunck of
927: memory is used.

929: */
930: static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1)
931: {
932:   Mat_MPIAIJ     *c = (Mat_MPIAIJ*)C->data;
933:   Mat            A  = c->A,B = c->B;
934:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data;
936:   PetscInt       rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k;
937:   PetscInt       row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end;
938:   PetscInt       val,max1,max2,m,no_malloc =0,*tmp,new_estimate,ctr;
939:   PetscInt       *rbuf_i,kmax,rbuf_0;
940:   PetscBT        xtable;

943:   m      = C->rmap->N;
944:   rstart = C->rmap->rstart;
945:   cstart = C->cmap->rstart;
946:   ai     = a->i;
947:   aj     = a->j;
948:   bi     = b->i;
949:   bj     = b->j;
950:   garray = c->garray;


953:   for (i=0,ct=0,total_sz=0; i<nrqr; ++i) {
954:     rbuf_i =  rbuf[i];
955:     rbuf_0 =  rbuf_i[0];
956:     ct    += rbuf_0;
957:     for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j];
958:   }

960:   if (C->rmap->n) max1 = ct*(a->nz + b->nz)/C->rmap->n;
961:   else max1 = 1;
962:   mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1);
963:   PetscMalloc1(mem_estimate,&xdata[0]);
964:   ++no_malloc;
965:   PetscBTCreate(m,&xtable);
966:   PetscArrayzero(isz1,nrqr);

968:   ct3 = 0;
969:   for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */
970:     rbuf_i =  rbuf[i];
971:     rbuf_0 =  rbuf_i[0];
972:     ct1    =  2*rbuf_0+1;
973:     ct2    =  ct1;
974:     ct3   += ct1;
975:     for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/
976:       PetscBTMemzero(m,xtable);
977:       oct2 = ct2;
978:       kmax = rbuf_i[2*j];
979:       for (k=0; k<kmax; k++,ct1++) {
980:         row = rbuf_i[ct1];
981:         if (!PetscBTLookupSet(xtable,row)) {
982:           if (!(ct3 < mem_estimate)) {
983:             new_estimate = (PetscInt)(1.5*mem_estimate)+1;
984:             PetscMalloc1(new_estimate,&tmp);
985:             PetscArraycpy(tmp,xdata[0],mem_estimate);
986:             PetscFree(xdata[0]);
987:             xdata[0]     = tmp;
988:             mem_estimate = new_estimate; ++no_malloc;
989:             for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
990:           }
991:           xdata[i][ct2++] = row;
992:           ct3++;
993:         }
994:       }
995:       for (k=oct2,max2=ct2; k<max2; k++) {
996:         row   = xdata[i][k] - rstart;
997:         start = ai[row];
998:         end   = ai[row+1];
999:         for (l=start; l<end; l++) {
1000:           val = aj[l] + cstart;
1001:           if (!PetscBTLookupSet(xtable,val)) {
1002:             if (!(ct3 < mem_estimate)) {
1003:               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
1004:               PetscMalloc1(new_estimate,&tmp);
1005:               PetscArraycpy(tmp,xdata[0],mem_estimate);
1006:               PetscFree(xdata[0]);
1007:               xdata[0]     = tmp;
1008:               mem_estimate = new_estimate; ++no_malloc;
1009:               for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
1010:             }
1011:             xdata[i][ct2++] = val;
1012:             ct3++;
1013:           }
1014:         }
1015:         start = bi[row];
1016:         end   = bi[row+1];
1017:         for (l=start; l<end; l++) {
1018:           val = garray[bj[l]];
1019:           if (!PetscBTLookupSet(xtable,val)) {
1020:             if (!(ct3 < mem_estimate)) {
1021:               new_estimate = (PetscInt)(1.5*mem_estimate)+1;
1022:               PetscMalloc1(new_estimate,&tmp);
1023:               PetscArraycpy(tmp,xdata[0],mem_estimate);
1024:               PetscFree(xdata[0]);
1025:               xdata[0]     = tmp;
1026:               mem_estimate = new_estimate; ++no_malloc;
1027:               for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];
1028:             }
1029:             xdata[i][ct2++] = val;
1030:             ct3++;
1031:           }
1032:         }
1033:       }
1034:       /* Update the header*/
1035:       xdata[i][2*j]   = ct2 - oct2; /* Undo the vector isz1 and use only a var*/
1036:       xdata[i][2*j-1] = rbuf_i[2*j-1];
1037:     }
1038:     xdata[i][0] = rbuf_0;
1039:     xdata[i+1]  = xdata[i] + ct2;
1040:     isz1[i]     = ct2; /* size of each message */
1041:   }
1042:   PetscBTDestroy(&xtable);
1043:   PetscInfo3(C,"Allocated %D bytes, required %D bytes, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);
1044:   return(0);
1045: }
1046: /* -------------------------------------------------------------------------*/
1047: extern PetscErrorCode MatCreateSubMatrices_MPIAIJ_Local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*);
1048: /*
1049:     Every processor gets the entire matrix
1050: */
1051: PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat A,MatCreateSubMatrixOption flag,MatReuse scall,Mat *Bin[])
1052: {
1053:   Mat            B;
1054:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
1055:   Mat_SeqAIJ     *b,*ad = (Mat_SeqAIJ*)a->A->data,*bd = (Mat_SeqAIJ*)a->B->data;
1057:   PetscMPIInt    size,rank,*recvcounts = 0,*displs = 0;
1058:   PetscInt       sendcount,i,*rstarts = A->rmap->range,n,cnt,j;
1059:   PetscInt       m,*b_sendj,*garray = a->garray,*lens,*jsendbuf,*a_jsendbuf,*b_jsendbuf;
1060:   MatScalar      *sendbuf,*recvbuf,*a_sendbuf,*b_sendbuf;

1063:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
1064:   MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);

1066:   if (scall == MAT_INITIAL_MATRIX) {
1067:     /* ----------------------------------------------------------------
1068:          Tell every processor the number of nonzeros per row
1069:     */
1070:     PetscMalloc1(A->rmap->N,&lens);
1071:     for (i=A->rmap->rstart; i<A->rmap->rend; i++) {
1072:       lens[i] = ad->i[i-A->rmap->rstart+1] - ad->i[i-A->rmap->rstart] + bd->i[i-A->rmap->rstart+1] - bd->i[i-A->rmap->rstart];
1073:     }
1074:     PetscMalloc2(size,&recvcounts,size,&displs);
1075:     for (i=0; i<size; i++) {
1076:       recvcounts[i] = A->rmap->range[i+1] - A->rmap->range[i];
1077:       displs[i]     = A->rmap->range[i];
1078:     }
1079: #if defined(PETSC_HAVE_MPI_IN_PLACE)
1080:     MPI_Allgatherv(MPI_IN_PLACE,0,MPI_DATATYPE_NULL,lens,recvcounts,displs,MPIU_INT,PetscObjectComm((PetscObject)A));
1081: #else
1082:     sendcount = A->rmap->rend - A->rmap->rstart;
1083:     MPI_Allgatherv(lens+A->rmap->rstart,sendcount,MPIU_INT,lens,recvcounts,displs,MPIU_INT,PetscObjectComm((PetscObject)A));
1084: #endif
1085:     /* ---------------------------------------------------------------
1086:          Create the sequential matrix of the same type as the local block diagonal
1087:     */
1088:     MatCreate(PETSC_COMM_SELF,&B);
1089:     MatSetSizes(B,A->rmap->N,A->cmap->N,PETSC_DETERMINE,PETSC_DETERMINE);
1090:     MatSetBlockSizesFromMats(B,A,A);
1091:     MatSetType(B,((PetscObject)a->A)->type_name);
1092:     MatSeqAIJSetPreallocation(B,0,lens);
1093:     PetscCalloc1(2,Bin);
1094:     **Bin = B;
1095:     b     = (Mat_SeqAIJ*)B->data;

1097:     /*--------------------------------------------------------------------
1098:        Copy my part of matrix column indices over
1099:     */
1100:     sendcount  = ad->nz + bd->nz;
1101:     jsendbuf   = b->j + b->i[rstarts[rank]];
1102:     a_jsendbuf = ad->j;
1103:     b_jsendbuf = bd->j;
1104:     n          = A->rmap->rend - A->rmap->rstart;
1105:     cnt        = 0;
1106:     for (i=0; i<n; i++) {

1108:       /* put in lower diagonal portion */
1109:       m = bd->i[i+1] - bd->i[i];
1110:       while (m > 0) {
1111:         /* is it above diagonal (in bd (compressed) numbering) */
1112:         if (garray[*b_jsendbuf] > A->rmap->rstart + i) break;
1113:         jsendbuf[cnt++] = garray[*b_jsendbuf++];
1114:         m--;
1115:       }

1117:       /* put in diagonal portion */
1118:       for (j=ad->i[i]; j<ad->i[i+1]; j++) {
1119:         jsendbuf[cnt++] = A->rmap->rstart + *a_jsendbuf++;
1120:       }

1122:       /* put in upper diagonal portion */
1123:       while (m-- > 0) {
1124:         jsendbuf[cnt++] = garray[*b_jsendbuf++];
1125:       }
1126:     }
1127:     if (cnt != sendcount) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupted PETSc matrix: nz given %D actual nz %D",sendcount,cnt);

1129:     /*--------------------------------------------------------------------
1130:        Gather all column indices to all processors
1131:     */
1132:     for (i=0; i<size; i++) {
1133:       recvcounts[i] = 0;
1134:       for (j=A->rmap->range[i]; j<A->rmap->range[i+1]; j++) {
1135:         recvcounts[i] += lens[j];
1136:       }
1137:     }
1138:     displs[0] = 0;
1139:     for (i=1; i<size; i++) {
1140:       displs[i] = displs[i-1] + recvcounts[i-1];
1141:     }
1142: #if defined(PETSC_HAVE_MPI_IN_PLACE)
1143:     MPI_Allgatherv(MPI_IN_PLACE,0,MPI_DATATYPE_NULL,b->j,recvcounts,displs,MPIU_INT,PetscObjectComm((PetscObject)A));
1144: #else
1145:     MPI_Allgatherv(jsendbuf,sendcount,MPIU_INT,b->j,recvcounts,displs,MPIU_INT,PetscObjectComm((PetscObject)A));
1146: #endif
1147:     /*--------------------------------------------------------------------
1148:         Assemble the matrix into useable form (note numerical values not yet set)
1149:     */
1150:     /* set the b->ilen (length of each row) values */
1151:     PetscArraycpy(b->ilen,lens,A->rmap->N);
1152:     /* set the b->i indices */
1153:     b->i[0] = 0;
1154:     for (i=1; i<=A->rmap->N; i++) {
1155:       b->i[i] = b->i[i-1] + lens[i-1];
1156:     }
1157:     PetscFree(lens);
1158:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1159:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

1161:   } else {
1162:     B = **Bin;
1163:     b = (Mat_SeqAIJ*)B->data;
1164:   }

1166:   /*--------------------------------------------------------------------
1167:        Copy my part of matrix numerical values into the values location
1168:   */
1169:   if (flag == MAT_GET_VALUES) {
1170:     sendcount = ad->nz + bd->nz;
1171:     sendbuf   = b->a + b->i[rstarts[rank]];
1172:     a_sendbuf = ad->a;
1173:     b_sendbuf = bd->a;
1174:     b_sendj   = bd->j;
1175:     n         = A->rmap->rend - A->rmap->rstart;
1176:     cnt       = 0;
1177:     for (i=0; i<n; i++) {

1179:       /* put in lower diagonal portion */
1180:       m = bd->i[i+1] - bd->i[i];
1181:       while (m > 0) {
1182:         /* is it above diagonal (in bd (compressed) numbering) */
1183:         if (garray[*b_sendj] > A->rmap->rstart + i) break;
1184:         sendbuf[cnt++] = *b_sendbuf++;
1185:         m--;
1186:         b_sendj++;
1187:       }

1189:       /* put in diagonal portion */
1190:       for (j=ad->i[i]; j<ad->i[i+1]; j++) {
1191:         sendbuf[cnt++] = *a_sendbuf++;
1192:       }

1194:       /* put in upper diagonal portion */
1195:       while (m-- > 0) {
1196:         sendbuf[cnt++] = *b_sendbuf++;
1197:         b_sendj++;
1198:       }
1199:     }
1200:     if (cnt != sendcount) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupted PETSc matrix: nz given %D actual nz %D",sendcount,cnt);

1202:     /* -----------------------------------------------------------------
1203:        Gather all numerical values to all processors
1204:     */
1205:     if (!recvcounts) {
1206:       PetscMalloc2(size,&recvcounts,size,&displs);
1207:     }
1208:     for (i=0; i<size; i++) {
1209:       recvcounts[i] = b->i[rstarts[i+1]] - b->i[rstarts[i]];
1210:     }
1211:     displs[0] = 0;
1212:     for (i=1; i<size; i++) {
1213:       displs[i] = displs[i-1] + recvcounts[i-1];
1214:     }
1215:     recvbuf = b->a;
1216: #if defined(PETSC_HAVE_MPI_IN_PLACE)
1217:     MPI_Allgatherv(MPI_IN_PLACE,0,MPI_DATATYPE_NULL,recvbuf,recvcounts,displs,MPIU_SCALAR,PetscObjectComm((PetscObject)A));
1218: #else
1219:     MPI_Allgatherv(sendbuf,sendcount,MPIU_SCALAR,recvbuf,recvcounts,displs,MPIU_SCALAR,PetscObjectComm((PetscObject)A));
1220: #endif
1221:   }  /* endof (flag == MAT_GET_VALUES) */
1222:   PetscFree2(recvcounts,displs);

1224:   if (A->symmetric) {
1225:     MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1226:   } else if (A->hermitian) {
1227:     MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
1228:   } else if (A->structurally_symmetric) {
1229:     MatSetOption(B,MAT_STRUCTURALLY_SYMMETRIC,PETSC_TRUE);
1230:   }
1231:   return(0);
1232: }

1234: PetscErrorCode MatCreateSubMatrices_MPIAIJ_SingleIS_Local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool allcolumns,Mat *submats)
1235: {
1236:   Mat_MPIAIJ     *c = (Mat_MPIAIJ*)C->data;
1237:   Mat            submat,A = c->A,B = c->B;
1238:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data,*subc;
1239:   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,nzA,nzB;
1240:   PetscInt       cstart = C->cmap->rstart,cend = C->cmap->rend,rstart = C->rmap->rstart,*bmap = c->garray;
1241:   const PetscInt *icol,*irow;
1242:   PetscInt       nrow,ncol,start;
1244:   PetscMPIInt    rank,size,tag1,tag2,tag3,tag4,*w1,*w2,nrqr;
1245:   PetscInt       **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,ct3,**rbuf1,row,proc;
1246:   PetscInt       nrqs=0,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol,*iptr;
1247:   PetscInt       **rbuf3,*req_source1,*req_source2,**sbuf_aj,**rbuf2,max1,nnz;
1248:   PetscInt       *lens,rmax,ncols,*cols,Crow;
1249: #if defined(PETSC_USE_CTABLE)
1250:   PetscTable     cmap,rmap;
1251:   PetscInt       *cmap_loc,*rmap_loc;
1252: #else
1253:   PetscInt       *cmap,*rmap;
1254: #endif
1255:   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*sbuf1_i,*rbuf2_i,*rbuf3_i;
1256:   PetscInt       *cworkB,lwrite,*subcols,*row2proc;
1257:   PetscScalar    *vworkA,*vworkB,*a_a = a->a,*b_a = b->a,*subvals=NULL;
1258:   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
1259:   MPI_Request    *r_waits4,*s_waits3 = NULL,*s_waits4;
1260:   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3 = NULL,*s_status2;
1261:   MPI_Status     *r_status3 = NULL,*r_status4,*s_status4;
1262:   MPI_Comm       comm;
1263:   PetscScalar    **rbuf4,**sbuf_aa,*vals,*sbuf_aa_i,*rbuf4_i;
1264:   PetscMPIInt    *onodes1,*olengths1,idex,end;
1265:   Mat_SubSppt    *smatis1;
1266:   PetscBool      isrowsorted,iscolsorted;

1269:   if (ismax != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This routine only works when all processes have ismax=1");

1271:   PetscObjectGetComm((PetscObject)C,&comm);
1272:   size = c->size;
1273:   rank = c->rank;

1275:   ISSorted(iscol[0],&iscolsorted);
1276:   ISSorted(isrow[0],&isrowsorted);
1277:   ISGetIndices(isrow[0],&irow);
1278:   ISGetLocalSize(isrow[0],&nrow);
1279:   if (allcolumns) {
1280:     icol = NULL;
1281:     ncol = C->cmap->N;
1282:   } else {
1283:     ISGetIndices(iscol[0],&icol);
1284:     ISGetLocalSize(iscol[0],&ncol);
1285:   }

1287:   if (scall == MAT_INITIAL_MATRIX) {
1288:     PetscInt *sbuf2_i,*cworkA,lwrite,ctmp;

1290:     /* Get some new tags to keep the communication clean */
1291:     tag1 = ((PetscObject)C)->tag;
1292:     PetscObjectGetNewTag((PetscObject)C,&tag2);
1293:     PetscObjectGetNewTag((PetscObject)C,&tag3);

1295:     /* evaluate communication - mesg to who, length of mesg, and buffer space
1296:      required. Based on this, buffers are allocated, and data copied into them */
1297:     PetscCalloc2(size,&w1,size,&w2);
1298:     PetscMalloc1(nrow,&row2proc);

1300:     /* w1[proc] = num of rows owned by proc -- to be requested */
1301:     proc = 0;
1302:     nrqs = 0; /* num of outgoing messages */
1303:     for (j=0; j<nrow; j++) {
1304:       row  = irow[j];
1305:       if (!isrowsorted) proc = 0;
1306:       while (row >= C->rmap->range[proc+1]) proc++;
1307:       w1[proc]++;
1308:       row2proc[j] = proc; /* map row index to proc */

1310:       if (proc != rank && !w2[proc]) {
1311:         w2[proc] = 1; nrqs++;
1312:       }
1313:     }
1314:     w1[rank] = 0;  /* rows owned by self will not be requested */

1316:     PetscMalloc1(nrqs+1,&pa); /*(proc -array)*/
1317:     for (proc=0,j=0; proc<size; proc++) {
1318:       if (w1[proc]) { pa[j++] = proc;}
1319:     }

1321:     /* Each message would have a header = 1 + 2*(num of IS) + data (here,num of IS = 1) */
1322:     msz = 0;              /* total mesg length (for all procs) */
1323:     for (i=0; i<nrqs; i++) {
1324:       proc      = pa[i];
1325:       w1[proc] += 3;
1326:       msz      += w1[proc];
1327:     }
1328:     PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);

1330:     /* Determine nrqr, the number of messages to expect, their lengths, from from-ids */
1331:     /* if w2[proc]=1, a message of length w1[proc] will be sent to proc; */
1332:     PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);

1334:     /* Input: nrqs: nsend; nrqr: nrecv; w1: msg length to be sent;
1335:        Output: onodes1: recv node-ids; olengths1: corresponding recv message length */
1336:     PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);

1338:     /* Now post the Irecvs corresponding to these messages */
1339:     PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);

1341:     PetscFree(onodes1);
1342:     PetscFree(olengths1);

1344:     /* Allocate Memory for outgoing messages */
1345:     PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);
1346:     PetscArrayzero(sbuf1,size);
1347:     PetscArrayzero(ptr,size);

1349:     /* subf1[pa[0]] = tmp, subf1[pa[i]] = subf1[pa[i-1]] + w1[pa[i-1]] */
1350:     iptr = tmp;
1351:     for (i=0; i<nrqs; i++) {
1352:       proc        = pa[i];
1353:       sbuf1[proc] = iptr;
1354:       iptr       += w1[proc];
1355:     }

1357:     /* Form the outgoing messages */
1358:     /* Initialize the header space */
1359:     for (i=0; i<nrqs; i++) {
1360:       proc      = pa[i];
1361:       PetscArrayzero(sbuf1[proc],3);
1362:       ptr[proc] = sbuf1[proc] + 3;
1363:     }

1365:     /* Parse the isrow and copy data into outbuf */
1366:     PetscArrayzero(ctr,size);
1367:     for (j=0; j<nrow; j++) {  /* parse the indices of each IS */
1368:       proc = row2proc[j];
1369:       if (proc != rank) { /* copy to the outgoing buf*/
1370:         *ptr[proc] = irow[j];
1371:         ctr[proc]++; ptr[proc]++;
1372:       }
1373:     }

1375:     /* Update the headers for the current IS */
1376:     for (j=0; j<size; j++) { /* Can Optimise this loop too */
1377:       if ((ctr_j = ctr[j])) {
1378:         sbuf1_j        = sbuf1[j];
1379:         k              = ++sbuf1_j[0];
1380:         sbuf1_j[2*k]   = ctr_j;
1381:         sbuf1_j[2*k-1] = 0;
1382:       }
1383:     }

1385:     /* Now post the sends */
1386:     PetscMalloc1(nrqs+1,&s_waits1);
1387:     for (i=0; i<nrqs; ++i) {
1388:       proc = pa[i];
1389:       MPI_Isend(sbuf1[proc],w1[proc],MPIU_INT,proc,tag1,comm,s_waits1+i);
1390:     }

1392:     /* Post Receives to capture the buffer size */
1393:     PetscMalloc4(nrqs+1,&r_status2,nrqr+1,&s_waits2,nrqs+1,&r_waits2,nrqr+1,&s_status2);
1394:     PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);

1396:     rbuf2[0] = tmp + msz;
1397:     for (i=1; i<nrqs; ++i) rbuf2[i] = rbuf2[i-1] + w1[pa[i-1]];

1399:     for (i=0; i<nrqs; ++i) {
1400:       proc = pa[i];
1401:       MPI_Irecv(rbuf2[i],w1[proc],MPIU_INT,proc,tag2,comm,r_waits2+i);
1402:     }

1404:     PetscFree2(w1,w2);

1406:     /* Send to other procs the buf size they should allocate */
1407:     /* Receive messages*/
1408:     PetscMalloc1(nrqr+1,&r_status1);
1409:     PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);

1411:     MPI_Waitall(nrqr,r_waits1,r_status1);
1412:     for (i=0; i<nrqr; ++i) {
1413:       req_size[i] = 0;
1414:       rbuf1_i        = rbuf1[i];
1415:       start          = 2*rbuf1_i[0] + 1;
1416:       MPI_Get_count(r_status1+i,MPIU_INT,&end);
1417:       PetscMalloc1(end+1,&sbuf2[i]);
1418:       sbuf2_i        = sbuf2[i];
1419:       for (j=start; j<end; j++) {
1420:         k            = rbuf1_i[j] - rstart;
1421:         ncols        = ai[k+1] - ai[k] + bi[k+1] - bi[k];
1422:         sbuf2_i[j]   = ncols;
1423:         req_size[i] += ncols;
1424:       }
1425:       req_source1[i] = r_status1[i].MPI_SOURCE;

1427:       /* form the header */
1428:       sbuf2_i[0] = req_size[i];
1429:       for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];

1431:       MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);
1432:     }

1434:     PetscFree(r_status1);
1435:     PetscFree(r_waits1);

1437:     /* rbuf2 is received, Post recv column indices a->j */
1438:     MPI_Waitall(nrqs,r_waits2,r_status2);

1440:     PetscMalloc4(nrqs+1,&r_waits3,nrqr+1,&s_waits3,nrqs+1,&r_status3,nrqr+1,&s_status3);
1441:     for (i=0; i<nrqs; ++i) {
1442:       PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);
1443:       req_source2[i] = r_status2[i].MPI_SOURCE;
1444:       MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);
1445:     }

1447:     /* Wait on sends1 and sends2 */
1448:     PetscMalloc1(nrqs+1,&s_status1);
1449:     MPI_Waitall(nrqs,s_waits1,s_status1);
1450:     PetscFree(s_waits1);
1451:     PetscFree(s_status1);

1453:     MPI_Waitall(nrqr,s_waits2,s_status2);
1454:     PetscFree4(r_status2,s_waits2,r_waits2,s_status2);

1456:     /* Now allocate sending buffers for a->j, and send them off */
1457:     PetscMalloc1(nrqr+1,&sbuf_aj);
1458:     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1459:     PetscMalloc1(j+1,&sbuf_aj[0]);
1460:     for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];

1462:     for (i=0; i<nrqr; i++) { /* for each requested message */
1463:       rbuf1_i   = rbuf1[i];
1464:       sbuf_aj_i = sbuf_aj[i];
1465:       ct1       = 2*rbuf1_i[0] + 1;
1466:       ct2       = 0;
1467:       /* max1=rbuf1_i[0]; if (max1 != 1) SETERRQ1(PETSC_COMM_SELF,0,"max1 %d != 1",max1); */

1469:       kmax = rbuf1[i][2];
1470:       for (k=0; k<kmax; k++,ct1++) { /* for each row */
1471:         row    = rbuf1_i[ct1] - rstart;
1472:         nzA    = ai[row+1] - ai[row];
1473:         nzB    = bi[row+1] - bi[row];
1474:         ncols  = nzA + nzB;
1475:         cworkA = aj + ai[row]; cworkB = bj + bi[row];

1477:         /* load the column indices for this row into cols*/
1478:         cols = sbuf_aj_i + ct2;

1480:         lwrite = 0;
1481:         for (l=0; l<nzB; l++) {
1482:           if ((ctmp = bmap[cworkB[l]]) < cstart) cols[lwrite++] = ctmp;
1483:         }
1484:         for (l=0; l<nzA; l++) cols[lwrite++] = cstart + cworkA[l];
1485:         for (l=0; l<nzB; l++) {
1486:           if ((ctmp = bmap[cworkB[l]]) >= cend) cols[lwrite++] = ctmp;
1487:         }

1489:         ct2 += ncols;
1490:       }
1491:       MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);
1492:     }

1494:     /* create column map (cmap): global col of C -> local col of submat */
1495: #if defined(PETSC_USE_CTABLE)
1496:     if (!allcolumns) {
1497:       PetscTableCreate(ncol+1,C->cmap->N+1,&cmap);
1498:       PetscCalloc1(C->cmap->n,&cmap_loc);
1499:       for (j=0; j<ncol; j++) { /* use array cmap_loc[] for local col indices */
1500:         if (icol[j] >= cstart && icol[j] <cend) {
1501:           cmap_loc[icol[j] - cstart] = j+1;
1502:         } else { /* use PetscTable for non-local col indices */
1503:           PetscTableAdd(cmap,icol[j]+1,j+1,INSERT_VALUES);
1504:         }
1505:       }
1506:     } else {
1507:       cmap     = NULL;
1508:       cmap_loc = NULL;
1509:     }
1510:     PetscCalloc1(C->rmap->n,&rmap_loc);
1511: #else
1512:     if (!allcolumns) {
1513:       PetscCalloc1(C->cmap->N,&cmap);
1514:       for (j=0; j<ncol; j++) cmap[icol[j]] = j+1;
1515:     } else {
1516:       cmap = NULL;
1517:     }
1518: #endif

1520:     /* Create lens for MatSeqAIJSetPreallocation() */
1521:     PetscCalloc1(nrow,&lens);

1523:     /* Compute lens from local part of C */
1524:     for (j=0; j<nrow; j++) {
1525:       row  = irow[j];
1526:       proc = row2proc[j];
1527:       if (proc == rank) {
1528:         /* diagonal part A = c->A */
1529:         ncols = ai[row-rstart+1] - ai[row-rstart];
1530:         cols  = aj + ai[row-rstart];
1531:         if (!allcolumns) {
1532:           for (k=0; k<ncols; k++) {
1533: #if defined(PETSC_USE_CTABLE)
1534:             tcol = cmap_loc[cols[k]];
1535: #else
1536:             tcol = cmap[cols[k]+cstart];
1537: #endif
1538:             if (tcol) lens[j]++;
1539:           }
1540:         } else { /* allcolumns */
1541:           lens[j] = ncols;
1542:         }

1544:         /* off-diagonal part B = c->B */
1545:         ncols = bi[row-rstart+1] - bi[row-rstart];
1546:         cols  = bj + bi[row-rstart];
1547:         if (!allcolumns) {
1548:           for (k=0; k<ncols; k++) {
1549: #if defined(PETSC_USE_CTABLE)
1550:             PetscTableFind(cmap,bmap[cols[k]]+1,&tcol);
1551: #else
1552:             tcol = cmap[bmap[cols[k]]];
1553: #endif
1554:             if (tcol) lens[j]++;
1555:           }
1556:         } else { /* allcolumns */
1557:           lens[j] += ncols;
1558:         }
1559:       }
1560:     }

1562:     /* Create row map (rmap): global row of C -> local row of submat */
1563: #if defined(PETSC_USE_CTABLE)
1564:     PetscTableCreate(nrow+1,C->rmap->N+1,&rmap);
1565:     for (j=0; j<nrow; j++) {
1566:       row  = irow[j];
1567:       proc = row2proc[j];
1568:       if (proc == rank) { /* a local row */
1569:         rmap_loc[row - rstart] = j;
1570:       } else {
1571:         PetscTableAdd(rmap,irow[j]+1,j+1,INSERT_VALUES);
1572:       }
1573:     }
1574: #else
1575:     PetscCalloc1(C->rmap->N,&rmap);
1576:     for (j=0; j<nrow; j++) {
1577:       rmap[irow[j]] = j;
1578:     }
1579: #endif

1581:     /* Update lens from offproc data */
1582:     /* recv a->j is done */
1583:     MPI_Waitall(nrqs,r_waits3,r_status3);
1584:     for (i=0; i<nrqs; i++) {
1585:       proc    = pa[i];
1586:       sbuf1_i = sbuf1[proc];
1587:       /* jmax    = sbuf1_i[0]; if (jmax != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"jmax !=1"); */
1588:       ct1     = 2 + 1;
1589:       ct2     = 0;
1590:       rbuf2_i = rbuf2[i]; /* received length of C->j */
1591:       rbuf3_i = rbuf3[i]; /* received C->j */

1593:       /* is_no  = sbuf1_i[2*j-1]; if (is_no != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_no !=0"); */
1594:       max1   = sbuf1_i[2];
1595:       for (k=0; k<max1; k++,ct1++) {
1596: #if defined(PETSC_USE_CTABLE)
1597:         PetscTableFind(rmap,sbuf1_i[ct1]+1,&row);
1598:         row--;
1599:         if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
1600: #else
1601:         row = rmap[sbuf1_i[ct1]]; /* the row index in submat */
1602: #endif
1603:         /* Now, store row index of submat in sbuf1_i[ct1] */
1604:         sbuf1_i[ct1] = row;

1606:         nnz = rbuf2_i[ct1];
1607:         if (!allcolumns) {
1608:           for (l=0; l<nnz; l++,ct2++) {
1609: #if defined(PETSC_USE_CTABLE)
1610:             if (rbuf3_i[ct2] >= cstart && rbuf3_i[ct2] <cend) {
1611:               tcol = cmap_loc[rbuf3_i[ct2] - cstart];
1612:             } else {
1613:               PetscTableFind(cmap,rbuf3_i[ct2]+1,&tcol);
1614:             }
1615: #else
1616:             tcol = cmap[rbuf3_i[ct2]]; /* column index in submat */
1617: #endif
1618:             if (tcol) lens[row]++;
1619:           }
1620:         } else { /* allcolumns */
1621:           lens[row] += nnz;
1622:         }
1623:       }
1624:     }
1625:     MPI_Waitall(nrqr,s_waits3,s_status3);
1626:     PetscFree4(r_waits3,s_waits3,r_status3,s_status3);

1628:     /* Create the submatrices */
1629:     MatCreate(PETSC_COMM_SELF,&submat);
1630:     MatSetSizes(submat,nrow,ncol,PETSC_DETERMINE,PETSC_DETERMINE);

1632:     ISGetBlockSize(isrow[0],&i);
1633:     ISGetBlockSize(iscol[0],&j);
1634:     MatSetBlockSizes(submat,i,j);
1635:     MatSetType(submat,((PetscObject)A)->type_name);
1636:     MatSeqAIJSetPreallocation(submat,0,lens);

1638:     /* create struct Mat_SubSppt and attached it to submat */
1639:     PetscNew(&smatis1);
1640:     subc = (Mat_SeqAIJ*)submat->data;
1641:     subc->submatis1 = smatis1;

1643:     smatis1->id          = 0;
1644:     smatis1->nrqs        = nrqs;
1645:     smatis1->nrqr        = nrqr;
1646:     smatis1->rbuf1       = rbuf1;
1647:     smatis1->rbuf2       = rbuf2;
1648:     smatis1->rbuf3       = rbuf3;
1649:     smatis1->sbuf2       = sbuf2;
1650:     smatis1->req_source2 = req_source2;

1652:     smatis1->sbuf1       = sbuf1;
1653:     smatis1->ptr         = ptr;
1654:     smatis1->tmp         = tmp;
1655:     smatis1->ctr         = ctr;

1657:     smatis1->pa           = pa;
1658:     smatis1->req_size     = req_size;
1659:     smatis1->req_source1  = req_source1;

1661:     smatis1->allcolumns  = allcolumns;
1662:     smatis1->singleis    = PETSC_TRUE;
1663:     smatis1->row2proc    = row2proc;
1664:     smatis1->rmap        = rmap;
1665:     smatis1->cmap        = cmap;
1666: #if defined(PETSC_USE_CTABLE)
1667:     smatis1->rmap_loc    = rmap_loc;
1668:     smatis1->cmap_loc    = cmap_loc;
1669: #endif

1671:     smatis1->destroy     = submat->ops->destroy;
1672:     submat->ops->destroy = MatDestroySubMatrix_SeqAIJ;
1673:     submat->factortype   = C->factortype;

1675:     /* compute rmax */
1676:     rmax = 0;
1677:     for (i=0; i<nrow; i++) rmax = PetscMax(rmax,lens[i]);

1679:   } else { /* scall == MAT_REUSE_MATRIX */
1680:     submat = submats[0];
1681:     if (submat->rmap->n != nrow || submat->cmap->n != ncol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");

1683:     subc    = (Mat_SeqAIJ*)submat->data;
1684:     rmax    = subc->rmax;
1685:     smatis1 = subc->submatis1;
1686:     nrqs        = smatis1->nrqs;
1687:     nrqr        = smatis1->nrqr;
1688:     rbuf1       = smatis1->rbuf1;
1689:     rbuf2       = smatis1->rbuf2;
1690:     rbuf3       = smatis1->rbuf3;
1691:     req_source2 = smatis1->req_source2;

1693:     sbuf1     = smatis1->sbuf1;
1694:     sbuf2     = smatis1->sbuf2;
1695:     ptr       = smatis1->ptr;
1696:     tmp       = smatis1->tmp;
1697:     ctr       = smatis1->ctr;

1699:     pa         = smatis1->pa;
1700:     req_size   = smatis1->req_size;
1701:     req_source1 = smatis1->req_source1;

1703:     allcolumns = smatis1->allcolumns;
1704:     row2proc   = smatis1->row2proc;
1705:     rmap       = smatis1->rmap;
1706:     cmap       = smatis1->cmap;
1707: #if defined(PETSC_USE_CTABLE)
1708:     rmap_loc   = smatis1->rmap_loc;
1709:     cmap_loc   = smatis1->cmap_loc;
1710: #endif
1711:   }

1713:   /* Post recv matrix values */
1714:   PetscMalloc3(nrqs+1,&rbuf4, rmax,&subcols, rmax,&subvals);
1715:   PetscMalloc4(nrqs+1,&r_waits4,nrqr+1,&s_waits4,nrqs+1,&r_status4,nrqr+1,&s_status4);
1716:   PetscObjectGetNewTag((PetscObject)C,&tag4);
1717:   for (i=0; i<nrqs; ++i) {
1718:     PetscMalloc1(rbuf2[i][0]+1,&rbuf4[i]);
1719:     MPI_Irecv(rbuf4[i],rbuf2[i][0],MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);
1720:   }

1722:   /* Allocate sending buffers for a->a, and send them off */
1723:   PetscMalloc1(nrqr+1,&sbuf_aa);
1724:   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
1725:   PetscMalloc1(j+1,&sbuf_aa[0]);
1726:   for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1];

1728:   for (i=0; i<nrqr; i++) {
1729:     rbuf1_i   = rbuf1[i];
1730:     sbuf_aa_i = sbuf_aa[i];
1731:     ct1       = 2*rbuf1_i[0]+1;
1732:     ct2       = 0;
1733:     /* max1=rbuf1_i[0]; if (max1 != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"max1 !=1"); */

1735:     kmax = rbuf1_i[2];
1736:     for (k=0; k<kmax; k++,ct1++) {
1737:       row = rbuf1_i[ct1] - rstart;
1738:       nzA = ai[row+1] - ai[row];
1739:       nzB = bi[row+1] - bi[row];
1740:       ncols  = nzA + nzB;
1741:       cworkB = bj + bi[row];
1742:       vworkA = a_a + ai[row];
1743:       vworkB = b_a + bi[row];

1745:       /* load the column values for this row into vals*/
1746:       vals = sbuf_aa_i + ct2;

1748:       lwrite = 0;
1749:       for (l=0; l<nzB; l++) {
1750:         if ((bmap[cworkB[l]]) < cstart) vals[lwrite++] = vworkB[l];
1751:       }
1752:       for (l=0; l<nzA; l++) vals[lwrite++] = vworkA[l];
1753:       for (l=0; l<nzB; l++) {
1754:         if ((bmap[cworkB[l]]) >= cend) vals[lwrite++] = vworkB[l];
1755:       }

1757:       ct2 += ncols;
1758:     }
1759:     MPI_Isend(sbuf_aa_i,req_size[i],MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);
1760:   }

1762:   /* Assemble submat */
1763:   /* First assemble the local rows */
1764:   for (j=0; j<nrow; j++) {
1765:     row  = irow[j];
1766:     proc = row2proc[j];
1767:     if (proc == rank) {
1768:       Crow = row - rstart;  /* local row index of C */
1769: #if defined(PETSC_USE_CTABLE)
1770:       row = rmap_loc[Crow]; /* row index of submat */
1771: #else
1772:       row = rmap[row];
1773: #endif

1775:       if (allcolumns) {
1776:         /* diagonal part A = c->A */
1777:         ncols = ai[Crow+1] - ai[Crow];
1778:         cols  = aj + ai[Crow];
1779:         vals  = a->a + ai[Crow];
1780:         i     = 0;
1781:         for (k=0; k<ncols; k++) {
1782:           subcols[i]   = cols[k] + cstart;
1783:           subvals[i++] = vals[k];
1784:         }

1786:         /* off-diagonal part B = c->B */
1787:         ncols = bi[Crow+1] - bi[Crow];
1788:         cols  = bj + bi[Crow];
1789:         vals  = b->a + bi[Crow];
1790:         for (k=0; k<ncols; k++) {
1791:           subcols[i]   = bmap[cols[k]];
1792:           subvals[i++] = vals[k];
1793:         }

1795:         MatSetValues_SeqAIJ(submat,1,&row,i,subcols,subvals,INSERT_VALUES);

1797:       } else { /* !allcolumns */
1798: #if defined(PETSC_USE_CTABLE)
1799:         /* diagonal part A = c->A */
1800:         ncols = ai[Crow+1] - ai[Crow];
1801:         cols  = aj + ai[Crow];
1802:         vals  = a->a + ai[Crow];
1803:         i     = 0;
1804:         for (k=0; k<ncols; k++) {
1805:           tcol = cmap_loc[cols[k]];
1806:           if (tcol) {
1807:             subcols[i]   = --tcol;
1808:             subvals[i++] = vals[k];
1809:           }
1810:         }

1812:         /* off-diagonal part B = c->B */
1813:         ncols = bi[Crow+1] - bi[Crow];
1814:         cols  = bj + bi[Crow];
1815:         vals  = b->a + bi[Crow];
1816:         for (k=0; k<ncols; k++) {
1817:           PetscTableFind(cmap,bmap[cols[k]]+1,&tcol);
1818:           if (tcol) {
1819:             subcols[i]   = --tcol;
1820:             subvals[i++] = vals[k];
1821:           }
1822:         }
1823: #else
1824:         /* diagonal part A = c->A */
1825:         ncols = ai[Crow+1] - ai[Crow];
1826:         cols  = aj + ai[Crow];
1827:         vals  = a->a + ai[Crow];
1828:         i     = 0;
1829:         for (k=0; k<ncols; k++) {
1830:           tcol = cmap[cols[k]+cstart];
1831:           if (tcol) {
1832:             subcols[i]   = --tcol;
1833:             subvals[i++] = vals[k];
1834:           }
1835:         }

1837:         /* off-diagonal part B = c->B */
1838:         ncols = bi[Crow+1] - bi[Crow];
1839:         cols  = bj + bi[Crow];
1840:         vals  = b->a + bi[Crow];
1841:         for (k=0; k<ncols; k++) {
1842:           tcol = cmap[bmap[cols[k]]];
1843:           if (tcol) {
1844:             subcols[i]   = --tcol;
1845:             subvals[i++] = vals[k];
1846:           }
1847:         }
1848: #endif
1849:         MatSetValues_SeqAIJ(submat,1,&row,i,subcols,subvals,INSERT_VALUES);
1850:       }
1851:     }
1852:   }

1854:   /* Now assemble the off-proc rows */
1855:   for (i=0; i<nrqs; i++) { /* for each requested message */
1856:     /* recv values from other processes */
1857:     MPI_Waitany(nrqs,r_waits4,&idex,r_status4+i);
1858:     proc    = pa[idex];
1859:     sbuf1_i = sbuf1[proc];
1860:     /* jmax    = sbuf1_i[0]; if (jmax != 1)SETERRQ1(PETSC_COMM_SELF,0,"jmax %d != 1",jmax); */
1861:     ct1     = 2 + 1;
1862:     ct2     = 0; /* count of received C->j */
1863:     ct3     = 0; /* count of received C->j that will be inserted into submat */
1864:     rbuf2_i = rbuf2[idex]; /* int** received length of C->j from other processes */
1865:     rbuf3_i = rbuf3[idex]; /* int** received C->j from other processes */
1866:     rbuf4_i = rbuf4[idex]; /* scalar** received C->a from other processes */

1868:     /* is_no = sbuf1_i[2*j-1]; if (is_no != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_no !=0"); */
1869:     max1 = sbuf1_i[2];             /* num of rows */
1870:     for (k=0; k<max1; k++,ct1++) { /* for each recved row */
1871:       row = sbuf1_i[ct1]; /* row index of submat */
1872:       if (!allcolumns) {
1873:         idex = 0;
1874:         if (scall == MAT_INITIAL_MATRIX || !iscolsorted) {
1875:           nnz  = rbuf2_i[ct1]; /* num of C entries in this row */
1876:           for (l=0; l<nnz; l++,ct2++) { /* for each recved column */
1877: #if defined(PETSC_USE_CTABLE)
1878:             if (rbuf3_i[ct2] >= cstart && rbuf3_i[ct2] <cend) {
1879:               tcol = cmap_loc[rbuf3_i[ct2] - cstart];
1880:             } else {
1881:               PetscTableFind(cmap,rbuf3_i[ct2]+1,&tcol);
1882:             }
1883: #else
1884:             tcol = cmap[rbuf3_i[ct2]];
1885: #endif
1886:             if (tcol) {
1887:               subcols[idex]   = --tcol; /* may not be sorted */
1888:               subvals[idex++] = rbuf4_i[ct2];

1890:               /* We receive an entire column of C, but a subset of it needs to be inserted into submat.
1891:                For reuse, we replace received C->j with index that should be inserted to submat */
1892:               if (iscolsorted) rbuf3_i[ct3++] = ct2;
1893:             }
1894:           }
1895:           MatSetValues_SeqAIJ(submat,1,&row,idex,subcols,subvals,INSERT_VALUES);
1896:         } else { /* scall == MAT_REUSE_MATRIX */
1897:           submat = submats[0];
1898:           subc   = (Mat_SeqAIJ*)submat->data;

1900:           nnz = subc->i[row+1] - subc->i[row]; /* num of submat entries in this row */
1901:           for (l=0; l<nnz; l++) {
1902:             ct2 = rbuf3_i[ct3++]; /* index of rbuf4_i[] which needs to be inserted into submat */
1903:             subvals[idex++] = rbuf4_i[ct2];
1904:           }

1906:           bj = subc->j + subc->i[row]; /* sorted column indices */
1907:           MatSetValues_SeqAIJ(submat,1,&row,nnz,bj,subvals,INSERT_VALUES);
1908:         }
1909:       } else { /* allcolumns */
1910:         nnz  = rbuf2_i[ct1]; /* num of C entries in this row */
1911:         MatSetValues_SeqAIJ(submat,1,&row,nnz,rbuf3_i+ct2,rbuf4_i+ct2,INSERT_VALUES);
1912:         ct2 += nnz;
1913:       }
1914:     }
1915:   }

1917:   /* sending a->a are done */
1918:   MPI_Waitall(nrqr,s_waits4,s_status4);
1919:   PetscFree4(r_waits4,s_waits4,r_status4,s_status4);

1921:   MatAssemblyBegin(submat,MAT_FINAL_ASSEMBLY);
1922:   MatAssemblyEnd(submat,MAT_FINAL_ASSEMBLY);
1923:   submats[0] = submat;

1925:   /* Restore the indices */
1926:   ISRestoreIndices(isrow[0],&irow);
1927:   if (!allcolumns) {
1928:     ISRestoreIndices(iscol[0],&icol);
1929:   }

1931:   /* Destroy allocated memory */
1932:   for (i=0; i<nrqs; ++i) {
1933:     PetscFree3(rbuf4[i],subcols,subvals);
1934:   }
1935:   PetscFree3(rbuf4,subcols,subvals);
1936:   PetscFree(sbuf_aa[0]);
1937:   PetscFree(sbuf_aa);

1939:   if (scall == MAT_INITIAL_MATRIX) {
1940:     PetscFree(lens);
1941:     PetscFree(sbuf_aj[0]);
1942:     PetscFree(sbuf_aj);
1943:   }
1944:   return(0);
1945: }

1947: PetscErrorCode MatCreateSubMatrices_MPIAIJ_SingleIS(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
1948: {
1950:   PetscInt       ncol;
1951:   PetscBool      colflag,allcolumns=PETSC_FALSE;

1954:   /* Allocate memory to hold all the submatrices */
1955:   if (scall == MAT_INITIAL_MATRIX) {
1956:     PetscCalloc1(2,submat);
1957:   }

1959:   /* Check for special case: each processor gets entire matrix columns */
1960:   ISIdentity(iscol[0],&colflag);
1961:   ISGetLocalSize(iscol[0],&ncol);
1962:   if (colflag && ncol == C->cmap->N) allcolumns = PETSC_TRUE;

1964:   MatCreateSubMatrices_MPIAIJ_SingleIS_Local(C,ismax,isrow,iscol,scall,allcolumns,*submat);
1965:   return(0);
1966: }

1968: PetscErrorCode MatCreateSubMatrices_MPIAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
1969: {
1971:   PetscInt       nmax,nstages=0,i,pos,max_no,nrow,ncol,in[2],out[2];
1972:   PetscBool      rowflag,colflag,wantallmatrix=PETSC_FALSE;
1973:   Mat_SeqAIJ     *subc;
1974:   Mat_SubSppt    *smat;

1977:   /* Check for special case: each processor has a single IS */
1978:   if (C->submat_singleis) { /* flag is set in PCSetUp_ASM() to skip MPIU_Allreduce() */
1979:     MatCreateSubMatrices_MPIAIJ_SingleIS(C,ismax,isrow,iscol,scall,submat);
1980:     C->submat_singleis = PETSC_FALSE; /* resume its default value in case C will be used for non-single IS */
1981:     return(0);
1982:   }

1984:   /* Collect global wantallmatrix and nstages */
1985:   if (!C->cmap->N) nmax=20*1000000/sizeof(PetscInt);
1986:   else nmax = 20*1000000 / (C->cmap->N * sizeof(PetscInt));
1987:   if (!nmax) nmax = 1;

1989:   if (scall == MAT_INITIAL_MATRIX) {
1990:     /* Collect global wantallmatrix and nstages */
1991:     if (ismax == 1 && C->rmap->N == C->cmap->N) {
1992:       ISIdentity(*isrow,&rowflag);
1993:       ISIdentity(*iscol,&colflag);
1994:       ISGetLocalSize(*isrow,&nrow);
1995:       ISGetLocalSize(*iscol,&ncol);
1996:       if (rowflag && colflag && nrow == C->rmap->N && ncol == C->cmap->N) {
1997:         wantallmatrix = PETSC_TRUE;

1999:         PetscOptionsGetBool(((PetscObject)C)->options,((PetscObject)C)->prefix,"-use_fast_submatrix",&wantallmatrix,NULL);
2000:       }
2001:     }

2003:     /* Determine the number of stages through which submatrices are done
2004:        Each stage will extract nmax submatrices.
2005:        nmax is determined by the matrix column dimension.
2006:        If the original matrix has 20M columns, only one submatrix per stage is allowed, etc.
2007:     */
2008:     nstages = ismax/nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */

2010:     in[0] = -1*(PetscInt)wantallmatrix;
2011:     in[1] = nstages;
2012:     MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));
2013:     wantallmatrix = (PetscBool)(-out[0]);
2014:     nstages       = out[1]; /* Make sure every processor loops through the global nstages */

2016:   } else { /* MAT_REUSE_MATRIX */
2017:     if (ismax) {
2018:       subc = (Mat_SeqAIJ*)(*submat)[0]->data;
2019:       smat = subc->submatis1;
2020:     } else { /* (*submat)[0] is a dummy matrix */
2021:       smat = (Mat_SubSppt*)(*submat)[0]->data;
2022:     }
2023:     if (!smat) {
2024:       /* smat is not generated by MatCreateSubMatrix_MPIAIJ_All(...,MAT_INITIAL_MATRIX,...) */
2025:       wantallmatrix = PETSC_TRUE;
2026:     } else if (smat->singleis) {
2027:       MatCreateSubMatrices_MPIAIJ_SingleIS(C,ismax,isrow,iscol,scall,submat);
2028:       return(0);
2029:     } else {
2030:       nstages = smat->nstages;
2031:     }
2032:   }

2034:   if (wantallmatrix) {
2035:     MatCreateSubMatrix_MPIAIJ_All(C,MAT_GET_VALUES,scall,submat);
2036:     return(0);
2037:   }

2039:   /* Allocate memory to hold all the submatrices and dummy submatrices */
2040:   if (scall == MAT_INITIAL_MATRIX) {
2041:     PetscCalloc1(ismax+nstages,submat);
2042:   }

2044:   for (i=0,pos=0; i<nstages; i++) {
2045:     if (pos+nmax <= ismax) max_no = nmax;
2046:     else if (pos >= ismax) max_no = 0;
2047:     else                   max_no = ismax-pos;

2049:     MatCreateSubMatrices_MPIAIJ_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);
2050:     if (!max_no && scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */
2051:       smat = (Mat_SubSppt*)(*submat)[pos]->data; pos++;
2052:       smat->nstages = nstages;
2053:     }
2054:     pos += max_no;
2055:   }

2057:   if (ismax && scall == MAT_INITIAL_MATRIX) {
2058:     /* save nstages for reuse */
2059:     subc = (Mat_SeqAIJ*)(*submat)[0]->data;
2060:     smat = subc->submatis1;
2061:     smat->nstages = nstages;
2062:   }
2063:   return(0);
2064: }

2066: /* -------------------------------------------------------------------------*/
2067: PetscErrorCode MatCreateSubMatrices_MPIAIJ_Local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
2068: {
2069:   Mat_MPIAIJ     *c = (Mat_MPIAIJ*)C->data;
2070:   Mat            A  = c->A;
2071:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)c->B->data,*subc;
2072:   const PetscInt **icol,**irow;
2073:   PetscInt       *nrow,*ncol,start;
2075:   PetscMPIInt    rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr;
2076:   PetscInt       **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1;
2077:   PetscInt       nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol;
2078:   PetscInt       **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2;
2079:   PetscInt       **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax;
2080: #if defined(PETSC_USE_CTABLE)
2081:   PetscTable     *cmap,cmap_i=NULL,*rmap,rmap_i;
2082: #else
2083:   PetscInt       **cmap,*cmap_i=NULL,**rmap,*rmap_i;
2084: #endif
2085:   const PetscInt *irow_i;
2086:   PetscInt       ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i;
2087:   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
2088:   MPI_Request    *r_waits4,*s_waits3,*s_waits4;
2089:   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
2090:   MPI_Status     *r_status3,*r_status4,*s_status4;
2091:   MPI_Comm       comm;
2092:   PetscScalar    **rbuf4,*rbuf4_i,**sbuf_aa,*vals,*mat_a,*imat_a,*sbuf_aa_i;
2093:   PetscMPIInt    *onodes1,*olengths1,end;
2094:   PetscInt       **row2proc,*row2proc_i,ilen_row,*imat_ilen,*imat_j,*imat_i,old_row;
2095:   Mat_SubSppt    *smat_i;
2096:   PetscBool      *issorted,*allcolumns,colflag,iscsorted=PETSC_TRUE;
2097:   PetscInt       *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen;

2100:   PetscObjectGetComm((PetscObject)C,&comm);
2101:   size = c->size;
2102:   rank = c->rank;

2104:   PetscMalloc4(ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax+1,&allcolumns);
2105:   PetscMalloc5(ismax,(PetscInt***)&irow,ismax,(PetscInt***)&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);

2107:   for (i=0; i<ismax; i++) {
2108:     ISSorted(iscol[i],&issorted[i]);
2109:     if (!issorted[i]) iscsorted = issorted[i];

2111:     ISSorted(isrow[i],&issorted[i]);

2113:     ISGetIndices(isrow[i],&irow[i]);
2114:     ISGetLocalSize(isrow[i],&nrow[i]);

2116:     /* Check for special case: allcolumn */
2117:     ISIdentity(iscol[i],&colflag);
2118:     ISGetLocalSize(iscol[i],&ncol[i]);
2119:     if (colflag && ncol[i] == C->cmap->N) {
2120:       allcolumns[i] = PETSC_TRUE;
2121:       icol[i] = NULL;
2122:     } else {
2123:       allcolumns[i] = PETSC_FALSE;
2124:       ISGetIndices(iscol[i],&icol[i]);
2125:     }
2126:   }

2128:   if (scall == MAT_REUSE_MATRIX) {
2129:     /* Assumes new rows are same length as the old rows */
2130:     for (i=0; i<ismax; i++) {
2131:       if (!submats[i]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats[%D] is null, cannot reuse",i);
2132:       subc = (Mat_SeqAIJ*)submats[i]->data;
2133:       if ((submats[i]->rmap->n != nrow[i]) || (submats[i]->cmap->n != ncol[i])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");

2135:       /* Initial matrix as if empty */
2136:       PetscArrayzero(subc->ilen,submats[i]->rmap->n);

2138:       smat_i   = subc->submatis1;

2140:       nrqs        = smat_i->nrqs;
2141:       nrqr        = smat_i->nrqr;
2142:       rbuf1       = smat_i->rbuf1;
2143:       rbuf2       = smat_i->rbuf2;
2144:       rbuf3       = smat_i->rbuf3;
2145:       req_source2 = smat_i->req_source2;

2147:       sbuf1     = smat_i->sbuf1;
2148:       sbuf2     = smat_i->sbuf2;
2149:       ptr       = smat_i->ptr;
2150:       tmp       = smat_i->tmp;
2151:       ctr       = smat_i->ctr;

2153:       pa          = smat_i->pa;
2154:       req_size    = smat_i->req_size;
2155:       req_source1 = smat_i->req_source1;

2157:       allcolumns[i] = smat_i->allcolumns;
2158:       row2proc[i]   = smat_i->row2proc;
2159:       rmap[i]       = smat_i->rmap;
2160:       cmap[i]       = smat_i->cmap;
2161:     }

2163:     if (!ismax){ /* Get dummy submatrices and retrieve struct submatis1 */
2164:       if (!submats[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"submats are null, cannot reuse");
2165:       smat_i = (Mat_SubSppt*)submats[0]->data;

2167:       nrqs        = smat_i->nrqs;
2168:       nrqr        = smat_i->nrqr;
2169:       rbuf1       = smat_i->rbuf1;
2170:       rbuf2       = smat_i->rbuf2;
2171:       rbuf3       = smat_i->rbuf3;
2172:       req_source2 = smat_i->req_source2;

2174:       sbuf1       = smat_i->sbuf1;
2175:       sbuf2       = smat_i->sbuf2;
2176:       ptr         = smat_i->ptr;
2177:       tmp         = smat_i->tmp;
2178:       ctr         = smat_i->ctr;

2180:       pa          = smat_i->pa;
2181:       req_size    = smat_i->req_size;
2182:       req_source1 = smat_i->req_source1;

2184:       allcolumns[0] = PETSC_FALSE;
2185:     }
2186:   } else { /* scall == MAT_INITIAL_MATRIX */
2187:     /* Get some new tags to keep the communication clean */
2188:     PetscObjectGetNewTag((PetscObject)C,&tag2);
2189:     PetscObjectGetNewTag((PetscObject)C,&tag3);

2191:     /* evaluate communication - mesg to who, length of mesg, and buffer space
2192:      required. Based on this, buffers are allocated, and data copied into them*/
2193:     PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);   /* mesg size, initialize work vectors */

2195:     for (i=0; i<ismax; i++) {
2196:       jmax   = nrow[i];
2197:       irow_i = irow[i];

2199:       PetscMalloc1(jmax,&row2proc_i);
2200:       row2proc[i] = row2proc_i;

2202:       if (issorted[i]) proc = 0;
2203:       for (j=0; j<jmax; j++) {
2204:         if (!issorted[i]) proc = 0;
2205:         row = irow_i[j];
2206:         while (row >= C->rmap->range[proc+1]) proc++;
2207:         w4[proc]++;
2208:         row2proc_i[j] = proc; /* map row index to proc */
2209:       }
2210:       for (j=0; j<size; j++) {
2211:         if (w4[j]) { w1[j] += w4[j];  w3[j]++; w4[j] = 0;}
2212:       }
2213:     }

2215:     nrqs     = 0;              /* no of outgoing messages */
2216:     msz      = 0;              /* total mesg length (for all procs) */
2217:     w1[rank] = 0;              /* no mesg sent to self */
2218:     w3[rank] = 0;
2219:     for (i=0; i<size; i++) {
2220:       if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
2221:     }
2222:     PetscMalloc1(nrqs+1,&pa); /*(proc -array)*/
2223:     for (i=0,j=0; i<size; i++) {
2224:       if (w1[i]) { pa[j] = i; j++; }
2225:     }

2227:     /* Each message would have a header = 1 + 2*(no of IS) + data */
2228:     for (i=0; i<nrqs; i++) {
2229:       j      = pa[i];
2230:       w1[j] += w2[j] + 2* w3[j];
2231:       msz   += w1[j];
2232:     }
2233:     PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);

2235:     /* Determine the number of messages to expect, their lengths, from from-ids */
2236:     PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
2237:     PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);

2239:     /* Now post the Irecvs corresponding to these messages */
2240:     tag0 = ((PetscObject)C)->tag;
2241:     PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);

2243:     PetscFree(onodes1);
2244:     PetscFree(olengths1);

2246:     /* Allocate Memory for outgoing messages */
2247:     PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);
2248:     PetscArrayzero(sbuf1,size);
2249:     PetscArrayzero(ptr,size);

2251:     {
2252:       PetscInt *iptr = tmp;
2253:       k    = 0;
2254:       for (i=0; i<nrqs; i++) {
2255:         j        = pa[i];
2256:         iptr    += k;
2257:         sbuf1[j] = iptr;
2258:         k        = w1[j];
2259:       }
2260:     }

2262:     /* Form the outgoing messages. Initialize the header space */
2263:     for (i=0; i<nrqs; i++) {
2264:       j           = pa[i];
2265:       sbuf1[j][0] = 0;
2266:       PetscArrayzero(sbuf1[j]+1,2*w3[j]);
2267:       ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
2268:     }

2270:     /* Parse the isrow and copy data into outbuf */
2271:     for (i=0; i<ismax; i++) {
2272:       row2proc_i = row2proc[i];
2273:       PetscArrayzero(ctr,size);
2274:       irow_i = irow[i];
2275:       jmax   = nrow[i];
2276:       for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
2277:         proc = row2proc_i[j];
2278:         if (proc != rank) { /* copy to the outgoing buf*/
2279:           ctr[proc]++;
2280:           *ptr[proc] = irow_i[j];
2281:           ptr[proc]++;
2282:         }
2283:       }
2284:       /* Update the headers for the current IS */
2285:       for (j=0; j<size; j++) { /* Can Optimise this loop too */
2286:         if ((ctr_j = ctr[j])) {
2287:           sbuf1_j        = sbuf1[j];
2288:           k              = ++sbuf1_j[0];
2289:           sbuf1_j[2*k]   = ctr_j;
2290:           sbuf1_j[2*k-1] = i;
2291:         }
2292:       }
2293:     }

2295:     /*  Now  post the sends */
2296:     PetscMalloc1(nrqs+1,&s_waits1);
2297:     for (i=0; i<nrqs; ++i) {
2298:       j    = pa[i];
2299:       MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);
2300:     }

2302:     /* Post Receives to capture the buffer size */
2303:     PetscMalloc1(nrqs+1,&r_waits2);
2304:     PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);
2305:     rbuf2[0] = tmp + msz;
2306:     for (i=1; i<nrqs; ++i) {
2307:       rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
2308:     }
2309:     for (i=0; i<nrqs; ++i) {
2310:       j    = pa[i];
2311:       MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);
2312:     }

2314:     /* Send to other procs the buf size they should allocate */
2315:     /* Receive messages*/
2316:     PetscMalloc1(nrqr+1,&s_waits2);
2317:     PetscMalloc1(nrqr+1,&r_status1);
2318:     PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);
2319:     {
2320:       PetscInt   *sAi = a->i,*sBi = b->i,id,rstart = C->rmap->rstart;
2321:       PetscInt   *sbuf2_i;

2323:       MPI_Waitall(nrqr,r_waits1,r_status1);
2324:       for (i=0; i<nrqr; ++i) {
2325:         req_size[i] = 0;
2326:         rbuf1_i        = rbuf1[i];
2327:         start          = 2*rbuf1_i[0] + 1;
2328:         MPI_Get_count(r_status1+i,MPIU_INT,&end);
2329:         PetscMalloc1(end+1,&sbuf2[i]);
2330:         sbuf2_i        = sbuf2[i];
2331:         for (j=start; j<end; j++) {
2332:           id              = rbuf1_i[j] - rstart;
2333:           ncols           = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id];
2334:           sbuf2_i[j]      = ncols;
2335:           req_size[i] += ncols;
2336:         }
2337:         req_source1[i] = r_status1[i].MPI_SOURCE;
2338:         /* form the header */
2339:         sbuf2_i[0] = req_size[i];
2340:         for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j];

2342:         MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);
2343:       }
2344:     }
2345:     PetscFree(r_status1);
2346:     PetscFree(r_waits1);
2347:     PetscFree4(w1,w2,w3,w4);

2349:     /* Receive messages*/
2350:     PetscMalloc1(nrqs+1,&r_waits3);
2351:     PetscMalloc1(nrqs+1,&r_status2);

2353:     MPI_Waitall(nrqs,r_waits2,r_status2);
2354:     for (i=0; i<nrqs; ++i) {
2355:       PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);
2356:       req_source2[i] = r_status2[i].MPI_SOURCE;
2357:       MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);
2358:     }
2359:     PetscFree(r_status2);
2360:     PetscFree(r_waits2);

2362:     /* Wait on sends1 and sends2 */
2363:     PetscMalloc1(nrqs+1,&s_status1);
2364:     PetscMalloc1(nrqr+1,&s_status2);

2366:     if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
2367:     if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
2368:     PetscFree(s_status1);
2369:     PetscFree(s_status2);
2370:     PetscFree(s_waits1);
2371:     PetscFree(s_waits2);

2373:     /* Now allocate sending buffers for a->j, and send them off */
2374:     PetscMalloc1(nrqr+1,&sbuf_aj);
2375:     for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2376:     PetscMalloc1(j+1,&sbuf_aj[0]);
2377:     for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];

2379:     PetscMalloc1(nrqr+1,&s_waits3);
2380:     {
2381:       PetscInt nzA,nzB,*a_i = a->i,*b_i = b->i,lwrite;
2382:       PetscInt *cworkA,*cworkB,cstart = C->cmap->rstart,rstart = C->rmap->rstart,*bmap = c->garray;
2383:       PetscInt cend = C->cmap->rend;
2384:       PetscInt *a_j = a->j,*b_j = b->j,ctmp;

2386:       for (i=0; i<nrqr; i++) {
2387:         rbuf1_i   = rbuf1[i];
2388:         sbuf_aj_i = sbuf_aj[i];
2389:         ct1       = 2*rbuf1_i[0] + 1;
2390:         ct2       = 0;
2391:         for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2392:           kmax = rbuf1[i][2*j];
2393:           for (k=0; k<kmax; k++,ct1++) {
2394:             row    = rbuf1_i[ct1] - rstart;
2395:             nzA    = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row];
2396:             ncols  = nzA + nzB;
2397:             cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row];

2399:             /* load the column indices for this row into cols */
2400:             cols = sbuf_aj_i + ct2;

2402:             lwrite = 0;
2403:             for (l=0; l<nzB; l++) {
2404:               if ((ctmp = bmap[cworkB[l]]) < cstart) cols[lwrite++] = ctmp;
2405:             }
2406:             for (l=0; l<nzA; l++) cols[lwrite++] = cstart + cworkA[l];
2407:             for (l=0; l<nzB; l++) {
2408:               if ((ctmp = bmap[cworkB[l]]) >= cend) cols[lwrite++] = ctmp;
2409:             }

2411:             ct2 += ncols;
2412:           }
2413:         }
2414:         MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);
2415:       }
2416:     }
2417:     PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);

2419:     /* create col map: global col of C -> local col of submatrices */
2420:     {
2421:       const PetscInt *icol_i;
2422: #if defined(PETSC_USE_CTABLE)
2423:       for (i=0; i<ismax; i++) {
2424:         if (!allcolumns[i]) {
2425:           PetscTableCreate(ncol[i]+1,C->cmap->N+1,&cmap[i]);

2427:           jmax   = ncol[i];
2428:           icol_i = icol[i];
2429:           cmap_i = cmap[i];
2430:           for (j=0; j<jmax; j++) {
2431:             PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);
2432:           }
2433:         } else cmap[i] = NULL;
2434:       }
2435: #else
2436:       for (i=0; i<ismax; i++) {
2437:         if (!allcolumns[i]) {
2438:           PetscCalloc1(C->cmap->N,&cmap[i]);
2439:           jmax   = ncol[i];
2440:           icol_i = icol[i];
2441:           cmap_i = cmap[i];
2442:           for (j=0; j<jmax; j++) {
2443:             cmap_i[icol_i[j]] = j+1;
2444:           }
2445:         } else cmap[i] = NULL;
2446:       }
2447: #endif
2448:     }

2450:     /* Create lens which is required for MatCreate... */
2451:     for (i=0,j=0; i<ismax; i++) j += nrow[i];
2452:     PetscMalloc1(ismax,&lens);

2454:     if (ismax) {
2455:       PetscCalloc1(j,&lens[0]);
2456:     }
2457:     for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1];

2459:     /* Update lens from local data */
2460:     for (i=0; i<ismax; i++) {
2461:       row2proc_i = row2proc[i];
2462:       jmax = nrow[i];
2463:       if (!allcolumns[i]) cmap_i = cmap[i];
2464:       irow_i = irow[i];
2465:       lens_i = lens[i];
2466:       for (j=0; j<jmax; j++) {
2467:         row = irow_i[j];
2468:         proc = row2proc_i[j];
2469:         if (proc == rank) {
2470:           MatGetRow_MPIAIJ(C,row,&ncols,&cols,0);
2471:           if (!allcolumns[i]) {
2472:             for (k=0; k<ncols; k++) {
2473: #if defined(PETSC_USE_CTABLE)
2474:               PetscTableFind(cmap_i,cols[k]+1,&tcol);
2475: #else
2476:               tcol = cmap_i[cols[k]];
2477: #endif
2478:               if (tcol) lens_i[j]++;
2479:             }
2480:           } else { /* allcolumns */
2481:             lens_i[j] = ncols;
2482:           }
2483:           MatRestoreRow_MPIAIJ(C,row,&ncols,&cols,0);
2484:         }
2485:       }
2486:     }

2488:     /* Create row map: global row of C -> local row of submatrices */
2489: #if defined(PETSC_USE_CTABLE)
2490:     for (i=0; i<ismax; i++) {
2491:       PetscTableCreate(nrow[i]+1,C->rmap->N+1,&rmap[i]);
2492:       irow_i = irow[i];
2493:       jmax   = nrow[i];
2494:       for (j=0; j<jmax; j++) {
2495:       PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);
2496:       }
2497:     }
2498: #else
2499:     for (i=0; i<ismax; i++) {
2500:       PetscCalloc1(C->rmap->N,&rmap[i]);
2501:       rmap_i = rmap[i];
2502:       irow_i = irow[i];
2503:       jmax   = nrow[i];
2504:       for (j=0; j<jmax; j++) {
2505:         rmap_i[irow_i[j]] = j;
2506:       }
2507:     }
2508: #endif

2510:     /* Update lens from offproc data */
2511:     {
2512:       PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i;

2514:       MPI_Waitall(nrqs,r_waits3,r_status3);
2515:       for (tmp2=0; tmp2<nrqs; tmp2++) {
2516:         sbuf1_i = sbuf1[pa[tmp2]];
2517:         jmax    = sbuf1_i[0];
2518:         ct1     = 2*jmax+1;
2519:         ct2     = 0;
2520:         rbuf2_i = rbuf2[tmp2];
2521:         rbuf3_i = rbuf3[tmp2];
2522:         for (j=1; j<=jmax; j++) {
2523:           is_no  = sbuf1_i[2*j-1];
2524:           max1   = sbuf1_i[2*j];
2525:           lens_i = lens[is_no];
2526:           if (!allcolumns[is_no]) cmap_i = cmap[is_no];
2527:           rmap_i = rmap[is_no];
2528:           for (k=0; k<max1; k++,ct1++) {
2529: #if defined(PETSC_USE_CTABLE)
2530:             PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);
2531:             row--;
2532:             if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table");
2533: #else
2534:             row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
2535: #endif
2536:             max2 = rbuf2_i[ct1];
2537:             for (l=0; l<max2; l++,ct2++) {
2538:               if (!allcolumns[is_no]) {
2539: #if defined(PETSC_USE_CTABLE)
2540:                 PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);
2541: #else
2542:                 tcol = cmap_i[rbuf3_i[ct2]];
2543: #endif
2544:                 if (tcol) lens_i[row]++;
2545:               } else { /* allcolumns */
2546:                 lens_i[row]++; /* lens_i[row] += max2 ? */
2547:               }
2548:             }
2549:           }
2550:         }
2551:       }
2552:     }
2553:     PetscFree(r_waits3);
2554:     if (nrqr) {MPI_Waitall(nrqr,s_waits3,s_status3);}
2555:     PetscFree2(r_status3,s_status3);
2556:     PetscFree(s_waits3);

2558:     /* Create the submatrices */
2559:     for (i=0; i<ismax; i++) {
2560:       PetscInt    rbs,cbs;

2562:       ISGetBlockSize(isrow[i],&rbs);
2563:       ISGetBlockSize(iscol[i],&cbs);

2565:       MatCreate(PETSC_COMM_SELF,submats+i);
2566:       MatSetSizes(submats[i],nrow[i],ncol[i],PETSC_DETERMINE,PETSC_DETERMINE);

2568:       MatSetBlockSizes(submats[i],rbs,cbs);
2569:       MatSetType(submats[i],((PetscObject)A)->type_name);
2570:       MatSeqAIJSetPreallocation(submats[i],0,lens[i]);

2572:       /* create struct Mat_SubSppt and attached it to submat */
2573:       PetscNew(&smat_i);
2574:       subc = (Mat_SeqAIJ*)submats[i]->data;
2575:       subc->submatis1 = smat_i;

2577:       smat_i->destroy          = submats[i]->ops->destroy;
2578:       submats[i]->ops->destroy = MatDestroySubMatrix_SeqAIJ;
2579:       submats[i]->factortype   = C->factortype;

2581:       smat_i->id          = i;
2582:       smat_i->nrqs        = nrqs;
2583:       smat_i->nrqr        = nrqr;
2584:       smat_i->rbuf1       = rbuf1;
2585:       smat_i->rbuf2       = rbuf2;
2586:       smat_i->rbuf3       = rbuf3;
2587:       smat_i->sbuf2       = sbuf2;
2588:       smat_i->req_source2 = req_source2;

2590:       smat_i->sbuf1       = sbuf1;
2591:       smat_i->ptr         = ptr;
2592:       smat_i->tmp         = tmp;
2593:       smat_i->ctr         = ctr;

2595:       smat_i->pa           = pa;
2596:       smat_i->req_size     = req_size;
2597:       smat_i->req_source1  = req_source1;

2599:       smat_i->allcolumns  = allcolumns[i];
2600:       smat_i->singleis    = PETSC_FALSE;
2601:       smat_i->row2proc    = row2proc[i];
2602:       smat_i->rmap        = rmap[i];
2603:       smat_i->cmap        = cmap[i];
2604:     }

2606:     if (!ismax) { /* Create dummy submats[0] for reuse struct subc */
2607:       MatCreate(PETSC_COMM_SELF,&submats[0]);
2608:       MatSetSizes(submats[0],0,0,PETSC_DETERMINE,PETSC_DETERMINE);
2609:       MatSetType(submats[0],MATDUMMY);

2611:       /* create struct Mat_SubSppt and attached it to submat */
2612:       PetscNewLog(submats[0],&smat_i);
2613:       submats[0]->data = (void*)smat_i;

2615:       smat_i->destroy          = submats[0]->ops->destroy;
2616:       submats[0]->ops->destroy = MatDestroySubMatrix_Dummy;
2617:       submats[0]->factortype   = C->factortype;

2619:       smat_i->id          = 0;
2620:       smat_i->nrqs        = nrqs;
2621:       smat_i->nrqr        = nrqr;
2622:       smat_i->rbuf1       = rbuf1;
2623:       smat_i->rbuf2       = rbuf2;
2624:       smat_i->rbuf3       = rbuf3;
2625:       smat_i->sbuf2       = sbuf2;
2626:       smat_i->req_source2 = req_source2;

2628:       smat_i->sbuf1       = sbuf1;
2629:       smat_i->ptr         = ptr;
2630:       smat_i->tmp         = tmp;
2631:       smat_i->ctr         = ctr;

2633:       smat_i->pa           = pa;
2634:       smat_i->req_size     = req_size;
2635:       smat_i->req_source1  = req_source1;

2637:       smat_i->allcolumns  = PETSC_FALSE;
2638:       smat_i->singleis    = PETSC_FALSE;
2639:       smat_i->row2proc    = NULL;
2640:       smat_i->rmap        = NULL;
2641:       smat_i->cmap        = NULL;
2642:     }

2644:     if (ismax) {PetscFree(lens[0]);}
2645:     PetscFree(lens);
2646:     PetscFree(sbuf_aj[0]);
2647:     PetscFree(sbuf_aj);

2649:   } /* endof scall == MAT_INITIAL_MATRIX */

2651:   /* Post recv matrix values */
2652:   PetscObjectGetNewTag((PetscObject)C,&tag4);
2653:   PetscMalloc1(nrqs+1,&rbuf4);
2654:   PetscMalloc1(nrqs+1,&r_waits4);
2655:   PetscMalloc1(nrqs+1,&r_status4);
2656:   PetscMalloc1(nrqr+1,&s_status4);
2657:   for (i=0; i<nrqs; ++i) {
2658:     PetscMalloc1(rbuf2[i][0]+1,&rbuf4[i]);
2659:     MPI_Irecv(rbuf4[i],rbuf2[i][0],MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);
2660:   }

2662:   /* Allocate sending buffers for a->a, and send them off */
2663:   PetscMalloc1(nrqr+1,&sbuf_aa);
2664:   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2665:   PetscMalloc1(j+1,&sbuf_aa[0]);
2666:   for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1];

2668:   PetscMalloc1(nrqr+1,&s_waits4);
2669:   {
2670:     PetscInt    nzA,nzB,*a_i = a->i,*b_i = b->i, *cworkB,lwrite;
2671:     PetscInt    cstart = C->cmap->rstart,rstart = C->rmap->rstart,*bmap = c->garray;
2672:     PetscInt    cend   = C->cmap->rend;
2673:     PetscInt    *b_j   = b->j;
2674:     PetscScalar *vworkA,*vworkB,*a_a = a->a,*b_a = b->a;

2676:     for (i=0; i<nrqr; i++) {
2677:       rbuf1_i   = rbuf1[i];
2678:       sbuf_aa_i = sbuf_aa[i];
2679:       ct1       = 2*rbuf1_i[0]+1;
2680:       ct2       = 0;
2681:       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2682:         kmax = rbuf1_i[2*j];
2683:         for (k=0; k<kmax; k++,ct1++) {
2684:           row    = rbuf1_i[ct1] - rstart;
2685:           nzA    = a_i[row+1] - a_i[row];     nzB = b_i[row+1] - b_i[row];
2686:           ncols  = nzA + nzB;
2687:           cworkB = b_j + b_i[row];
2688:           vworkA = a_a + a_i[row];
2689:           vworkB = b_a + b_i[row];

2691:           /* load the column values for this row into vals*/
2692:           vals = sbuf_aa_i+ct2;

2694:           lwrite = 0;
2695:           for (l=0; l<nzB; l++) {
2696:             if ((bmap[cworkB[l]]) < cstart) vals[lwrite++] = vworkB[l];
2697:           }
2698:           for (l=0; l<nzA; l++) vals[lwrite++] = vworkA[l];
2699:           for (l=0; l<nzB; l++) {
2700:             if ((bmap[cworkB[l]]) >= cend) vals[lwrite++] = vworkB[l];
2701:           }

2703:           ct2 += ncols;
2704:         }
2705:       }
2706:       MPI_Isend(sbuf_aa_i,req_size[i],MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);
2707:     }
2708:   }

2710:   /* Assemble the matrices */
2711:   /* First assemble the local rows */
2712:   for (i=0; i<ismax; i++) {
2713:     row2proc_i = row2proc[i];
2714:     subc      = (Mat_SeqAIJ*)submats[i]->data;
2715:     imat_ilen = subc->ilen;
2716:     imat_j    = subc->j;
2717:     imat_i    = subc->i;
2718:     imat_a    = subc->a;

2720:     if (!allcolumns[i]) cmap_i = cmap[i];
2721:     rmap_i = rmap[i];
2722:     irow_i = irow[i];
2723:     jmax   = nrow[i];
2724:     for (j=0; j<jmax; j++) {
2725:       row  = irow_i[j];
2726:       proc = row2proc_i[j];
2727:       if (proc == rank) {
2728:         old_row = row;
2729: #if defined(PETSC_USE_CTABLE)
2730:         PetscTableFind(rmap_i,row+1,&row);
2731:         row--;
2732: #else
2733:         row = rmap_i[row];
2734: #endif
2735:         ilen_row = imat_ilen[row];
2736:         MatGetRow_MPIAIJ(C,old_row,&ncols,&cols,&vals);
2737:         mat_i    = imat_i[row];
2738:         mat_a    = imat_a + mat_i;
2739:         mat_j    = imat_j + mat_i;
2740:         if (!allcolumns[i]) {
2741:           for (k=0; k<ncols; k++) {
2742: #if defined(PETSC_USE_CTABLE)
2743:             PetscTableFind(cmap_i,cols[k]+1,&tcol);
2744: #else
2745:             tcol = cmap_i[cols[k]];
2746: #endif
2747:             if (tcol) {
2748:               *mat_j++ = tcol - 1;
2749:               *mat_a++ = vals[k];
2750:               ilen_row++;
2751:             }
2752:           }
2753:         } else { /* allcolumns */
2754:           for (k=0; k<ncols; k++) {
2755:             *mat_j++ = cols[k];  /* global col index! */
2756:             *mat_a++ = vals[k];
2757:             ilen_row++;
2758:           }
2759:         }
2760:         MatRestoreRow_MPIAIJ(C,old_row,&ncols,&cols,&vals);

2762:         imat_ilen[row] = ilen_row;
2763:       }
2764:     }
2765:   }

2767:   /* Now assemble the off proc rows */
2768:   MPI_Waitall(nrqs,r_waits4,r_status4);
2769:   for (tmp2=0; tmp2<nrqs; tmp2++) {
2770:     sbuf1_i = sbuf1[pa[tmp2]];
2771:     jmax    = sbuf1_i[0];
2772:     ct1     = 2*jmax + 1;
2773:     ct2     = 0;
2774:     rbuf2_i = rbuf2[tmp2];
2775:     rbuf3_i = rbuf3[tmp2];
2776:     rbuf4_i = rbuf4[tmp2];
2777:     for (j=1; j<=jmax; j++) {
2778:       is_no     = sbuf1_i[2*j-1];
2779:       rmap_i    = rmap[is_no];
2780:       if (!allcolumns[is_no]) cmap_i = cmap[is_no];
2781:       subc      = (Mat_SeqAIJ*)submats[is_no]->data;
2782:       imat_ilen = subc->ilen;
2783:       imat_j    = subc->j;
2784:       imat_i    = subc->i;
2785:       imat_a    = subc->a;
2786:       max1      = sbuf1_i[2*j];
2787:       for (k=0; k<max1; k++,ct1++) {
2788:         row = sbuf1_i[ct1];
2789: #if defined(PETSC_USE_CTABLE)
2790:         PetscTableFind(rmap_i,row+1,&row);
2791:         row--;
2792: #else
2793:         row = rmap_i[row];
2794: #endif
2795:         ilen  = imat_ilen[row];
2796:         mat_i = imat_i[row];
2797:         mat_a = imat_a + mat_i;
2798:         mat_j = imat_j + mat_i;
2799:         max2  = rbuf2_i[ct1];
2800:         if (!allcolumns[is_no]) {
2801:           for (l=0; l<max2; l++,ct2++) {
2802: #if defined(PETSC_USE_CTABLE)
2803:             PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);
2804: #else
2805:             tcol = cmap_i[rbuf3_i[ct2]];
2806: #endif
2807:             if (tcol) {
2808:               *mat_j++ = tcol - 1;
2809:               *mat_a++ = rbuf4_i[ct2];
2810:               ilen++;
2811:             }
2812:           }
2813:         } else { /* allcolumns */
2814:           for (l=0; l<max2; l++,ct2++) {
2815:             *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */
2816:             *mat_a++ = rbuf4_i[ct2];
2817:             ilen++;
2818:           }
2819:         }
2820:         imat_ilen[row] = ilen;
2821:       }
2822:     }
2823:   }

2825:   if (!iscsorted) { /* sort column indices of the rows */
2826:     for (i=0; i<ismax; i++) {
2827:       subc      = (Mat_SeqAIJ*)submats[i]->data;
2828:       imat_j    = subc->j;
2829:       imat_i    = subc->i;
2830:       imat_a    = subc->a;
2831:       imat_ilen = subc->ilen;

2833:       if (allcolumns[i]) continue;
2834:       jmax = nrow[i];
2835:       for (j=0; j<jmax; j++) {
2836:         mat_i = imat_i[j];
2837:         mat_a = imat_a + mat_i;
2838:         mat_j = imat_j + mat_i;
2839:         PetscSortIntWithScalarArray(imat_ilen[j],mat_j,mat_a);
2840:       }
2841:     }
2842:   }

2844:   PetscFree(r_status4);
2845:   PetscFree(r_waits4);
2846:   if (nrqr) {MPI_Waitall(nrqr,s_waits4,s_status4);}
2847:   PetscFree(s_waits4);
2848:   PetscFree(s_status4);

2850:   /* Restore the indices */
2851:   for (i=0; i<ismax; i++) {
2852:     ISRestoreIndices(isrow[i],irow+i);
2853:     if (!allcolumns[i]) {
2854:       ISRestoreIndices(iscol[i],icol+i);
2855:     }
2856:   }

2858:   for (i=0; i<ismax; i++) {
2859:     MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);
2860:     MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);
2861:   }

2863:   /* Destroy allocated memory */
2864:   PetscFree(sbuf_aa[0]);
2865:   PetscFree(sbuf_aa);
2866:   PetscFree5(*(PetscInt***)&irow,*(PetscInt***)&icol,nrow,ncol,issorted);

2868:   for (i=0; i<nrqs; ++i) {
2869:     PetscFree(rbuf4[i]);
2870:   }
2871:   PetscFree(rbuf4);

2873:   PetscFree4(row2proc,cmap,rmap,allcolumns);
2874:   return(0);
2875: }

2877: /*
2878:  Permute A & B into C's *local* index space using rowemb,dcolemb for A and rowemb,ocolemb for B.
2879:  Embeddings are supposed to be injections and the above implies that the range of rowemb is a subset
2880:  of [0,m), dcolemb is in [0,n) and ocolemb is in [N-n).
2881:  If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A&B.
2882:  After that B's columns are mapped into C's global column space, so that C is in the "disassembled"
2883:  state, and needs to be "assembled" later by compressing B's column space.

2885:  This function may be called in lieu of preallocation, so C should not be expected to be preallocated.
2886:  Following this call, C->A & C->B have been created, even if empty.
2887:  */
2888: PetscErrorCode MatSetSeqMats_MPIAIJ(Mat C,IS rowemb,IS dcolemb,IS ocolemb,MatStructure pattern,Mat A,Mat B)
2889: {
2890:   /* If making this function public, change the error returned in this function away from _PLIB. */
2892:   Mat_MPIAIJ     *aij;
2893:   Mat_SeqAIJ     *Baij;
2894:   PetscBool      seqaij,Bdisassembled;
2895:   PetscInt       m,n,*nz,i,j,ngcol,col,rstart,rend,shift,count;
2896:   PetscScalar    v;
2897:   const PetscInt *rowindices,*colindices;

2900:   /* Check to make sure the component matrices (and embeddings) are compatible with C. */
2901:   if (A) {
2902:     PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij);
2903:     if (!seqaij) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diagonal matrix is of wrong type");
2904:     if (rowemb) {
2905:       ISGetLocalSize(rowemb,&m);
2906:       if (m != A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Row IS of size %D is incompatible with diag matrix row size %D",m,A->rmap->n);
2907:     } else {
2908:       if (C->rmap->n != A->rmap->n) {
2909:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diag seq matrix is row-incompatible with the MPIAIJ matrix");
2910:       }
2911:     }
2912:     if (dcolemb) {
2913:       ISGetLocalSize(dcolemb,&n);
2914:       if (n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diag col IS of size %D is incompatible with diag matrix col size %D",n,A->cmap->n);
2915:     } else {
2916:       if (C->cmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diag seq matrix is col-incompatible with the MPIAIJ matrix");
2917:     }
2918:   }
2919:   if (B) {
2920:     PetscObjectBaseTypeCompare((PetscObject)B,MATSEQAIJ,&seqaij);
2921:     if (!seqaij) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Off-diagonal matrix is of wrong type");
2922:     if (rowemb) {
2923:       ISGetLocalSize(rowemb,&m);
2924:       if (m != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Row IS of size %D is incompatible with off-diag matrix row size %D",m,A->rmap->n);
2925:     } else {
2926:       if (C->rmap->n != B->rmap->n) {
2927:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Off-diag seq matrix is row-incompatible with the MPIAIJ matrix");
2928:       }
2929:     }
2930:     if (ocolemb) {
2931:       ISGetLocalSize(ocolemb,&n);
2932:       if (n != B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Off-diag col IS of size %D is incompatible with off-diag matrix col size %D",n,B->cmap->n);
2933:     } else {
2934:       if (C->cmap->N - C->cmap->n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Off-diag seq matrix is col-incompatible with the MPIAIJ matrix");
2935:     }
2936:   }

2938:   aij = (Mat_MPIAIJ*)C->data;
2939:   if (!aij->A) {
2940:     /* Mimic parts of MatMPIAIJSetPreallocation() */
2941:     MatCreate(PETSC_COMM_SELF,&aij->A);
2942:     MatSetSizes(aij->A,C->rmap->n,C->cmap->n,C->rmap->n,C->cmap->n);
2943:     MatSetBlockSizesFromMats(aij->A,C,C);
2944:     MatSetType(aij->A,MATSEQAIJ);
2945:     PetscLogObjectParent((PetscObject)C,(PetscObject)aij->A);
2946:   }
2947:   if (A) {
2948:     MatSetSeqMat_SeqAIJ(aij->A,rowemb,dcolemb,pattern,A);
2949:   } else {
2950:     MatSetUp(aij->A);
2951:   }
2952:   if (B) { /* Destroy the old matrix or the column map, depending on the sparsity pattern. */
2953:     /*
2954:       If pattern == DIFFERENT_NONZERO_PATTERN, we reallocate B and
2955:       need to "disassemble" B -- convert it to using C's global indices.
2956:       To insert the values we take the safer, albeit more expensive, route of MatSetValues().

2958:       If pattern == SUBSET_NONZERO_PATTERN, we do not "disassemble" B and do not reallocate;
2959:       we MatZeroValues(B) first, so there may be a bunch of zeros that, perhaps, could be compacted out.

2961:       TODO: Put B's values into aij->B's aij structure in place using the embedding ISs?
2962:       At least avoid calling MatSetValues() and the implied searches?
2963:     */

2965:     if (B && pattern == DIFFERENT_NONZERO_PATTERN) {
2966: #if defined(PETSC_USE_CTABLE)
2967:       PetscTableDestroy(&aij->colmap);
2968: #else
2969:       PetscFree(aij->colmap);
2970:       /* A bit of a HACK: ideally we should deal with case aij->B all in one code block below. */
2971:       if (aij->B) {
2972:         PetscLogObjectMemory((PetscObject)C,-aij->B->cmap->n*sizeof(PetscInt));
2973:       }
2974: #endif
2975:       ngcol = 0;
2976:       if (aij->lvec) {
2977:         VecGetSize(aij->lvec,&ngcol);
2978:       }
2979:       if (aij->garray) {
2980:         PetscFree(aij->garray);
2981:         PetscLogObjectMemory((PetscObject)C,-ngcol*sizeof(PetscInt));
2982:       }
2983:       VecDestroy(&aij->lvec);
2984:       VecScatterDestroy(&aij->Mvctx);
2985:     }
2986:     if (aij->B && B && pattern == DIFFERENT_NONZERO_PATTERN) {
2987:       MatDestroy(&aij->B);
2988:     }
2989:     if (aij->B && B && pattern == SUBSET_NONZERO_PATTERN) {
2990:       MatZeroEntries(aij->B);
2991:     }
2992:   }
2993:   Bdisassembled = PETSC_FALSE;
2994:   if (!aij->B) {
2995:     MatCreate(PETSC_COMM_SELF,&aij->B);
2996:     PetscLogObjectParent((PetscObject)C,(PetscObject)aij->B);
2997:     MatSetSizes(aij->B,C->rmap->n,C->cmap->N,C->rmap->n,C->cmap->N);
2998:     MatSetBlockSizesFromMats(aij->B,B,B);
2999:     MatSetType(aij->B,MATSEQAIJ);
3000:     Bdisassembled = PETSC_TRUE;
3001:   }
3002:   if (B) {
3003:     Baij = (Mat_SeqAIJ*)B->data;
3004:     if (pattern == DIFFERENT_NONZERO_PATTERN) {
3005:       PetscMalloc1(B->rmap->n,&nz);
3006:       for (i=0; i<B->rmap->n; i++) {
3007:         nz[i] = Baij->i[i+1] - Baij->i[i];
3008:       }
3009:       MatSeqAIJSetPreallocation(aij->B,0,nz);
3010:       PetscFree(nz);
3011:     }

3013:     PetscLayoutGetRange(C->rmap,&rstart,&rend);
3014:     shift = rend-rstart;
3015:     count = 0;
3016:     rowindices = NULL;
3017:     colindices = NULL;
3018:     if (rowemb) {
3019:       ISGetIndices(rowemb,&rowindices);
3020:     }
3021:     if (ocolemb) {
3022:       ISGetIndices(ocolemb,&colindices);
3023:     }
3024:     for (i=0; i<B->rmap->n; i++) {
3025:       PetscInt row;
3026:       row = i;
3027:       if (rowindices) row = rowindices[i];
3028:       for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
3029:         col  = Baij->j[count];
3030:         if (colindices) col = colindices[col];
3031:         if (Bdisassembled && col>=rstart) col += shift;
3032:         v    = Baij->a[count];
3033:         MatSetValues(aij->B,1,&row,1,&col,&v,INSERT_VALUES);
3034:         ++count;
3035:       }
3036:     }
3037:     /* No assembly for aij->B is necessary. */
3038:     /* FIXME: set aij->B's nonzerostate correctly. */
3039:   } else {
3040:     MatSetUp(aij->B);
3041:   }
3042:   C->preallocated  = PETSC_TRUE;
3043:   C->was_assembled = PETSC_FALSE;
3044:   C->assembled     = PETSC_FALSE;
3045:    /*
3046:       C will need to be assembled so that aij->B can be compressed into local form in MatSetUpMultiply_MPIAIJ().
3047:       Furthermore, its nonzerostate will need to be based on that of aij->A's and aij->B's.
3048:    */
3049:   return(0);
3050: }

3052: /*
3053:   B uses local indices with column indices ranging between 0 and N-n; they  must be interpreted using garray.
3054:  */
3055: PetscErrorCode MatGetSeqMats_MPIAIJ(Mat C,Mat *A,Mat *B)
3056: {
3057:   Mat_MPIAIJ *aij = (Mat_MPIAIJ*)C->data;

3062:   /* FIXME: make sure C is assembled */
3063:   *A = aij->A;
3064:   *B = aij->B;
3065:   /* Note that we don't incref *A and *B, so be careful! */
3066:   return(0);
3067: }

3069: /*
3070:   Extract MPI submatrices encoded by pairs of IS that may live on subcomms of C.
3071:   NOT SCALABLE due to the use of ISGetNonlocalIS() (see below).
3072: */
3073: PetscErrorCode MatCreateSubMatricesMPI_MPIXAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[],
3074:                                                  PetscErrorCode(*getsubmats_seq)(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat**),
3075:                                                  PetscErrorCode(*getlocalmats)(Mat,Mat*,Mat*),
3076:                                                  PetscErrorCode(*setseqmat)(Mat,IS,IS,MatStructure,Mat),
3077:                                                  PetscErrorCode(*setseqmats)(Mat,IS,IS,IS,MatStructure,Mat,Mat))
3078: {
3080:   PetscMPIInt    isize,flag;
3081:   PetscInt       i,ii,cismax,ispar;
3082:   Mat            *A,*B;
3083:   IS             *isrow_p,*iscol_p,*cisrow,*ciscol,*ciscol_p;

3086:   if (!ismax) return(0);

3088:   for (i = 0, cismax = 0; i < ismax; ++i) {
3089:     PetscMPIInt isize;
3090:     MPI_Comm_compare(((PetscObject)isrow[i])->comm,((PetscObject)iscol[i])->comm,&flag);
3091:     if (flag != MPI_IDENT) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Row and column index sets must have the same communicator");
3092:     MPI_Comm_size(((PetscObject)isrow[i])->comm, &isize);
3093:     if (isize > 1) ++cismax;
3094:   }

3096:   /*
3097:      If cismax is zero on all C's ranks, then and only then can we use purely sequential matrix extraction.
3098:      ispar counts the number of parallel ISs across C's comm.
3099:   */
3100:   MPIU_Allreduce(&cismax,&ispar,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));
3101:   if (!ispar) { /* Sequential ISs only across C's comm, so can call the sequential matrix extraction subroutine. */
3102:     (*getsubmats_seq)(C,ismax,isrow,iscol,scall,submat);
3103:     return(0);
3104:   }

3106:   /* if (ispar) */
3107:   /*
3108:     Construct the "complements" -- the off-processor indices -- of the iscol ISs for parallel ISs only.
3109:     These are used to extract the off-diag portion of the resulting parallel matrix.
3110:     The row IS for the off-diag portion is the same as for the diag portion,
3111:     so we merely alias (without increfing) the row IS, while skipping those that are sequential.
3112:   */
3113:   PetscMalloc2(cismax,&cisrow,cismax,&ciscol);
3114:   PetscMalloc1(cismax,&ciscol_p);
3115:   for (i = 0, ii = 0; i < ismax; ++i) {
3116:     MPI_Comm_size(((PetscObject)isrow[i])->comm,&isize);
3117:     if (isize > 1) {
3118:       /*
3119:          TODO: This is the part that's ***NOT SCALABLE***.
3120:          To fix this we need to extract just the indices of C's nonzero columns
3121:          that lie on the intersection of isrow[i] and ciscol[ii] -- the nonlocal
3122:          part of iscol[i] -- without actually computing ciscol[ii]. This also has
3123:          to be done without serializing on the IS list, so, most likely, it is best
3124:          done by rewriting MatCreateSubMatrices_MPIAIJ() directly.
3125:       */
3126:       ISGetNonlocalIS(iscol[i],&(ciscol[ii]));
3127:       /* Now we have to
3128:          (a) make sure ciscol[ii] is sorted, since, even if the off-proc indices
3129:              were sorted on each rank, concatenated they might no longer be sorted;
3130:          (b) Use ISSortPermutation() to construct ciscol_p, the mapping from the
3131:              indices in the nondecreasing order to the original index positions.
3132:          If ciscol[ii] is strictly increasing, the permutation IS is NULL.
3133:       */
3134:       ISSortPermutation(ciscol[ii],PETSC_FALSE,ciscol_p+ii);
3135:       ISSort(ciscol[ii]);
3136:       ++ii;
3137:     }
3138:   }
3139:   PetscMalloc2(ismax,&isrow_p,ismax,&iscol_p);
3140:   for (i = 0, ii = 0; i < ismax; ++i) {
3141:     PetscInt       j,issize;
3142:     const PetscInt *indices;

3144:     /*
3145:        Permute the indices into a nondecreasing order. Reject row and col indices with duplicates.
3146:      */
3147:     ISSortPermutation(isrow[i],PETSC_FALSE,isrow_p+i);
3148:     ISSort(isrow[i]);
3149:     ISGetLocalSize(isrow[i],&issize);
3150:     ISGetIndices(isrow[i],&indices);
3151:     for (j = 1; j < issize; ++j) {
3152:       if (indices[j] == indices[j-1]) {
3153:         SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Repeated indices in row IS %D: indices at %D and %D are both %D",i,j-1,j,indices[j]);
3154:       }
3155:     }
3156:     ISRestoreIndices(isrow[i],&indices);


3159:     ISSortPermutation(iscol[i],PETSC_FALSE,iscol_p+i);
3160:     ISSort(iscol[i]);
3161:     ISGetLocalSize(iscol[i],&issize);
3162:     ISGetIndices(iscol[i],&indices);
3163:     for (j = 1; j < issize; ++j) {
3164:       if (indices[j-1] == indices[j]) {
3165:         SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Repeated indices in col IS %D: indices at %D and %D are both %D",i,j-1,j,indices[j]);
3166:       }
3167:     }
3168:     ISRestoreIndices(iscol[i],&indices);
3169:     MPI_Comm_size(((PetscObject)isrow[i])->comm,&isize);
3170:     if (isize > 1) {
3171:       cisrow[ii] = isrow[i];
3172:       ++ii;
3173:     }
3174:   }
3175:   /*
3176:     Allocate the necessary arrays to hold the resulting parallel matrices as well as the intermediate
3177:     array of sequential matrices underlying the resulting parallel matrices.
3178:     Which arrays to allocate is based on the value of MatReuse scall and whether ISs are sorted and/or
3179:     contain duplicates.

3181:     There are as many diag matrices as there are original index sets. There are only as many parallel
3182:     and off-diag matrices, as there are parallel (comm size > 1) index sets.

3184:     ARRAYS that can hold Seq matrices get allocated in any event -- either here or by getsubmats_seq():
3185:     - If the array of MPI matrices already exists and is being reused, we need to allocate the array
3186:       and extract the underlying seq matrices into it to serve as placeholders, into which getsubmats_seq
3187:       will deposite the extracted diag and off-diag parts. Thus, we allocate the A&B arrays and fill them
3188:       with A[i] and B[ii] extracted from the corresponding MPI submat.
3189:     - However, if the rows, A's column indices or B's column indices are not sorted, the extracted A[i] & B[ii]
3190:       will have a different order from what getsubmats_seq expects.  To handle this case -- indicated
3191:       by a nonzero isrow_p[i], iscol_p[i], or ciscol_p[ii] -- we duplicate A[i] --> AA[i], B[ii] --> BB[ii]
3192:       (retrieve composed AA[i] or BB[ii]) and reuse them here. AA[i] and BB[ii] are then used to permute its
3193:       values into A[i] and B[ii] sitting inside the corresponding submat.
3194:     - If no reuse is taking place then getsubmats_seq will allocate the A&B arrays and create the corresponding
3195:       A[i], B[ii], AA[i] or BB[ii] matrices.
3196:   */
3197:   /* Parallel matrix array is allocated here only if no reuse is taking place. If reused, it is passed in by the caller. */
3198:   if (scall == MAT_INITIAL_MATRIX) {
3199:     PetscMalloc1(ismax,submat);
3200:   }

3202:   /* Now obtain the sequential A and B submatrices separately. */
3203:   /* scall=MAT_REUSE_MATRIX is not handled yet, because getsubmats_seq() requires reuse of A and B */
3204:   (*getsubmats_seq)(C,ismax,isrow,iscol,MAT_INITIAL_MATRIX,&A);
3205:   (*getsubmats_seq)(C,cismax,cisrow,ciscol,MAT_INITIAL_MATRIX,&B);

3207:   /*
3208:     If scall == MAT_REUSE_MATRIX AND the permutations are NULL, we are done, since the sequential
3209:     matrices A & B have been extracted directly into the parallel matrices containing them, or
3210:     simply into the sequential matrix identical with the corresponding A (if isize == 1).
3211:     Note that in that case colmap doesn't need to be rebuilt, since the matrices are expected
3212:     to have the same sparsity pattern.
3213:     Otherwise, A and/or B have to be properly embedded into C's index spaces and the correct colmap
3214:     must be constructed for C. This is done by setseqmat(s).
3215:   */
3216:   for (i = 0, ii = 0; i < ismax; ++i) {
3217:     /*
3218:        TODO: cache ciscol, permutation ISs and maybe cisrow? What about isrow & iscol?
3219:        That way we can avoid sorting and computing permutations when reusing.
3220:        To this end:
3221:         - remove the old cache, if it exists, when extracting submatrices with MAT_INITIAL_MATRIX
3222:         - if caching arrays to hold the ISs, make and compose a container for them so that it can
3223:           be destroyed upon destruction of C (use PetscContainerUserDestroy() to clear out the contents).
3224:     */
3225:     MatStructure pattern;
3226:     pattern = DIFFERENT_NONZERO_PATTERN;

3228:     MPI_Comm_size(((PetscObject)isrow[i])->comm,&isize);
3229:     /* Construct submat[i] from the Seq pieces A (and B, if necessary). */
3230:     if (isize > 1) {
3231:       if (scall == MAT_INITIAL_MATRIX) {
3232:         MatCreate(((PetscObject)isrow[i])->comm,(*submat)+i);
3233:         MatSetSizes((*submat)[i],A[i]->rmap->n,A[i]->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
3234:         MatSetType((*submat)[i],MATMPIAIJ);
3235:         PetscLayoutSetUp((*submat)[i]->rmap);
3236:         PetscLayoutSetUp((*submat)[i]->cmap);
3237:       }
3238:       /*
3239:         For each parallel isrow[i], insert the extracted sequential matrices into the parallel matrix.
3240:       */
3241:       {
3242:         Mat AA,BB;
3243:         AA = A[i];
3244:         BB = B[ii];
3245:         if (AA || BB) {
3246:           setseqmats((*submat)[i],isrow_p[i],iscol_p[i],ciscol_p[ii],pattern,AA,BB);
3247:           MatAssemblyBegin((*submat)[i],MAT_FINAL_ASSEMBLY);
3248:           MatAssemblyEnd((*submat)[i],MAT_FINAL_ASSEMBLY);
3249:         }

3251:         MatDestroy(&AA);
3252:       }
3253:       ISDestroy(ciscol+ii);
3254:       ISDestroy(ciscol_p+ii);
3255:       ++ii;
3256:     } else { /* if (isize == 1) */
3257:       if (scall == MAT_REUSE_MATRIX) {
3258:         MatDestroy(&(*submat)[i]);
3259:       }
3260:       if (isrow_p[i] || iscol_p[i]) {
3261:         MatDuplicate(A[i],MAT_DO_NOT_COPY_VALUES,(*submat)+i);
3262:         setseqmat((*submat)[i],isrow_p[i],iscol_p[i],pattern,A[i]);
3263:         /* Otherwise A is extracted straight into (*submats)[i]. */
3264:         /* TODO: Compose A[i] on (*submat([i] for future use, if ((isrow_p[i] || iscol_p[i]) && MAT_INITIAL_MATRIX). */
3265:         MatDestroy(A+i);
3266:       } else (*submat)[i] = A[i];
3267:     }
3268:     ISDestroy(&isrow_p[i]);
3269:     ISDestroy(&iscol_p[i]);
3270:   }
3271:   PetscFree2(cisrow,ciscol);
3272:   PetscFree2(isrow_p,iscol_p);
3273:   PetscFree(ciscol_p);
3274:   PetscFree(A);
3275:   MatDestroySubMatrices(cismax,&B);
3276:   return(0);
3277: }

3279: PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
3280: {

3284:   MatCreateSubMatricesMPI_MPIXAIJ(C,ismax,isrow,iscol,scall,submat,MatCreateSubMatrices_MPIAIJ,MatGetSeqMats_MPIAIJ,MatSetSeqMat_SeqAIJ,MatSetSeqMats_MPIAIJ);
3285:   return(0);
3286: }